30#define YAC_RAD (0.01745329251994329576923690768489)
32int main(
int argc,
char *argv[]) {
35 PUT_ERR(
"wrong number of arguments\n");
43 PUT_ERR(
"invalid argument (has to be either \"src\" or \"tgt\")\n");
49 xt_initialize(MPI_COMM_WORLD);
51 int comm_rank, comm_size;
52 MPI_Comm_rank(MPI_COMM_WORLD, &comm_rank);
53 MPI_Comm_size(MPI_COMM_WORLD, &comm_size);
54 MPI_Barrier(MPI_COMM_WORLD);
57 PUT_ERR(
"ERROR: wrong number of processes");
63 unsigned is_source = comm_rank > 0;
64 unsigned is_target = comm_rank < 2;
109 char const source_grid_name[] =
"source_grid";
110 char const target_grid_name[] =
"target_grid";
112 double src_coordinates_x[] = {0,1,2,3,4,5};
113 double src_coordinates_y[] = {0,1,2,3,4};
114 size_t src_global_num_cells[] = {5,4};
115 int src_with_halo = 1;
116 size_t src_local_start[2][2] = {{0,0}, {0,2}};
117 size_t src_local_count[2] = {5,2};
119 for (
size_t i = 0; i <= src_global_num_cells[0]; ++i)
120 src_coordinates_x[i] *=
YAC_RAD;
121 for (
size_t i = 0; i <= src_global_num_cells[1]; ++i)
122 src_coordinates_y[i] *=
YAC_RAD;
129 src_coordinates_x, src_coordinates_y, src_global_num_cells,
130 src_local_start[comm_rank-1], src_local_count, src_with_halo)):
133 double tgt_coordinates_x[] = {0.7,1.7,2.7,3.7,4.7,5.7,6.7};
134 double tgt_coordinates_y[] = {0.7,1.7,2.7,3.7};
135 size_t tgt_global_num_cells[] = {6,3};
136 int tgt_with_halo = 1;
137 size_t tgt_local_start[2][2] = {{3,0}, {0,0}};
138 size_t tgt_local_count[2] = {3,3};
140 for (
size_t i = 0; i <= tgt_global_num_cells[0]; ++i)
141 tgt_coordinates_x[i] *=
YAC_RAD;
142 for (
size_t i = 0; i <= tgt_global_num_cells[1]; ++i)
143 tgt_coordinates_y[i] *=
YAC_RAD;
150 tgt_coordinates_x, tgt_coordinates_y, tgt_global_num_cells,
151 tgt_local_start[comm_rank], tgt_local_count, tgt_with_halo)):
161 size_t num_src_fields =
sizeof(src_fields) /
sizeof(src_fields[0]);
163 {.location =
YAC_LOC_CORNER, .coordinates_idx = SIZE_MAX, .masks_idx = SIZE_MAX};
170 char const * sum_weight_file_name =
171 "test_interpolation_parallel5_weight_file.nc";
172 if (comm_rank == 0) {
173 int src_indices[4*5*4] =
174 {0,1,6,7, 1,2,7,8, 2,3,8,9, 3,4,9,10, 4,5,10,11,
175 6,7,12,13, 7,8,13,14, 8,9,14,15, 9,10,15,16, 10,11,16,17,
176 12,13,18,19, 13,14,19,20, 14,15,20,21, 15,16,21,22, 16,17,22,23,
177 18,19,24,25, 19,20,25,26, 20,21,26,27, 21,22,27,28, 22,23,28,29};
178 int tgt_indices[4*5*4] =
179 {0,0,0,0, 1,1,1,1, 2,2,2,2, 3,3,3,3, 4,4,4,4,
180 7,7,7,7, 8,8,8,8, 9,9,9,9, 10,10,10,10, 11,11,11,11,
181 14,14,14,14, 15,15,15,15, 16,16,16,16, 17,17,17,17, 18,18,18,18,
182 21,21,21,21, 22,22,22,22, 23,23,23,23, 24,24,24,24, 25,25,25,25};
183 double weights[4*5*4];
184 for (
size_t i = 0; i < 4*5*4; ++i) weights[i] = 1.0;
185 size_t num_links = 4*5*4;
189 int num_links_per_field[1] = {num_links};
190 int * tgt_id_fixed = NULL;
191 size_t num_fixed_tgt = 0;
192 double * fixed_values = NULL;
193 int * num_tgt_per_fixed_value = NULL;
194 size_t num_fixed_values = 0;
197 sum_weight_file_name, src_indices, tgt_indices, weights, num_links,
199 num_fixed_tgt, fixed_values, num_tgt_per_fixed_value,
200 num_fixed_values, tgt_location,
201 source_grid_name, target_grid_name);
203 MPI_Barrier(MPI_COMM_WORLD);
205 size_t const num_method_stacks = 3;
213 .max_search_distance =
223 for (
size_t method_stack_index = 0;
224 method_stack_index < num_method_stacks; ++method_stack_index) {
228 method_stacks[method_stack_index], interp_grid);
234 weights, reorder_type, 1,
238 "test_interpolation_parallel5_weight_file_out.nc";
249 double * source_data_field =
252 double * source_data_pointset[1] = {source_data_field};
253 double ** source_data[1] = {source_data_pointset};
257 double * target_data_field =
260 double * target_data[1] = {target_data_field};
267 for (
size_t i = 0; i < src_grid_data->
num_vertices; ++i)
268 source_data_field[i] =
270 ((double)(src_grid_data->
vertex_ids[i])):(-1.0);
273 for (
int use_put_get = 0; use_put_get < 2; ++use_put_get) {
279 for (
size_t i = 0; i < tgt_grid_data->
num_vertices; ++i)
280 target_data_field[i] = -1.0;
298 double ref_target_data[3][4*7] =
299 {{3.5,4.5,5.5,6.5,7.5,1337,1337,
300 9.5,10.5,11.5,12.5,13.5,1337,1337,
301 15.5,16.5,17.5,18.5,19.5,1337,1337,
302 21.5,22.5,23.5,24.5,25.5,1337,1337},
304 13,14,15,16,17,17,17,
305 19,20,21,22,23,23,23,
306 25,26,27,28,29,29,29},
307 {14,18,22,26,30,1337,1337,
308 38,42,46,50,54,1337,1337,
309 62,66,70,74,78,1337,1337,
310 86,90,94,98,102,1337,1337}};
312 for (
size_t i = 0; i < tgt_grid_data->
num_vertices; ++i) {
317 [method_stack_index][tgt_grid_data->
vertex_ids[i]]))
318 PUT_ERR(
"error in interpolated data on target side\n");
320 if (target_data[0][i] != -1.0)
321 PUT_ERR(
"error in interpolated data on target side\n");
327 if (is_target) free(target_data_field);
328 if (is_source) free(source_data_field);
340 if (comm_rank == 0) unlink(sum_weight_file_name);
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)
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)
struct yac_basic_grid_data yac_generate_basic_grid_data_reg2d(double const *global_coords_x, double const *global_coords_y, size_t const num_global_cells_[2], size_t const local_start[2], size_t const local_count[2], int with_halo)
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, enum yac_interp_file_on_missing_file on_missing_file, enum yac_interp_file_on_success on_success)
#define YAC_INTERP_FILE_ON_SUCCESS_DEFAULT
#define YAC_INTERP_FILE_ON_MISSING_FILE_DEFAULT
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
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, int is_source, int is_target)
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, enum yac_weight_file_on_existing on_existing)
yac_interp_weights_reorder_type
@ YAC_MAPPING_ON_TGT
weights will be applied at target processes
@ YAC_MAPPING_ON_SRC
weights will be applied at source processes
@ YAC_WEIGHT_FILE_ERROR
error when weight file existis already
void yac_interpolation_execute(struct yac_interpolation *interp, double ***src_fields, double **tgt_field)
Execute interpolation synchronously and write results to the target field.
void yac_interpolation_delete(struct yac_interpolation *interp)
Free an interpolation object and release all resources.
void yac_interpolation_execute_get(struct yac_interpolation *interp, double **tgt_field)
Complete interpolation and write results to the target field (get phase).
void yac_interpolation_execute_put(struct yac_interpolation *interp, double ***src_fields)
Provide source field data and start asynchronous execution of interpolation (put phase).
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[]
int double_are_unequal(double a, double b)
char const weight_file_out[]
void write_weight_file(char const *file_name, int const *src_id, int const *tgt_id, double const *weights, unsigned num_links, enum yac_location const *src_locations, unsigned num_src_fields, int const *num_links_per_src_field, int *tgt_id_fixed, unsigned num_fixed_tgt, double *fixed_values, int *num_tgt_per_fixed_value, unsigned num_fixed_values, enum yac_location tgt_location, char const *src_grid_name, char const *tgt_grid_name)