34 double const tgt_coords[3],
int src_cell_id,
size_t src_cell_idx,
35 int const ** global_results_points,
double ** result_weights,
38 double const tgt_coords[3],
int src_cell_id,
size_t src_cell_idx,
39 int const ** global_results_points,
double ** result_weights,
42 double const tgt_coords[3],
int src_cell_id,
size_t src_cell_idx,
43 int const ** global_results_points,
double ** result_weights,
46 double const tgt_coords[3],
int src_cell_id,
size_t src_cell_idx,
47 int const ** global_results_points,
double ** result_weights,
49static char const *
grid_names[2] = {
"src_grid",
"tgt_grid"};
56 int comm_rank, comm_size;
57 MPI_Comm_rank(MPI_COMM_WORLD, &comm_rank);
58 MPI_Comm_size(MPI_COMM_WORLD, &comm_size);
59 MPI_Barrier(MPI_COMM_WORLD);
62 PUT_ERR(
"ERROR: wrong number of processes");
78 enum {NUM_CALLBACKS =
sizeof(callbacks) /
sizeof(callbacks[0])};
80 for (
size_t i = 0;
i < NUM_CALLBACKS; ++
i)
82 callbacks[i].func, (
void*)&(callbacks[
i].user_data), callbacks[i].
key);
85 callbacks[1].func, (
void*)&(callbacks[1].
user_data), callbacks[1].
key);
87 for (
size_t i = 0;
i < NUM_CALLBACKS; ++
i) {
95 if (func != callbacks[i].func)
97 "ERROR in yac_interp_method_callback_get_compute_weights_callback");
100 "ERROR in yac_interp_method_callback_get_compute_weights_callback");
120 int is_tgt = comm_rank == 2;
121 double coordinates_x[2][5] = {{0.0,1.0,2.0,3.0}, {-0.5,0.5,1.5,2.5,3.5}};
122 double coordinates_y[2][4] = {{0.0,1.0,2.0}, {-0.5,0.5,1.5,2.5}};
123 size_t const num_cells[2][2] = {{3,2}, {4,3}};
124 size_t local_start[2][2][2] = {{{0,0},{1,0}}, {{0,0}}};
125 size_t local_count[2][2][2] = {{{1,2},{2,2}}, {{4,3}}};
135 local_start[is_tgt][(is_tgt)?0:comm_rank],
136 local_count[is_tgt][(is_tgt)?0:comm_rank],
with_halo);
160 {.location =
YAC_LOC_CELL, .coordinates_idx = SIZE_MAX, .masks_idx = SIZE_MAX}};
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};
182 {{-1.0, -1.0, -1.0, -1.0, -1.0,
183 -1.0, -1.0, -1.0, -1.0, -1.0,
184 -1.0, -1.0, -1.0, -1.0, -1.0,
185 -1.0, -1.0, -1.0, -1.0, -1.0},
186 {-1.0, -1.0, -1.0, -1.0, -1.0,
187 -1.0, 0.5*0+0.125*(0+1+4+5), 0.5*1+0.125*(1+2+5+6), 0.5*2+0.125*(2+3+6+7), -1.0,
188 -1.0, 0.5*3+0.125*(4+5+8+9), 0.5*4+0.125*(5+6+9+10), 0.5*5+0.125*(6+7+10+11), -1.0,
189 -1.0, -1.0, -1.0, -1.0, -1.0},
190 {-1.0, -1.0, -1.0, -1.0, -1.0,
191 -1.0, -1.0, 0.5*1+0.125*(1+2+5+6), -1.0, -1.0,
192 -1.0, 0.5*3+0.125*(4+5+8+9), -1.0, 0.5*5+0.125*(6+7+10+11), -1.0,
193 -1.0, -1.0, -1.0, -1.0, -1.0}};
194 double collection_factor[] = {0.0, 10.0, 10.0};
212 for (
size_t reorder_idx = 0; reorder_idx < 2; ++reorder_idx) {
221 double *** src_data = NULL;
222 double ** tgt_data = NULL;
228 tgt_data[collection_idx] =
230 for (
size_t i = 0;
i <
grid_data.num_vertices; ++
i)
231 tgt_data[collection_idx][i] = -1337.0;
237 src_data[collection_idx] =
xmalloc(2 *
sizeof(**src_data));
238 src_data[collection_idx][0] =
240 for (
size_t i = 0;
i <
grid_data.num_vertices; ++
i)
241 src_data[collection_idx][0][i] =
242 (
double)(vertex_ids[
i] + 10 * collection_idx);
243 src_data[collection_idx][1] =
246 src_data[collection_idx][1][i] =
247 (
double)(cell_ids[
i] + 10 * collection_idx);
258 for (
size_t i = 0;
i <
grid_data.num_vertices; ++
i)
259 if (ref_tgt_data[test_idx][vertex_ids[i]] == -1.0) {
260 if (tgt_data[collection_idx][i] != -1.0)
263 if (fabs(tgt_data[collection_idx][i] -
264 (ref_tgt_data[test_idx][vertex_ids[i]] +
265 (
double)collection_idx *
266 collection_factor[test_idx])) > 1e-9)
272 ++collection_idx) free(tgt_data[collection_idx]);
277 for (
size_t src_field_idx = 0; src_field_idx < num_src_fields;
279 free(src_data[collection_idx][src_field_idx]);
280 free(src_data[collection_idx]);
307 double src_coordinates_x[51];
308 double src_coordinates_y[51];
309 size_t const src_num_cells[2] = {50,50};
310 size_t src_local_start[3][2] = {{20,10},{20,0},{0,0}};
311 size_t src_local_count[3][2] = {{30,40},{30,10},{20,50}};
313 for (
size_t i = 0;
i <= src_num_cells[0]; ++
i)
314 src_coordinates_x[i] = (
double)
i *
YAC_RAD;
315 for (
size_t i = 0;
i <= src_num_cells[1]; ++
i)
316 src_coordinates_y[i] = (
double)
i *
YAC_RAD;
320 src_coordinates_x, src_coordinates_y, src_num_cells,
321 src_local_start[comm_rank], src_local_count[comm_rank],
with_halo);
323 double tgt_coordinates_x[50];
324 double tgt_coordinates_y[50];
325 size_t const tgt_num_cells[2] = {49,49};
326 size_t tgt_local_start[2] = {0,0};
327 size_t tgt_local_count[2] = {49,49};
328 for (
size_t i = 0;
i <= tgt_num_cells[0]; ++
i)
329 tgt_coordinates_x[i] = (0.5 + (
double)i) *
YAC_RAD;
330 for (
size_t i = 0;
i <= tgt_num_cells[1]; ++
i)
331 tgt_coordinates_y[i] = (0.5 + (
double)i) *
YAC_RAD;
337 tgt_coordinates_x, tgt_coordinates_y, tgt_num_cells,
338 tgt_local_start, tgt_local_count,
with_halo);
353 size_t num_src_fields =
sizeof(src_fields) /
sizeof(src_fields[0]);
355 {.location =
YAC_LOC_CORNER, .coordinates_idx = SIZE_MAX, .masks_idx = SIZE_MAX};
375 double *** src_data = NULL;
376 double ** tgt_data = NULL;
378 if (comm_rank == 0) {
379 tgt_data =
xmalloc(1 *
sizeof(*tgt_data));
380 tgt_data[0] =
xmalloc(tgt_grid_data.num_vertices *
sizeof(**tgt_data));
381 for (
size_t i = 0;
i < tgt_grid_data.num_vertices; ++
i)
382 tgt_data[0][i] = -1.0;
385 src_data =
xmalloc(1 *
sizeof(*src_data));
386 src_data[0] =
xmalloc(1 *
sizeof(**src_data));
389 for (
size_t i = 0;
i < src_grid_data.
num_cells; ++
i)
390 src_data[0][0][i] = (
double)(src_grid_data.
cell_ids[
i]);
394 if (comm_rank == 0) {
397 for (
size_t i = 0;
i < tgt_grid_data.num_vertices; ++
i)
398 if (tgt_data[0][i] != (
double)(tgt_grid_data.vertex_ids[i]))
405 free(src_data[0][0]);
426 double const tgt_coords[3],
int src_cell_id,
size_t src_cell_idx,
427 int const ** global_results_points,
double ** result_weights,
428 size_t * result_count,
void *
user_data) {
434 size_t num_src_fields = *(
size_t*)
user_data;
436 for (
size_t i = 0; i < num_src_fields; ++i) {
437 global_results_points[i] = NULL;
438 result_weights[i] = NULL;
444 double const tgt_coords[3],
int src_cell_id,
size_t src_cell_idx,
445 int const ** global_results_points,
double ** result_weights,
446 size_t * result_count,
void *
user_data) {
453 if ((info->num_grid_cells <= src_cell_idx) ||
454 (info->cell_ids[src_cell_idx] != src_cell_id)) {
455 fputs(
"ERROR(test2_callback): inconsistent data\n", stderr);
459 static double vertex_weights[4] = {0.125, 0.125, 0.125, 0.125};
460 static double cell_weight[1] = {0.5};
462 static int vertex_result_points[4];
463 static int cell_result_points[1];
465 for (
size_t i = 0; i < 4; ++i)
466 vertex_result_points[i] =
467 (
int)(info->vertex_ids[info->cell_to_vertex[4*src_cell_idx + i]]);
468 cell_result_points[0] = (int)src_cell_id;
470 global_results_points[0] = vertex_result_points;
471 global_results_points[1] = cell_result_points;
472 result_weights[0] = vertex_weights;
473 result_weights[1] = cell_weight;
479 double const tgt_coords[3],
int src_cell_id,
size_t src_cell_idx,
480 int const ** global_results_points,
double ** result_weights,
481 size_t * result_count,
void *
user_data) {
488 static int result_points[1];
489 static double weight[1] = {1.0};
491 result_points[0] = (int)src_cell_id;
493 global_results_points[0] = result_points;
494 result_weights[0] = weight;
499 double const tgt_coords[3],
int src_cell_id,
size_t src_cell_idx,
500 int const ** global_results_points,
double ** result_weights,
501 size_t * result_count,
void *
user_data) {
508 if ((info->num_grid_cells <= src_cell_idx) ||
509 (info->cell_ids[src_cell_idx] != src_cell_id)) {
510 fputs(
"ERROR(test2_callback): inconsistent data\n", stderr);
514 static double vertex_weights[4] = {0.125, 0.125, 0.125, 0.125};
515 static double cell_weight[1] = {0.5};
517 static int vertex_result_points[4];
518 static int cell_result_points[1];
520 for (
size_t i = 0; i < 4; ++i)
521 vertex_result_points[i] =
522 (
int)(info->vertex_ids[info->cell_to_vertex[4*src_cell_idx + i]]);
523 cell_result_points[0] = (int)src_cell_id;
526 global_results_points[0] = vertex_result_points;
527 global_results_points[1] = cell_result_points;
528 result_weights[0] = vertex_weights;
529 result_weights[1] = cell_weight;
533 global_results_points[0] = NULL;
534 global_results_points[1] = NULL;
535 result_weights[0] = NULL;
536 result_weights[1] = NULL;
struct yac_basic_grid * yac_basic_grid_new(char const *name, struct yac_basic_grid_data grid_data)
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)
void yac_interp_method_callback_get_compute_weights_callback(char const *key, yac_func_compute_weights *compute_weights_callback, void **user_data)
struct interp_method * yac_interp_method_callback_new(yac_func_compute_weights compute_weights_callback, void *user_data)
void yac_interp_method_callback_add_compute_weights_callback(yac_func_compute_weights compute_weights_callback, void *user_data, char const *key)
void yac_interp_method_callback_buf_free()
void(* yac_func_compute_weights)(double const tgt_coords[3], int src_cell_id, size_t src_cell_idx, int const **global_results_points, double **result_weights, size_t *result_count, void *user_data)
struct interp_method * yac_interp_method_fixed_new(double value)
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)
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
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.
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[]
static char const * grid_names[2]
static void test2_callback(double const tgt_coords[3], int src_cell_id, size_t src_cell_idx, int const **global_results_points, double **result_weights, size_t *result_count, void *user_data)
static void test4_callback(double const tgt_coords[3], int src_cell_id, size_t src_cell_idx, int const **global_results_points, double **result_weights, size_t *result_count, void *user_data)
static void test3_callback(double const tgt_coords[3], int src_cell_id, size_t src_cell_idx, int const **global_results_points, double **result_weights, size_t *result_count, void *user_data)
static void test1_callback(double const tgt_coords[3], int src_cell_id, size_t src_cell_idx, int const **global_results_points, double **result_weights, size_t *result_count, void *user_data)
struct @1 test_functions[]
struct yac_basic_grid ** grids
void yac_yaxt_init(MPI_Comm comm)