A test for internal interpolation routines.
#include <stdlib.h>
#include <limits.h>
#include <math.h>
#include <mpi.h>
#include "tests.h"
#include "test_common.h"
#define UNSET (1337.0)
#define FRAC_MASK_VALUE (-1337.0)
#define SCALE(X) \
((X) * scaling_factors[scaling_factor_idx] + \
scaling_summand[scaling_summand_idx])
static int check_interpolation(
double * src_data_raw, size_t num_src_fields, size_t src_field_size,
double * ref_tgt_data_raw, size_t tgt_field_size, size_t collection_size);
static int check_interpolation_frac(
double * src_data_raw, double * src_frac_mask_raw,
size_t num_src_fields, size_t src_field_size,
double * ref_tgt_data_raw, size_t tgt_field_size, size_t collection_size);
MPI_Init(NULL, NULL);
xt_initialize(MPI_COMM_WORLD);
int comm_rank, comm_size;
MPI_Comm_rank(MPI_COMM_WORLD, &comm_rank);
MPI_Comm_size(MPI_COMM_WORLD, &comm_size);
MPI_Barrier(MPI_COMM_WORLD);
double scaling_factors[] = {1.0, 0.1, -1.0, 0.5, 10.0};
double scaling_summand[] = {0.0, -1.0, 1.0, 10.0};
enum {
num_procs = 4,
NUM_SCALING_FACTORS = sizeof(scaling_factors) / sizeof(scaling_factors[0]),
NUM_SCALING_SUMMAND = sizeof(scaling_summand) / sizeof(scaling_summand[0]),
};
if (comm_size != num_procs) {
PUT_ERR("ERROR: wrong number of processes");
xt_finalize();
MPI_Finalize();
return TEST_EXIT_CODE;
}
{
enum {num_fixed_values = 2};
double value[num_fixed_values] = {-1.0, 1.0};
size_t count[num_fixed_values] = {8, 8};
size_t pos[num_fixed_values][8] = {{0, 2, 4, 6, 8, 10, 12, 14},
{1, 3, 5, 7, 9, 11, 13, 15}};
for (int scaling_factor_idx = 0;
scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
for (int scaling_summand_idx = 0;
scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
scaling_factors[scaling_factor_idx],
scaling_summand[scaling_summand_idx]);
for (int i = 0; i < num_fixed_values; ++i)
interp,
value[i], count[i], pos[i]);
enum {num_src_fields = 1};
enum {src_field_size = 0};
double * src_data_raw = NULL;
enum {tgt_field_size = 16};
{{-1,1,-1,1,-1,1,-1,1,-1,1,-1,1,-1,1,-1,1}};
if (check_interpolation(
interp, src_data_raw, num_src_fields, src_field_size,
PUT_ERR("ERROR in yac_interpolation_add_fixed");
}
}
}
{
Xt_int src_indices[num_procs][4] =
{{0,1,2,3},{4,5,6,7},{8,9,10,11},{12,13,14,15}};
size_t num_src_indices[num_procs] = {4,4,0,0};
Xt_int dst_indices[num_procs][4] = {{-1},{0,2},{0,4,2,6},{-1}};
size_t num_dst_indices[num_procs] = {0,2,4,0};
Xt_idxlist src_idxlist =
xt_idxvec_new(src_indices[comm_rank], num_src_indices[comm_rank]);
Xt_idxlist dst_idxlist =
xt_idxvec_new(dst_indices[comm_rank], num_dst_indices[comm_rank]);
Xt_xmap xmap
= xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
Xt_redist redist = xt_redist_p2p_new(xmap, MPI_DOUBLE);
xt_xmap_delete(xmap);
xt_idxlist_delete(dst_idxlist);
xt_idxlist_delete(src_idxlist);
for (int scaling_factor_idx = 0;
scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
for (int scaling_summand_idx = 0;
scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
scaling_factors[scaling_factor_idx],
scaling_summand[scaling_summand_idx]);
enum {num_src_fields = 1};
size_t const src_field_size[num_procs] = {4,4,0,0};
double src_data_raw[num_procs][num_src_fields][4] =
{{{0,1,2,3}},{{4,5,6,7}},{{-1}},{{-1}}};
size_t const tgt_field_size[num_procs] = {0,4,4,0};
{{{UNSET}},
{{SCALE(0),SCALE(2),UNSET,UNSET}},
{{SCALE(0),SCALE(4),SCALE(2),SCALE(6)}},
{{UNSET}}};
if (check_interpolation(
interp,
(src_field_size[comm_rank] > 0)?
(&src_data_raw[comm_rank][0][0]):NULL,
num_src_fields, src_field_size[comm_rank],
(tgt_field_size[comm_rank] > 0)?
(&ref_tgt_data_raw[comm_rank][0][0]):NULL,
PUT_ERR("ERROR in yac_interpolation_add_direct");
}
}
xt_redist_delete(redist);
}
{
Xt_int src_indices[num_procs][4] =
{{0,1,2,3},{4,5,6,7},{8,9,10,11},{12,13,14,15}};
size_t num_src_indices[num_procs] = {4,4,0,0};
Xt_int dst_indices[num_procs][4] = {{-1},{0,2},{0,4,2,6},{-1}};
size_t num_dst_indices[num_procs] = {0,2,4,0};
Xt_idxlist src_idxlist =
xt_idxvec_new(src_indices[comm_rank], num_src_indices[comm_rank]);
Xt_idxlist dst_idxlist =
xt_idxvec_new(dst_indices[comm_rank], num_dst_indices[comm_rank]);
Xt_xmap xmap
= xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
Xt_redist redist = xt_redist_p2p_new(xmap, MPI_DOUBLE);
xt_xmap_delete(xmap);
xt_idxlist_delete(dst_idxlist);
xt_idxlist_delete(src_idxlist);
for (int scaling_factor_idx = 0;
scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
for (int scaling_summand_idx = 0;
scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
scaling_factors[scaling_factor_idx],
scaling_summand[scaling_summand_idx]);
enum {num_src_fields = 1};
size_t const src_field_size[num_procs] = {4,4,0,0};
double src_data_raw[num_procs][num_src_fields][4] =
{{{0,1,2,3}},{{4,5,6,7}},{{-1}},{{-1}}};
double src_frac_mask_raw[num_procs][num_src_fields][4] =
{{{1.0,1.0,0.0,1.0}},{{0.5,0.5,0.5,0.5}},{{-1}},{{-1}}};
size_t const tgt_field_size[num_procs] = {0,4,4,0};
{{{UNSET}},
{{SCALE(0),FRAC_MASK_VALUE,UNSET,UNSET}},
{{SCALE(0),SCALE(4),FRAC_MASK_VALUE,SCALE(6)}},
{{UNSET}}};
if (check_interpolation_frac(
interp,
(src_field_size[comm_rank] > 0)?
(&src_data_raw[comm_rank][0][0]):NULL,
(src_field_size[comm_rank] > 0)?
(&src_frac_mask_raw[comm_rank][0][0]):NULL,
num_src_fields, src_field_size[comm_rank],
(tgt_field_size[comm_rank] > 0)?
(&ref_tgt_data_raw[comm_rank][0][0]):NULL,
PUT_ERR("ERROR in yac_interpolation_add_direct");
}
}
xt_redist_delete(redist);
}
{
enum {num_src_fields = 2};
enum {num_src_indices = 4};
Xt_int src_indices[num_procs][num_src_indices] =
{{0,1,2,3},{4,5,6,7},{8,9,10,11},{12,13,14,15}};
Xt_int dst_indices[num_procs][2][2] = {{{-1},{-1}},
{{0,2},{1,3}},
{{1,4},{2,5}},
{{-1},{-1}}};
int src_offsets[num_procs] = {0,1,2,3};
int dst_offsets[num_procs][num_src_fields][2] = {{{-1},{-1}},
{{0,2},{3,1}},
{{0,1},{2,3}},
{{-1},{-1}}};
size_t num_dst_indices[num_procs][num_src_fields] =
{{0,0},{2,2},{2,2},{0,0}};
Xt_redist redists[num_src_fields];
for (int i = 0; i < num_src_fields; ++i) {
Xt_idxlist src_idxlist =
xt_idxvec_new(src_indices[comm_rank], num_src_indices);
Xt_idxlist dst_idxlist =
xt_idxvec_new(dst_indices[comm_rank][i], num_dst_indices[comm_rank][i]);
Xt_xmap xmap =
xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
redists[i] =
xt_redist_p2p_off_new(
xmap, src_offsets, dst_offsets[comm_rank][i], MPI_DOUBLE);
xt_xmap_delete(xmap);
xt_idxlist_delete(dst_idxlist);
xt_idxlist_delete(src_idxlist);
}
for (int scaling_factor_idx = 0;
scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
for (int scaling_summand_idx = 0;
scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
scaling_factors[scaling_factor_idx],
scaling_summand[scaling_summand_idx]);
size_t const src_field_size[num_procs] = {4,4,4,4};
double src_data_raw[num_procs][num_src_fields][num_src_indices] =
{{{0,1,2,3},{00,10,20,30}},
{{4,5,6,7},{40,50,60,70}},
{{8,9,10,11},{80,90,100,110}},
{{12,13,14,15},{120,130,140,150}}};
size_t const tgt_field_size[num_procs] = {0,4,4,0};
{{{UNSET}},
{{SCALE(0),SCALE(30),SCALE(2),SCALE(10)}},
{{SCALE(1),SCALE(4),SCALE(20),SCALE(50)}},
{{UNSET}}};
if (check_interpolation(
interp,
(src_field_size[comm_rank] > 0)?
(&src_data_raw[comm_rank][0][0]):NULL,
num_src_fields, src_field_size[comm_rank],
(tgt_field_size[comm_rank] > 0)?
(&ref_tgt_data_raw[comm_rank][0][0]):NULL,
PUT_ERR("ERROR in yac_interpolation_add_direct_mf");
}
}
xt_redist_delete(redists[1]);
xt_redist_delete(redists[0]);
}
{
enum {num_src_fields = 2};
enum {num_src_indices = 4};
Xt_int src_indices[num_procs][num_src_indices] =
{{0,1,2,3},{4,5,6,7},{8,9,10,11},{12,13,14,15}};
Xt_int dst_indices[num_procs][2][2] = {{{-1},{-1}},
{{0,2},{1,3}},
{{1,4},{2,5}},
{{-1},{-1}}};
int src_offsets[num_procs] = {0,1,2,3};
int dst_offsets[num_procs][num_src_fields][2] = {{{-1},{-1}},
{{0,2},{3,1}},
{{0,1},{2,3}},
{{-1},{-1}}};
size_t num_dst_indices[num_procs][num_src_fields] =
{{0,0},{2,2},{2,2},{0,0}};
Xt_redist redists[num_src_fields];
for (int i = 0; i < num_src_fields; ++i) {
Xt_idxlist src_idxlist =
xt_idxvec_new(src_indices[comm_rank], num_src_indices);
Xt_idxlist dst_idxlist =
xt_idxvec_new(dst_indices[comm_rank][i], num_dst_indices[comm_rank][i]);
Xt_xmap xmap =
xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
redists[i] =
xt_redist_p2p_off_new(
xmap, src_offsets, dst_offsets[comm_rank][i], MPI_DOUBLE);
xt_xmap_delete(xmap);
xt_idxlist_delete(dst_idxlist);
xt_idxlist_delete(src_idxlist);
}
for (int scaling_factor_idx = 0;
scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
for (int scaling_summand_idx = 0;
scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
scaling_factors[scaling_factor_idx],
scaling_summand[scaling_summand_idx]);
size_t const src_field_size[num_procs] = {4,4,4,4};
double src_data_raw[num_procs][num_src_fields][num_src_indices] =
{{{0,1,2,3},{00,10,20,30}},
{{4,5,6,7},{40,50,60,70}},
{{8,9,10,11},{80,90,100,110}},
{{12,13,14,15},{120,130,140,150}}};
double src_frac_mask_raw[num_procs][num_src_fields][num_src_indices] =
{{{1.0,0.0,1.0,1.0},{1.0,1.0,1.0,1.0}},
{{0.0,1.0,1.0,1.0},{1.0,0.0,1.0,1.0}},
{{1.0,1.0,1.0,1.0},{1.0,1.0,1.0,1.0}},
{{1.0,1.0,1.0,1.0},{1.0,1.0,1.0,1.0}}};
size_t const tgt_field_size[num_procs] = {0,4,4,0};
{{{UNSET}},
{{SCALE(0),SCALE(30),SCALE(2),SCALE(10)}},
{{FRAC_MASK_VALUE,FRAC_MASK_VALUE,SCALE(20),FRAC_MASK_VALUE}},
{{UNSET}}};
if (check_interpolation_frac(
interp,
(src_field_size[comm_rank] > 0)?
(&src_data_raw[comm_rank][0][0]):NULL,
(src_field_size[comm_rank] > 0)?
(&src_frac_mask_raw[comm_rank][0][0]):NULL,
num_src_fields, src_field_size[comm_rank],
(tgt_field_size[comm_rank] > 0)?
(&ref_tgt_data_raw[comm_rank][0][0]):NULL,
PUT_ERR("ERROR in yac_interpolation_add_direct_mf");
}
}
xt_redist_delete(redists[1]);
xt_redist_delete(redists[0]);
}
{
enum {num_src_fields = 2};
Xt_redist halo_redists[num_src_fields];
{
enum {num_src_indices = 4};
Xt_int src_indices[num_procs][num_src_fields][num_src_indices] =
{{{0,1,2,3},{16,17,18,19}},
{{4,5,6,7},{20,21,22,23}},
{{8,9,10,11},{24,25,26,27}},
{{12,13,14,15},{28,29,30,31}}};
Xt_int dst_indices[num_procs][num_src_fields][1] =
{{{7},{22}},
{{2},{19}},
{{-1},{-1}},
{{-1},{-1}}};
size_t num_dst_indices[num_procs][num_src_fields] =
{{1,1},{1,1},{0,0},{0,0}};
for (int field_id = 0; field_id < num_src_fields; ++field_id) {
Xt_idxlist src_idxlist =
xt_idxvec_new(src_indices[comm_rank][field_id], num_src_indices);
Xt_idxlist dst_idxlist =
xt_idxvec_new(
dst_indices[comm_rank][field_id],
num_dst_indices[comm_rank][field_id]);
Xt_xmap xmap =
xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
halo_redists[field_id] = xt_redist_p2p_new(xmap, MPI_DOUBLE);
xt_xmap_delete(xmap);
xt_idxlist_delete(dst_idxlist);
xt_idxlist_delete(src_idxlist);
}
}
size_t tgt_count[num_procs] = {2,5,0,0};
size_t num_src_per_tgt[num_procs][5] =
{{4,4},{4,4,1,2,2},{(size_t)-1},{(size_t)-1}};
double weights[num_procs][13] =
{{0.1,0.2,0.3,0.4,
1.4,1.5,1.6,1.7},
{0.5,0.6,0.7,0.8,
0.9,1.0,1.1,1.2,
1.3,
2.2,2.3,
2.4,2.5},{-1},{-1}};
size_t src_field_idx[num_procs][13] =
{{0,1,0,1,
0,1,SIZE_MAX,SIZE_MAX},
{SIZE_MAX,SIZE_MAX,1,0,
0,1,0,1,
1,
0,0,
1,1},{(size_t)-1},{(size_t)-1}};
size_t src_idx[num_procs][13] =
{{0,1,2,3,
2,3,1,0},
{0,1,0,1,
2,2,3,3,
3,
1,2,
1,2},{(size_t)-1},{(size_t)-1}};
Xt_redist result_redist;
{
Xt_int src_indices[num_procs][5] = {{4,8},{5,6,7,10,11},{-1},{-1}};
size_t num_src_indices[num_procs] = {2,5,0,0};
Xt_int dst_indices[num_procs][4] = {{-1},{4,5,6,7},{8,4,10,11},{-1}};
size_t num_dst_indices[num_procs] = {0,4,4,0};
Xt_idxlist src_idxlist =
xt_idxvec_new(src_indices[comm_rank], num_src_indices[comm_rank]);
Xt_idxlist dst_idxlist =
xt_idxvec_new(dst_indices[comm_rank], num_dst_indices[comm_rank]);
Xt_xmap xmap = xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
result_redist = xt_redist_p2p_new(xmap, MPI_DOUBLE);
xt_xmap_delete(xmap);
xt_idxlist_delete(dst_idxlist);
xt_idxlist_delete(src_idxlist);
}
for (int scaling_factor_idx = 0;
scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
for (int scaling_summand_idx = 0;
scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
for (int with_weights = 0; with_weights < 2; ++with_weights) {
scaling_factors[scaling_factor_idx],
scaling_summand[scaling_summand_idx]);
if (with_weights) {
interp, halo_redists, tgt_count[comm_rank],
num_src_per_tgt[comm_rank], weights[comm_rank],
src_field_idx[comm_rank], src_idx[comm_rank], num_src_fields,
result_redist);
} else {
interp, halo_redists, tgt_count[comm_rank],
num_src_per_tgt[comm_rank], src_field_idx[comm_rank],
src_idx[comm_rank], num_src_fields, result_redist);
}
enum {src_field_size = 4};
double src_data_raw[num_procs][2][src_field_size] =
{{{0,1,2,3}, {16,17,18,19}},
{{4,5,6,7}, {20,21,22,23}},
{{8,9,10,11}, {24,25,26,27}},
{{12,13,14,15},{28,29,30,31}}};
size_t const tgt_field_size[num_procs] = {0,4,4,0};
{{{{UNSET}},
{{SCALE(0+17+2+19),SCALE(2+19+20+5),SCALE(6+22+7+23),SCALE(23)}},
{{SCALE(2+19+22+7),SCALE(0+17+2+19),SCALE(5+6),SCALE(21+22)}},
{{UNSET}}},
{{{UNSET}},
{{SCALE(0.1*0+0.2*17+0.3*2+0.4*19),
SCALE(0.5*2+0.6*19+0.7*20+0.8*5),
SCALE(0.9*6+1.0*22+1.1*7+1.2*23),
SCALE(1.3*23)}},
{{SCALE(1.4*2+1.5*19+1.6*22+1.7*7),
SCALE(0.1*0+0.2*17+0.3*2+0.4*19),
SCALE(2.2*5+2.3*6),
SCALE(2.4*21+2.5*22)}},
{{UNSET}}}};
if (check_interpolation(
interp, &src_data_raw[comm_rank][0][0],
num_src_fields, src_field_size,
&ref_tgt_data_raw[with_weights][comm_rank][0][0],
PUT_ERR("ERROR in yac_interpolation_add_sum_at_src/"
"yac_interpolation_add_weight_sum_mvp_at_src");
}
}
}
xt_redist_delete(result_redist);
xt_redist_delete(halo_redists[1]);
xt_redist_delete(halo_redists[0]);
}
{
enum {num_src_fields = 2};
Xt_redist halo_redists[num_src_fields];
{
enum {num_src_indices = 4};
Xt_int src_indices[num_procs][num_src_fields][num_src_indices] =
{{{0,1,2,3},{16,17,18,19}},
{{4,5,6,7},{20,21,22,23}},
{{8,9,10,11},{24,25,26,27}},
{{12,13,14,15},{28,29,30,31}}};
Xt_int dst_indices[num_procs][num_src_fields][1] =
{{{7},{22}},
{{2},{19}},
{{-1},{-1}},
{{-1},{-1}}};
size_t num_dst_indices[num_procs][num_src_fields] =
{{1,1},{1,1},{0,0},{0,0}};
for (int field_id = 0; field_id < num_src_fields; ++field_id) {
Xt_idxlist src_idxlist =
xt_idxvec_new(src_indices[comm_rank][field_id], num_src_indices);
Xt_idxlist dst_idxlist =
xt_idxvec_new(
dst_indices[comm_rank][field_id],
num_dst_indices[comm_rank][field_id]);
Xt_xmap xmap =
xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
halo_redists[field_id] = xt_redist_p2p_new(xmap, MPI_DOUBLE);
xt_xmap_delete(xmap);
xt_idxlist_delete(dst_idxlist);
xt_idxlist_delete(src_idxlist);
}
}
size_t tgt_count[num_procs] = {2,5,0,0};
size_t num_src_per_tgt[num_procs][5] =
{{4,4},{4,4,1,2,2},{(size_t)-1},{(size_t)-1}};
double weights[num_procs][13] =
{{0.1,0.2,0.3,0.4,
1.4,1.5,1.6,1.7},
{0.5,0.6,0.7,0.8,
0.9,1.0,1.1,1.2,
1.3,
2.2,2.3,
2.4,2.5},{-1},{-1}};
size_t src_field_idx[num_procs][13] =
{{0,1,0,1,
0,1,SIZE_MAX,SIZE_MAX},
{SIZE_MAX,SIZE_MAX,1,0,
0,1,0,1,
1,
0,0,
1,1},{(size_t)-1},{(size_t)-1}};
size_t src_idx[num_procs][13] =
{{0,1,2,3,
2,3,1,0},
{0,1,0,1,
2,2,3,3,
3,
1,2,
1,2},{(size_t)-1},{(size_t)-1}};
Xt_redist result_redist;
{
Xt_int src_indices[num_procs][5] = {{4,8},{5,6,7,10,11},{-1},{-1}};
size_t num_src_indices[num_procs] = {2,5,0,0};
Xt_int dst_indices[num_procs][4] = {{-1},{4,5,6,7},{8,4,10,11},{-1}};
size_t num_dst_indices[num_procs] = {0,4,4,0};
Xt_idxlist src_idxlist =
xt_idxvec_new(src_indices[comm_rank], num_src_indices[comm_rank]);
Xt_idxlist dst_idxlist =
xt_idxvec_new(dst_indices[comm_rank], num_dst_indices[comm_rank]);
Xt_xmap xmap = xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
result_redist = xt_redist_p2p_new(xmap, MPI_DOUBLE);
xt_xmap_delete(xmap);
xt_idxlist_delete(dst_idxlist);
xt_idxlist_delete(src_idxlist);
}
for (int scaling_factor_idx = 0;
scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
for (int scaling_summand_idx = 0;
scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
for (int with_weights = 0; with_weights < 2; ++with_weights) {
scaling_factors[scaling_factor_idx],
scaling_summand[scaling_summand_idx]);
if (with_weights) {
interp, halo_redists, tgt_count[comm_rank],
num_src_per_tgt[comm_rank], weights[comm_rank],
src_field_idx[comm_rank], src_idx[comm_rank], num_src_fields,
result_redist);
} else {
interp, halo_redists, tgt_count[comm_rank],
num_src_per_tgt[comm_rank], src_field_idx[comm_rank],
src_idx[comm_rank], num_src_fields, result_redist);
}
enum {src_field_size = 4};
double src_data_raw[num_procs][2][src_field_size] =
{{{0,1,2,3}, {16,17,18,19}},
{{4,5,6,7}, {20,21,22,23}},
{{8,9,10,11}, {24,25,26,27}},
{{12,13,14,15},{28,29,30,31}}};
double src_frac_mask_raw[num_procs][2][src_field_size] =
{{{1.0, 1.0, 1.0, 1.0}, {1.0, 1.0, 1.0, 1.0}},
{{1.0, 0.5, 0.0, 0.0}, {0.5, 0.0, 0.0, 0.0}},
{{1.0, 1.0, 1.0, 1.0}, {1.0, 1.0, 1.0, 1.0}},
{{1.0, 1.0, 1.0, 1.0}, {1.0, 1.0, 1.0, 1.0}}};
size_t const tgt_field_size[num_procs] = {0,4,4,0};
{{{{UNSET}},
{{SCALE((1.0*0+1.0*17+1.0*2+1.0*19)/(1.0+1.0+1.0+1.0)),
SCALE((1.0*2+1.0*19+0.5*20+0.5*5)/(1.0+1.0+0.5+0.5)),
FRAC_MASK_VALUE,FRAC_MASK_VALUE}},
{{SCALE((1.0*2+1.0*19+0.0*22+0.0*7)/(1.0+1.0+0.0+0.0)),
SCALE((1.0*0+1.0*17+1.0*2+1.0*19)/(1.0+1.0+1.0+1.0)),
SCALE((0.5*5+0.0*6)/(0.5+0.0)),
FRAC_MASK_VALUE}},
{{UNSET}}},
{{{UNSET}},
{{SCALE((1.0*0.1*0+1.0*0.2*17+1.0*0.3*2+1.0*0.4*19)/
(1.0*0.1+1.0*0.2+1.0*0.3+1.0*0.4)),
SCALE((1.0*0.5*2+1.0*0.6*19+0.5*0.7*20+0.5*0.8*5)/
(1.0*0.5+1.0*0.6+0.5*0.7+0.5*0.8)),
FRAC_MASK_VALUE, FRAC_MASK_VALUE}},
{{SCALE((1.0*1.4*2+1.0*1.5*19+0.0*1.6*22+0.0*1.7*7)/
(1.0*1.4+1.0*1.5+0.0*1.6+0.0*1.7)),
SCALE((1.0*0.1*0+1.0*0.2*17+1.0*0.3*2+1.0*0.4*19)/
(1.0*0.1+1.0*0.2+1.0*0.3+1.0*0.4)),
SCALE((0.5*2.2*5+0.0*2.3*6)/(0.5*2.2+0.0*2.3)),
FRAC_MASK_VALUE}},
{{UNSET}}}};
if (check_interpolation_frac(
interp, &src_data_raw[comm_rank][0][0],
&src_frac_mask_raw[comm_rank][0][0],
num_src_fields, src_field_size,
&ref_tgt_data_raw[with_weights][comm_rank][0][0],
PUT_ERR("ERROR in yac_interpolation_add_sum_at_src/"
"yac_interpolation_add_weight_sum_mvp_at_src frac");
}
}
}
xt_redist_delete(result_redist);
xt_redist_delete(halo_redists[1]);
xt_redist_delete(halo_redists[0]);
}
{
enum {num_src_fields = 2};
Xt_redist src_redists[num_src_fields];
{
enum {num_src_indices = 4};
Xt_int src_indices[num_procs][num_src_fields][num_src_indices] =
{{{0,1,2,3},{16,17,18,19}},
{{4,5,6,7},{20,21,22,23}},
{{8,9,10,11},{24,25,26,27}},
{{12,13,14,15},{28,29,30,31}}};
Xt_int dst_indices[num_procs][num_src_fields][5] =
{{{-1},{-1}},
{{0,2},{17,19}},
{{0,2,5,6,7},{17,19,21,22}},
{{-1},{-1}}};
size_t num_dst_indices[num_procs][num_src_fields] =
{{0,0},{2,2},{5,4},{0,0}};
for (int field_id = 0; field_id < num_src_fields; ++field_id) {
Xt_idxlist src_idxlist =
xt_idxvec_new(src_indices[comm_rank][field_id], num_src_indices);
Xt_idxlist dst_idxlist =
xt_idxvec_new(
dst_indices[comm_rank][field_id],
num_dst_indices[comm_rank][field_id]);
Xt_xmap xmap =
xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
src_redists[field_id] = xt_redist_p2p_new(xmap, MPI_DOUBLE);
xt_xmap_delete(xmap);
xt_idxlist_delete(dst_idxlist);
xt_idxlist_delete(src_idxlist);
}
}
size_t tgt_pos[num_procs][4] =
{{(size_t)-1},{0,1,2,3},{0,1,2,3},{(size_t)-1}};
size_t tgt_count[num_procs] = {0,4,4,0};
size_t num_src_per_tgt[num_procs][4] =
{{(size_t)-1},{4,4,4,1},{4,4,2,2},{(size_t)-1}};
double weights[num_procs][13] =
{{-1},
{0.1,0.2,0.3,0.4,
0.5,0.6,0.7,0.8,
0.9,1.0,1.1,1.2,
1.3},
{1.4,1.5,1.6,1.7,
0.1,0.2,0.3,0.4,
2.2,2.3,
2.4,2.5},
{-1}};
size_t src_field_idx[num_procs][13] =
{{(size_t)-1},
{SIZE_MAX,SIZE_MAX,SIZE_MAX,SIZE_MAX,
SIZE_MAX,SIZE_MAX,1,0,
0,1,0,1,
1},
{SIZE_MAX,SIZE_MAX,SIZE_MAX,SIZE_MAX,
SIZE_MAX,SIZE_MAX,SIZE_MAX,SIZE_MAX,
SIZE_MAX,SIZE_MAX,
SIZE_MAX,SIZE_MAX},
{(size_t)-1}};
size_t src_idx[num_procs][13] =
{{(size_t)-1},
{0,2,1,3, 1,3,0,1, 2,2,3,3, 3},
{1,6,8,4, 0,5,1,6, 2,3, 7,8},
{(size_t)-1}};
for (int scaling_factor_idx = 0;
scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
for (int scaling_summand_idx = 0;
scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
for (int with_weights = 0; with_weights < 2; ++with_weights) {
scaling_factors[scaling_factor_idx],
scaling_summand[scaling_summand_idx]);
if (with_weights) {
interp, src_redists, tgt_pos[comm_rank], tgt_count[comm_rank],
num_src_per_tgt[comm_rank], weights[comm_rank],
src_field_idx[comm_rank], src_idx[comm_rank], num_src_fields);
} else {
interp, src_redists, tgt_pos[comm_rank], tgt_count[comm_rank],
num_src_per_tgt[comm_rank], src_field_idx[comm_rank],
src_idx[comm_rank], num_src_fields);
}
enum {src_field_size = 4};
double src_data_raw[num_procs][num_src_fields][src_field_size] =
{{{0,1,2,3}, {16,17,18,19}},
{{4,5,6,7}, {20,21,22,23}},
{{8,9,10,11}, {24,25,26,27}},
{{12,13,14,15},{28,29,30,31}}};
size_t const tgt_field_size[num_procs] = {0,4,4,0};
{{{{UNSET}},
{{SCALE(0+17+2+19),SCALE(2+19+20+5),SCALE(6+22+7+23),SCALE(23)}},
{{SCALE(2+19+22+7),SCALE(0+17+2+19),SCALE(5+6),SCALE(21+22)}},
{{UNSET}}},
{{{UNSET}},
{{SCALE(0.1*0+0.2*17+0.3*2+0.4*19),
SCALE(0.5*2+0.6*19+0.7*20+0.8*5),
SCALE(0.9*6+1.0*22+1.1*7+1.2*23),
SCALE(1.3*23)}},
{{SCALE(1.4*2+1.5*19+1.6*22+1.7*7),
SCALE(0.1*0+0.2*17+0.3*2+0.4*19),
SCALE(2.2*5+2.3*6),
SCALE(2.4*21+2.5*22)}},
{{UNSET}}}};
if (check_interpolation(
interp, &src_data_raw[comm_rank][0][0],
num_src_fields, src_field_size,
&ref_tgt_data_raw[with_weights][comm_rank][0][0],
PUT_ERR("ERROR in yac_interpolation_add_sum_at_tgt/"
"yac_interpolation_add_weight_sum_mvp_at_tgt");
}
}
}
xt_redist_delete(src_redists[1]);
xt_redist_delete(src_redists[0]);
}
{
enum {num_src_fields = 2};
Xt_redist src_redists[num_src_fields];
{
enum {num_src_indices = 4};
Xt_int src_indices[num_procs][num_src_fields][num_src_indices] =
{{{0,1,2,3},{16,17,18,19}},
{{4,5,6,7},{20,21,22,23}},
{{8,9,10,11},{24,25,26,27}},
{{12,13,14,15},{28,29,30,31}}};
Xt_int dst_indices[num_procs][num_src_fields][5] =
{{{-1},{-1}},
{{0,2},{17,19}},
{{0,2,5,6,7},{17,19,21,22}},
{{-1},{-1}}};
size_t num_dst_indices[num_procs][num_src_fields] =
{{0,0},{2,2},{5,4},{0,0}};
for (int field_id = 0; field_id < num_src_fields; ++field_id) {
Xt_idxlist src_idxlist =
xt_idxvec_new(src_indices[comm_rank][field_id], num_src_indices);
Xt_idxlist dst_idxlist =
xt_idxvec_new(
dst_indices[comm_rank][field_id],
num_dst_indices[comm_rank][field_id]);
Xt_xmap xmap =
xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
src_redists[field_id] = xt_redist_p2p_new(xmap, MPI_DOUBLE);
xt_xmap_delete(xmap);
xt_idxlist_delete(dst_idxlist);
xt_idxlist_delete(src_idxlist);
}
}
size_t tgt_pos[num_procs][4] =
{{(size_t)-1},{0,1,2,3},{0,1,2,3},{(size_t)-1}};
size_t tgt_count[num_procs] = {0,4,4,0};
size_t num_src_per_tgt[num_procs][4] =
{{(size_t)-1},{4,4,4,1},{4,4,2,2},{(size_t)-1}};
double weights[num_procs][13] =
{{-1},
{0.1,0.2,0.3,0.4,
0.5,0.6,0.7,0.8,
0.9,1.0,1.1,1.2,
1.3},
{1.4,1.5,1.6,1.7,
0.1,0.2,0.3,0.4,
2.2,2.3,
2.4,2.5},
{-1}};
size_t src_field_idx[num_procs][13] =
{{(size_t)-1},
{SIZE_MAX,SIZE_MAX,SIZE_MAX,SIZE_MAX,
SIZE_MAX,SIZE_MAX,1,0,
0,1,0,1,
1},
{SIZE_MAX,SIZE_MAX,SIZE_MAX,SIZE_MAX,
SIZE_MAX,SIZE_MAX,SIZE_MAX,SIZE_MAX,
SIZE_MAX,SIZE_MAX,
SIZE_MAX,SIZE_MAX},
{(size_t)-1}};
size_t src_idx[num_procs][13] =
{{(size_t)-1},
{0,2,1,3, 1,3,0,1, 2,2,3,3, 3},
{1,6,8,4, 0,5,1,6, 2,3, 7,8},
{(size_t)-1}};
for (int scaling_factor_idx = 0;
scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
for (int scaling_summand_idx = 0;
scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
for (int with_weights = 0; with_weights < 2; ++with_weights) {
scaling_factors[scaling_factor_idx],
scaling_summand[scaling_summand_idx]);
if (with_weights) {
interp, src_redists, tgt_pos[comm_rank], tgt_count[comm_rank],
num_src_per_tgt[comm_rank], weights[comm_rank],
src_field_idx[comm_rank], src_idx[comm_rank], num_src_fields);
} else {
interp, src_redists, tgt_pos[comm_rank], tgt_count[comm_rank],
num_src_per_tgt[comm_rank], src_field_idx[comm_rank],
src_idx[comm_rank], num_src_fields);
}
enum {src_field_size = 4};
double src_data_raw[num_procs][num_src_fields][src_field_size] =
{{{0,1,2,3}, {16,17,18,19}},
{{4,5,6,7}, {20,21,22,23}},
{{8,9,10,11}, {24,25,26,27}},
{{12,13,14,15},{28,29,30,31}}};
double src_frac_mask_raw[num_procs][2][src_field_size] =
{{{1.0, 1.0, 1.0, 1.0}, {1.0, 1.0, 1.0, 1.0}},
{{1.0, 0.5, 0.0, 0.0}, {0.5, 0.0, 0.0, 0.0}},
{{1.0, 1.0, 1.0, 1.0}, {1.0, 1.0, 1.0, 1.0}},
{{1.0, 1.0, 1.0, 1.0}, {1.0, 1.0, 1.0, 1.0}}};
size_t const tgt_field_size[num_procs] = {0,4,4,0};
{{{{UNSET}},
{{SCALE((1.0*0+1.0*17+1.0*2+1.0*19)/(1.0+1.0+1.0+1.0)),
SCALE((1.0*2+1.0*19+0.5*20+0.5*5)/(1.0+1.0+0.5+0.5)),
FRAC_MASK_VALUE,FRAC_MASK_VALUE}},
{{SCALE((1.0*2+1.0*19+0.0*22+0.0*7)/(1.0+1.0+0.0+0.0)),
SCALE((1.0*0+1.0*17+1.0*2+1.0*19)/(1.0+1.0+1.0+1.0)),
SCALE((0.5*5+0.0*6)/(0.5+0.0)),
FRAC_MASK_VALUE}},
{{UNSET}}},
{{{UNSET}},
{{SCALE((1.0*0.1*0+1.0*0.2*17+1.0*0.3*2+1.0*0.4*19)/
(1.0*0.1+1.0*0.2+1.0*0.3+1.0*0.4)),
SCALE((1.0*0.5*2+1.0*0.6*19+0.5*0.7*20+0.5*0.8*5)/
(1.0*0.5+1.0*0.6+0.5*0.7+0.5*0.8)),
FRAC_MASK_VALUE, FRAC_MASK_VALUE}},
{{SCALE((1.0*1.4*2+1.0*1.5*19+0.0*1.6*22+0.0*1.7*7)/
(1.0*1.4+1.0*1.5+0.0*1.6+0.0*1.7)),
SCALE((1.0*0.1*0+1.0*0.2*17+1.0*0.3*2+1.0*0.4*19)/
(1.0*0.1+1.0*0.2+1.0*0.3+1.0*0.4)),
SCALE((0.5*2.2*5+0.0*2.3*6)/(0.5*2.2+0.0*2.3)),
FRAC_MASK_VALUE}},
{{UNSET}}}};
if (check_interpolation_frac(
interp,
&src_data_raw[comm_rank][0][0],
&src_frac_mask_raw[comm_rank][0][0],
num_src_fields, src_field_size,
&ref_tgt_data_raw[with_weights][comm_rank][0][0],
PUT_ERR("ERROR in yac_interpolation_add_sum_at_tgt/"
"yac_interpolation_add_weight_sum_mvp_at_tgt frac");
}
}
}
xt_redist_delete(src_redists[1]);
xt_redist_delete(src_redists[0]);
}
xt_finalize();
MPI_Finalize();
return TEST_EXIT_CODE;
}
static void init_tgt_data(
if (tgt_data == NULL) return;
for (size_t j = 0; j < tgt_field_size; ++j)
tgt_data[i][j] = UNSET;
}
static int check_tgt_data(
double ** tgt_data, double * ref_tgt_data_raw,
int diff_count = 0;
if (tgt_data != NULL)
++i, ref_tgt_data_raw += tgt_field_size)
for (size_t j = 0; j < tgt_field_size; ++j)
if (fabs(tgt_data[i][j] - ref_tgt_data_raw[j]) > 1e-9)
++diff_count;
return diff_count;
}
static int check_interpolation_execute(
double * ref_tgt_data_raw,
size_t tgt_field_size,
size_t collection_size) {
int diff_count = 0;
PUT_ERR("ERROR in yac_interpolation_execute_put/get_test");
diff_count +=
check_tgt_data(tgt_data, ref_tgt_data_raw, tgt_field_size,
collection_size);
MPI_Barrier(MPI_COMM_WORLD);
PUT_ERR("ERROR in yac_interpolation_execute_put/get_test");
if (src_data != NULL) {
if (tgt_data == NULL)
}
diff_count +=
check_tgt_data(tgt_data, ref_tgt_data_raw, tgt_field_size,
collection_size);
if (tgt_data != NULL)
if (src_data != NULL)
diff_count +=
check_tgt_data(tgt_data, ref_tgt_data_raw, tgt_field_size,
collection_size);
if (tgt_data != NULL)
if (src_data != NULL)
check_tgt_data(tgt_data, ref_tgt_data_raw, tgt_field_size,
collection_size);
return diff_count;
}
static int check_interpolation_execute_frac(
double *** src_frac_mask, double ** tgt_data, double * ref_tgt_data_raw,
int diff_count = 0;
PUT_ERR("ERROR in yac_interpolation_execute_put/get_test");
diff_count +=
check_tgt_data(tgt_data, ref_tgt_data_raw, tgt_field_size,
collection_size);
MPI_Barrier(MPI_COMM_WORLD);
PUT_ERR("ERROR in yac_interpolation_execute_put/get_test");
if (src_data != NULL) {
if (tgt_data == NULL)
}
diff_count +=
check_tgt_data(tgt_data, ref_tgt_data_raw, tgt_field_size,
collection_size);
return diff_count;
}
static int check_interpolation(
double * src_data_raw, size_t num_src_fields, size_t src_field_size,
double * ref_tgt_data_raw,
size_t tgt_field_size,
size_t collection_size) {
double *** src_data = NULL;
if (src_data_raw != NULL) {
double ** src_field_data =
++i, src_field_data += num_src_fields) {
src_data[i] = src_field_data;
for (size_t j = 0; j < num_src_fields;
++j, src_data_raw += src_field_size)
src_data[i][j] = src_data_raw;
}
}
double ** tgt_data = NULL;
if (ref_tgt_data_raw != NULL) {
double * tgt_data_raw =
tgt_data_raw[i] = UNSET;
for (
size_t i = 0; i <
collection_size; ++i, tgt_data_raw += tgt_field_size)
tgt_data[i] = tgt_data_raw;
}
int diff_count =
check_interpolation_execute(
interp, src_data, tgt_data, ref_tgt_data_raw,
diff_count +=
check_interpolation_execute(
interp, src_data, tgt_data, ref_tgt_data_raw,
if (ref_tgt_data_raw != NULL) {
free(tgt_data[0]);
free(tgt_data);
}
if (src_data_raw != NULL) {
free(src_data[0]);
free(src_data);
}
return diff_count;
}
static int check_interpolation_frac(
double * src_data_raw, double * src_frac_mask_raw,
size_t num_src_fields, size_t src_field_size,
double * ref_tgt_data_raw,
size_t tgt_field_size,
size_t collection_size) {
double *** src_data = NULL;
if (src_data_raw != NULL) {
double ** src_field_data =
++i, src_field_data += num_src_fields) {
src_data[i] = src_field_data;
for (size_t j = 0; j < num_src_fields;
++j, src_data_raw += src_field_size)
src_data[i][j] = src_data_raw;
}
}
double *** src_frac_data = NULL;
if (src_frac_mask_raw != NULL) {
double ** src_frac_mask_data =
++i, src_frac_mask_data += num_src_fields) {
src_frac_data[i] = src_frac_mask_data;
for (size_t j = 0; j < num_src_fields;
++j, src_frac_mask_raw += src_field_size)
src_frac_data[i][j] = src_frac_mask_raw;
}
}
if ((src_data != NULL) && (src_frac_data != NULL))
for (size_t j = 0; j < num_src_fields; ++j)
for (size_t k = 0; k < src_field_size; ++k)
src_data[i][j][k] *= src_frac_data[i][j][k];
double ** tgt_data = NULL;
if (ref_tgt_data_raw != NULL) {
double * tgt_data_raw =
tgt_data_raw[i] = UNSET;
for (
size_t i = 0; i <
collection_size; ++i, tgt_data_raw += tgt_field_size)
tgt_data[i] = tgt_data_raw;
}
int diff_count =
check_interpolation_execute_frac(
interp, src_data, src_frac_data, tgt_data, ref_tgt_data_raw,
diff_count +=
check_interpolation_execute_frac(
interp, src_data, src_frac_data, tgt_data, ref_tgt_data_raw,
if (ref_tgt_data_raw != NULL) {
free(tgt_data[0]);
free(tgt_data);
}
if (src_frac_mask_raw != NULL) {
free(src_frac_data[0]);
free(src_frac_data);
}
if (src_data_raw != NULL) {
free(src_data[0]);
free(src_data);
}
return diff_count;
}
int main(int argc, char **argv)
void yac_interpolation_add_sum_at_src(struct yac_interpolation *interp, Xt_redist *halo_redists, size_t tgt_count, size_t *num_src_per_tgt, size_t *src_field_idx, size_t *src_idx, size_t num_src_fields, Xt_redist result_redist)
struct yac_interpolation * yac_interpolation_copy(struct yac_interpolation *interp)
void yac_interpolation_add_weight_sum_mvp_at_tgt(struct yac_interpolation *interp, Xt_redist *src_redists, size_t *tgt_pos, size_t tgt_count, size_t *num_src_per_tgt, double *weights, size_t *src_field_idx, size_t *src_idx, size_t num_src_fields)
struct yac_interpolation * yac_interpolation_new(size_t collection_size, double frac_mask_fallback_value, double scale_factor, double scale_summand)
void yac_interpolation_execute(struct yac_interpolation *interp, double ***src_fields, double **tgt_field)
void yac_interpolation_add_direct_mf(struct yac_interpolation *interp, Xt_redist *redists, size_t num_src_fields)
void yac_interpolation_delete(struct yac_interpolation *interp)
void yac_interpolation_add_sum_at_tgt(struct yac_interpolation *interp, Xt_redist *src_redists, size_t *tgt_pos, size_t tgt_count, size_t *num_src_per_tgt, size_t *src_field_idx, size_t *src_idx, size_t num_src_fields)
int yac_interpolation_execute_put_test(struct yac_interpolation *interp)
void yac_interpolation_add_fixed(struct yac_interpolation *interp, double value, size_t count, size_t *pos)
void yac_interpolation_execute_wait(struct yac_interpolation *interp)
void yac_interpolation_execute_frac(struct yac_interpolation *interp, double ***src_fields, double ***src_frac_masks, double **tgt_field)
void yac_interpolation_execute_get(struct yac_interpolation *interp, double **tgt_field)
void yac_interpolation_execute_get_async(struct yac_interpolation *interp, double **tgt_field)
void yac_interpolation_add_direct(struct yac_interpolation *interp, Xt_redist redist)
void yac_interpolation_execute_put_frac(struct yac_interpolation *interp, double ***src_fields, double ***src_frac_masks)
void yac_interpolation_execute_put(struct yac_interpolation *interp, double ***src_fields)
void yac_interpolation_add_weight_sum_mvp_at_src(struct yac_interpolation *interp, Xt_redist *halo_redists, size_t tgt_count, size_t *num_src_per_tgt, double *weights, size_t *src_field_idx, size_t *src_idx, size_t num_src_fields, Xt_redist result_redist)
int yac_interpolation_execute_get_test(struct yac_interpolation *interp)
double const YAC_FRAC_MASK_NO_VALUE