A test for internal interpolation routines.
#include <stdlib.h>
#include <stdio.h>
#include <math.h>
#include <unistd.h>
#include <string.h>
#include <mpi.h>
#include "tests.h"
#include "test_common.h"
#include "weight_file_common.h"
#include "dist_grid_utils.h"
#define YAC_RAD (0.01745329251994329576923690768489)
int main(
int argc,
char *argv[]) {
if (argc != 2) {
PUT_ERR("wrong number of arguments\n");
return TEST_EXIT_CODE;
}
PUT_ERR("invalid argument (has to be either \"src\" or \"tgt\")\n");
return TEST_EXIT_CODE;
}
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);
if (comm_size != 3) {
PUT_ERR("ERROR: wrong number of processes");
xt_finalize();
MPI_Finalize();
return TEST_EXIT_CODE;
}
unsigned is_source = comm_rank > 0;
unsigned is_target = comm_rank < 2;
char const source_grid_name[] = "source_grid";
char const target_grid_name[] = "target_grid";
double src_coordinates_x[] = {0,1,2,3,4,5};
double src_coordinates_y[] = {0,1,2,3,4};
size_t src_global_num_cells[] = {5,4};
int src_with_halo = 1;
size_t src_local_start[2][2] = {{0,0}, {0,2}};
size_t src_local_count[2] = {5,2};
for (size_t i = 0; i <= src_global_num_cells[0]; ++i)
for (size_t i = 0; i <= src_global_num_cells[1]; ++i)
(is_source)?
source_grid_name,
yac_generate_basic_grid_data_reg2d(
src_coordinates_x, src_coordinates_y, src_global_num_cells,
src_local_start[comm_rank-1], src_local_count, src_with_halo)):
double tgt_coordinates_x[] = {0.7,1.7,2.7,3.7,4.7,5.7,6.7};
double tgt_coordinates_y[] = {0.7,1.7,2.7,3.7};
size_t tgt_global_num_cells[] = {6,3};
int tgt_with_halo = 1;
size_t tgt_local_start[2][2] = {{3,0}, {0,0}};
size_t tgt_local_count[2] = {3,3};
for (size_t i = 0; i <= tgt_global_num_cells[0]; ++i)
for (size_t i = 0; i <= tgt_global_num_cells[1]; ++i)
(is_target)?
target_grid_name,
yac_generate_basic_grid_data_reg2d(
tgt_coordinates_x, tgt_coordinates_y, tgt_global_num_cells,
tgt_local_start[comm_rank], tgt_local_count, tgt_with_halo)):
{
size_t num_src_fields = sizeof(src_fields) / sizeof(src_fields[0]);
char const * sum_weight_file_name =
"test_interpolation_parallel5_weight_file.nc";
if (comm_rank == 0) {
int src_indices[4*5*4] =
{0,1,6,7, 1,2,7,8, 2,3,8,9, 3,4,9,10, 4,5,10,11,
6,7,12,13, 7,8,13,14, 8,9,14,15, 9,10,15,16, 10,11,16,17,
12,13,18,19, 13,14,19,20, 14,15,20,21, 15,16,21,22, 16,17,22,23,
18,19,24,25, 19,20,25,26, 20,21,26,27, 21,22,27,28, 22,23,28,29};
int tgt_indices[4*5*4] =
{0,0,0,0, 1,1,1,1, 2,2,2,2, 3,3,3,3, 4,4,4,4,
7,7,7,7, 8,8,8,8, 9,9,9,9, 10,10,10,10, 11,11,11,11,
14,14,14,14, 15,15,15,15, 16,16,16,16, 17,17,17,17, 18,18,18,18,
21,21,21,21, 22,22,22,22, 23,23,23,23, 24,24,24,24, 25,25,25,25};
double weights[4*5*4];
for (size_t i = 0; i < 4*5*4; ++i) weights[i] = 1.0;
size_t num_links = 4*5*4;
int num_links_per_field[1] = {num_links};
int * tgt_id_fixed = NULL;
size_t num_fixed_tgt = 0;
double * fixed_values = NULL;
int * num_tgt_per_fixed_value = NULL;
size_t num_fixed_values = 0;
write_weight_file(
sum_weight_file_name, src_indices, tgt_indices, weights, num_links,
num_fixed_tgt, fixed_values, num_tgt_per_fixed_value,
num_fixed_values, tgt_location,
source_grid_name, target_grid_name);
}
MPI_Barrier(MPI_COMM_WORLD);
size_t const num_method_stacks = 3;
NULL},
.max_search_distance =
NULL},
NULL}};
for (size_t method_stack_index = 0;
method_stack_index < num_method_stacks; ++method_stack_index) {
method_stacks[method_stack_index], interp_grid);
weights, reorder_type, 1,
char const weight_file_out[] =
"test_interpolation_parallel5_weight_file_out.nc";
weights, weight_file_out, source_grid_name, target_grid_name, 0, 0);
double * source_data_field =
is_source?
double * source_data_pointset[1] = {source_data_field};
double ** source_data[1] = {source_data_pointset};
double * target_data_field =
is_target?
double * target_data[1] = {target_data_field};
if (is_source) {
source_data_field[i] =
}
for (int use_put_get = 0; use_put_get < 2; ++use_put_get) {
if (is_target) {
target_data_field[i] = -1.0;
}
if (use_put_get) {
if (is_source)
if (is_target)
} else {
}
if (is_target) {
double ref_target_data[3][4*7] =
{{3.5,4.5,5.5,6.5,7.5,1337,1337,
9.5,10.5,11.5,12.5,13.5,1337,1337,
15.5,16.5,17.5,18.5,19.5,1337,1337,
21.5,22.5,23.5,24.5,25.5,1337,1337},
{7,8,9,10,11,11,11,
13,14,15,16,17,17,17,
19,20,21,22,23,23,23,
25,26,27,28,29,29,29},
{14,18,22,26,30,1337,1337,
38,42,46,50,54,1337,1337,
62,66,70,74,78,1337,1337,
86,90,94,98,102,1337,1337}};
if (double_are_unequal(
target_data[0][i],
ref_target_data
[method_stack_index][tgt_grid_data->
vertex_ids[i]]))
PUT_ERR("error in interpolated data on target side\n");
} else {
if (target_data[0][i] != -1.0)
PUT_ERR("error in interpolated data on target side\n");
}
}
}
}
if (is_target) free(target_data_field);
if (is_source) free(source_data_field);
if (comm_rank == 0) unlink(weight_file_out);
}
if (comm_rank == 0) unlink(sum_weight_file_name);
}
xt_finalize();
MPI_Finalize();
return TEST_EXIT_CODE;
}
struct yac_basic_grid * yac_basic_grid_new(char const *name, struct yac_basic_grid_data grid_data)
struct yac_basic_grid_data * yac_basic_grid_get_data(struct yac_basic_grid *grid)
struct yac_basic_grid * yac_basic_grid_empty_new(char const *name)
void yac_basic_grid_delete(struct yac_basic_grid *grid)
int main(int argc, char **argv)
void yac_dist_grid_pair_delete(struct yac_dist_grid_pair *grid_pair)
struct yac_dist_grid_pair * yac_dist_grid_pair_new(struct yac_basic_grid *grid_a, struct yac_basic_grid *grid_b, MPI_Comm comm)
void yac_interp_grid_delete(struct yac_interp_grid *interp_grid)
struct yac_interp_grid * yac_interp_grid_new(struct yac_dist_grid_pair *grid_pair, char const *src_grid_name, char const *tgt_grid_name, size_t num_src_fields, struct yac_interp_field const *src_fields, struct yac_interp_field const tgt_field)
void yac_interp_method_delete(struct interp_method **method)
struct yac_interp_weights * yac_interp_method_do_search(struct interp_method **method, struct yac_interp_grid *interp_grid)
struct interp_method * yac_interp_method_avg_new(enum yac_interp_avg_weight_type weight_type, int partial_coverage)
@ YAC_INTERP_AVG_ARITHMETIC
struct interp_method * yac_interp_method_file_new(char const *weight_file_name)
struct interp_method * yac_interp_method_fixed_new(double value)
struct interp_method * yac_interp_method_nnn_new(struct yac_nnn_config config)
@ YAC_INTERP_NNN_AVG
average of n source points
#define YAC_INTERP_NNN_MAX_SEARCH_DISTANCE_DEFAULT
void yac_interp_weights_delete(struct yac_interp_weights *weights)
void yac_interp_weights_write_to_file(struct yac_interp_weights *weights, char const *filename, char const *src_grid_name, char const *tgt_grid_name, size_t src_grid_size, size_t tgt_grid_size)
struct yac_interpolation * yac_interp_weights_get_interpolation(struct yac_interp_weights *weights, enum yac_interp_weights_reorder_type reorder, size_t collection_size, double frac_mask_fallback_value, double scaling_factor, double scaling_summand, char const *yaxt_exchanger_name)
yac_interp_weights_reorder_type
@ YAC_MAPPING_ON_TGT
weights will be applied at target processes
@ YAC_MAPPING_ON_SRC
weights will be appied at source processes
void yac_interpolation_execute(struct yac_interpolation *interp, double ***src_fields, double **tgt_field)
void yac_interpolation_delete(struct yac_interpolation *interp)
void yac_interpolation_execute_get(struct yac_interpolation *interp, double **tgt_field)
void yac_interpolation_execute_put(struct yac_interpolation *interp, double ***src_fields)
double const YAC_FRAC_MASK_NO_VALUE
enum yac_location location
struct yac_interp_field tgt_field
struct yac_dist_grid_pair * grid_pair
struct yac_interp_field src_fields[]