YAC 3.12.0
Yet Another Coupler
Loading...
Searching...
No Matches
test_interp_method_callback_parallel.c
Go to the documentation of this file.
1// Copyright (c) 2024 The YAC Authors
2//
3// SPDX-License-Identifier: BSD-3-Clause
4
5#include <stdlib.h>
6#include <string.h>
7
8#include "tests.h"
9#include "test_common.h"
10#include "geometry.h"
14#include "dist_grid_utils.h"
15#include "yac_mpi.h"
16
17#include <mpi.h>
18#include <yaxt.h>
19#include <netcdf.h>
20
26struct grid_info {
31};
32
33static void test1_callback(
34 double const tgt_coords[3], int src_cell_id, size_t src_cell_idx,
35 int const ** global_results_points, double ** result_weights,
36 size_t * result_count, void * user_data);
37static void test2_callback(
38 double const tgt_coords[3], int src_cell_id, size_t src_cell_idx,
39 int const ** global_results_points, double ** result_weights,
40 size_t * result_count, void * user_data);
41static void test3_callback(
42 double const tgt_coords[3], int src_cell_id, size_t src_cell_idx,
43 int const ** global_results_points, double ** result_weights,
44 size_t * result_count, void * user_data);
45static void test4_callback(
46 double const tgt_coords[3], int src_cell_id, size_t src_cell_idx,
47 int const ** global_results_points, double ** result_weights,
48 size_t * result_count, void * user_data);
49static char const * grid_names[2] = {"src_grid", "tgt_grid"};
50
51int main(void) {
52
54 yac_yaxt_init(MPI_COMM_WORLD);
55
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);
60
61 if (comm_size != 3) {
62 PUT_ERR("ERROR: wrong number of processes");
63 xt_finalize();
64 MPI_Finalize();
65 return TEST_EXIT_CODE;
66 }
67
68 { // test setting and getting callbacks
69
70 struct {
72 char const * key;
73 int user_data;
74 } callbacks[] =
75 {{.func = test1_callback, .key = "test1", .user_data = 1},
76 {.func = test2_callback, .key = "test2", .user_data = 2},
77 {.func = test3_callback, .key = "test3", .user_data = 3}};
78 enum {NUM_CALLBACKS = sizeof(callbacks) / sizeof(callbacks[0])};
79
80 for (size_t i = 0; i < NUM_CALLBACKS; ++i)
82 callbacks[i].func, (void*)&(callbacks[i].user_data), callbacks[i].key);
83 // add the second key again
85 callbacks[1].func, (void*)&(callbacks[1].user_data), callbacks[1].key);
86
87 for (size_t i = 0; i < NUM_CALLBACKS; ++i) {
88
90 void * user_data;
91
93 callbacks[i].key, &func, &user_data);
94
95 if (func != callbacks[i].func)
96 PUT_ERR(
97 "ERROR in yac_interp_method_callback_get_compute_weights_callback");
98 if (user_data != (void*)&(callbacks[i].user_data))
99 PUT_ERR(
100 "ERROR in yac_interp_method_callback_get_compute_weights_callback");
101 }
102 }
103
104 { // test with two source processes, one target process and two source fields
105 // the global source grid is a 4x2 grid:
106 // 08--14--09--15--10--16--11
107 // | | | |
108 // 08 03 10 04 12 05 13
109 // | | | |
110 // 04--07--05--09--06--11--07
111 // | | | |
112 // 01 00 03 01 05 02 06
113 // | | | |
114 // 00--00--01--02--02--04--03
115 //
116 //---------------
117 // setup
118 //---------------
119
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}}};
126 int with_halo = 0;
127 for (size_t i = 0; i <= num_cells[is_tgt][0]; ++i)
128 coordinates_x[is_tgt][i] *= YAC_RAD;
129 for (size_t i = 0; i <= num_cells[is_tgt][1]; ++i)
130 coordinates_y[is_tgt][i] *= YAC_RAD;
131
134 coordinates_x[is_tgt], coordinates_y[is_tgt], num_cells[is_tgt],
135 local_start[is_tgt][(is_tgt)?0:comm_rank],
136 local_count[is_tgt][(is_tgt)?0:comm_rank], with_halo);
137
138 yac_int * cell_ids =
139 xmalloc(grid_data.num_cells * sizeof(*cell_ids));
140 memcpy(cell_ids, grid_data.cell_ids,
141 grid_data.num_cells * sizeof(*cell_ids));
142 size_t * cell_to_vertex =
143 xmalloc(4 * grid_data.num_cells * sizeof(*cell_to_vertex));
145 4 * grid_data.num_cells * sizeof(*cell_to_vertex));
147 xmalloc(grid_data.num_vertices * sizeof(*vertex_ids));
148 memcpy(vertex_ids, grid_data.vertex_ids,
149 grid_data.num_vertices * sizeof(*vertex_ids));
150
151 struct yac_basic_grid * grids[2] =
154
155 struct yac_dist_grid_pair * grid_pair =
156 yac_dist_grid_pair_new(grids[0], grids[1], MPI_COMM_WORLD);
157
158 struct yac_interp_field src_fields[] =
159 {{.location = YAC_LOC_CORNER, .coordinates_idx = SIZE_MAX, .masks_idx = SIZE_MAX},
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};
164
165 struct yac_interp_grid * interp_grid =
168
169 struct grid_info grid_info = {
171 .vertex_ids = vertex_ids,
172 .cell_ids = cell_ids,
173 .num_grid_cells = grid_data.num_cells};
177 enum {
178 NUM_TEST_FUNCTIONS = sizeof(test_functions) / sizeof(test_functions[0])};
180 {(void*)&num_src_fields, (void*)&grid_info, (void*)&grid_info};
181 double ref_tgt_data[NUM_TEST_FUNCTIONS][20] =
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};
195
196 for (size_t test_idx = 0; test_idx < NUM_TEST_FUNCTIONS; ++test_idx) {
197
198 struct interp_method * method_stack[] =
200 test_functions[test_idx], user_data[test_idx]),
201 yac_interp_method_fixed_new(-1.0), NULL};
202
203 struct yac_interp_weights * weights =
204 yac_interp_method_do_search(method_stack, interp_grid);
205
206 enum yac_interp_weights_reorder_type reorder_type[2] =
208
209 for (size_t collection_size = 1; collection_size <= 16;
210 collection_size *= 2) {
211
212 for (size_t reorder_idx = 0; reorder_idx < 2; ++reorder_idx) {
213
214 struct yac_interpolation * interpolation =
216 weights, reorder_type[reorder_idx], collection_size,
217 YAC_FRAC_MASK_NO_VALUE, 1.0, 0.0, NULL, 1, 1);
218
219 // check generated interpolation
220 {
221 double *** src_data = NULL;
222 double ** tgt_data = NULL;
223
224 if (is_tgt) {
225 tgt_data = xmalloc(collection_size * sizeof(*tgt_data));
226 for (size_t collection_idx = 0; collection_idx < collection_size;
227 ++collection_idx) {
228 tgt_data[collection_idx] =
229 xmalloc(grid_data.num_vertices * sizeof(**tgt_data));
230 for (size_t i = 0; i < grid_data.num_vertices; ++i)
231 tgt_data[collection_idx][i] = -1337.0;
232 }
233 } else {
234 src_data = xmalloc(collection_size * sizeof(*src_data));
235 for (size_t collection_idx = 0; collection_idx < collection_size;
236 ++collection_idx) {
237 src_data[collection_idx] = xmalloc(2 * sizeof(**src_data));
238 src_data[collection_idx][0] =
239 xmalloc(grid_data.num_vertices * sizeof(***src_data));
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] =
244 xmalloc(grid_data.num_cells * sizeof(***src_data));
245 for (size_t i = 0; i < grid_data.num_cells; ++i)
246 src_data[collection_idx][1][i] =
247 (double)(cell_ids[i] + 10 * collection_idx);
248 }
249 }
250
251 yac_interpolation_execute(interpolation, src_data, tgt_data);
252
253 if (is_tgt) {
254
255 // check results
256 for (size_t collection_idx = 0; collection_idx < collection_size;
257 ++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)
261 PUT_ERR("wrong results");
262 } else {
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)
267 PUT_ERR("wrong results");
268 }
269 }
270
271 for (size_t collection_idx = 0; collection_idx < collection_size;
272 ++collection_idx) free(tgt_data[collection_idx]);
273 free(tgt_data);
274 } else {
275 for (size_t collection_idx = 0; collection_idx < collection_size;
276 ++collection_idx) {
277 for (size_t src_field_idx = 0; src_field_idx < num_src_fields;
278 ++src_field_idx)
279 free(src_data[collection_idx][src_field_idx]);
280 free(src_data[collection_idx]);
281 }
282 free(src_data);
283 }
284 }
285 yac_interpolation_delete(interpolation);
286 }
287 }
288
290 yac_interp_method_delete(method_stack);
291 }
292 yac_interp_grid_delete(interp_grid);
293 yac_dist_grid_pair_delete(grid_pair);
296 free(cell_to_vertex);
297 free(vertex_ids);
298 free(cell_ids);
299 }
300
301 {
302 // All three ranks are source processes. Rank 0 is also a target process.
303 //---------------
304 // setup
305 //---------------
306
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}};
312 int with_halo = 0;
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;
317
318 struct yac_basic_grid_data src_grid_data =
320 src_coordinates_x, src_coordinates_y, src_num_cells,
321 src_local_start[comm_rank], src_local_count[comm_rank], with_halo);
322
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;
332
333 struct yac_basic_grid_data tgt_grid_data;
334 if (comm_rank == 0)
335 tgt_grid_data =
337 tgt_coordinates_x, tgt_coordinates_y, tgt_num_cells,
338 tgt_local_start, tgt_local_count, with_halo);
339
340 struct yac_basic_grid * src_grid =
341 yac_basic_grid_new(grid_names[0], src_grid_data);
342 struct yac_basic_grid * tgt_grid = NULL;
343 if (comm_rank == 0)
344 tgt_grid = yac_basic_grid_new(grid_names[1], tgt_grid_data);
345 else
347
348 struct yac_dist_grid_pair * grid_pair =
349 yac_dist_grid_pair_new(src_grid, tgt_grid, MPI_COMM_WORLD);
350
351 struct yac_interp_field src_fields[] =
352 {{.location = YAC_LOC_CELL, .coordinates_idx = SIZE_MAX, .masks_idx = SIZE_MAX}};
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};
356
357 struct yac_interp_grid * interp_grid =
360
361 struct interp_method * method_stack[] =
363 yac_interp_method_fixed_new(-1.0), NULL};
364
365 struct yac_interp_weights * weights =
366 yac_interp_method_do_search(method_stack, interp_grid);
367
368 struct yac_interpolation * interpolation =
370 weights, YAC_MAPPING_ON_SRC, 1,
371 YAC_FRAC_MASK_NO_VALUE, 1.0, 0.0, NULL, 1, 1);
372
373 // check generated interpolation
374 {
375 double *** src_data = NULL;
376 double ** tgt_data = NULL;
377
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;
383 }
384
385 src_data = xmalloc(1 * sizeof(*src_data));
386 src_data[0] = xmalloc(1 * sizeof(**src_data));
387 src_data[0][0] =
388 xmalloc(src_grid_data.num_cells * 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]);
391
392 yac_interpolation_execute(interpolation, src_data, tgt_data);
393
394 if (comm_rank == 0) {
395
396 // check results
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]))
399 PUT_ERR("wrong results");
400
401 free(tgt_data[0]);
402 free(tgt_data);
403 }
404
405 free(src_data[0][0]);
406 free(src_data[0]);
407 free(src_data);
408 }
409 yac_interpolation_delete(interpolation);
411 yac_interp_method_delete(method_stack);
412 yac_interp_grid_delete(interp_grid);
413 yac_dist_grid_pair_delete(grid_pair);
414 yac_basic_grid_delete(tgt_grid);
415 yac_basic_grid_delete(src_grid);
416 }
417
418 // cleanup
421
422 return TEST_EXIT_CODE;
423}
424
425static void test1_callback(
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) {
429
430 UNUSED(tgt_coords);
431 UNUSED(src_cell_id);
432 UNUSED(src_cell_idx);
433
434 size_t num_src_fields = *(size_t*)user_data;
435
436 for (size_t i = 0; i < num_src_fields; ++i) {
437 global_results_points[i] = NULL;
438 result_weights[i] = NULL;
439 result_count[i] = 0;
440 }
441}
442
443static void test2_callback(
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) {
447
448 UNUSED(tgt_coords);
449
450 struct grid_info * info = (struct grid_info *)user_data;
451
452 // consistency check
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);
456 exit(EXIT_FAILURE);
457 }
458
459 static double vertex_weights[4] = {0.125, 0.125, 0.125, 0.125};
460 static double cell_weight[1] = {0.5};
461
462 static int vertex_result_points[4];
463 static int cell_result_points[1];
464
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;
469
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;
474 result_count[0] = 4;
475 result_count[1] = 1;
476}
477
478static void test3_callback(
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) {
482
483 UNUSED(tgt_coords);
484 UNUSED(src_cell_id);
485 UNUSED(src_cell_idx);
487
488 static int result_points[1];
489 static double weight[1] = {1.0};
490
491 result_points[0] = (int)src_cell_id;
492
493 global_results_points[0] = result_points;
494 result_weights[0] = weight;
495 result_count[0] = 1;
496}
497
498static void test4_callback(
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) {
502
503 UNUSED(tgt_coords);
504
505 struct grid_info * info = (struct grid_info *)user_data;
506
507 // consistency check
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);
511 exit(EXIT_FAILURE);
512 }
513
514 static double vertex_weights[4] = {0.125, 0.125, 0.125, 0.125};
515 static double cell_weight[1] = {0.5};
516
517 static int vertex_result_points[4];
518 static int cell_result_points[1];
519
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;
524
525 if (src_cell_id&1) {
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;
530 result_count[0] = 4;
531 result_count[1] = 1;
532 } else {
533 global_results_points[0] = NULL;
534 global_results_points[1] = NULL;
535 result_weights[0] = NULL;
536 result_weights[1] = NULL;
537 result_count[0] = 0;
538 result_count[1] = 0;
539 }
540}
struct yac_basic_grid * yac_basic_grid_new(char const *name, struct yac_basic_grid_data grid_data)
Definition basic_grid.c:50
struct yac_basic_grid * yac_basic_grid_empty_new(char const *name)
Definition basic_grid.c:63
void yac_basic_grid_delete(struct yac_basic_grid *grid)
Definition basic_grid.c:70
#define UNUSED(x)
Definition core.h:73
void yac_dist_grid_pair_delete(struct yac_dist_grid_pair *grid_pair)
Definition dist_grid.c:2313
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)
Definition dist_grid.c:2061
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)
#define YAC_RAD
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)
Definition interp_grid.c:31
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 * user_data
void yac_interp_method_callback_get_compute_weights_callback(char const *key, yac_func_compute_weights *compute_weights_callback, void **user_data)
char * key
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
@ YAC_LOC_CORNER
Definition location.h:15
@ YAC_LOC_CELL
Definition location.h:14
#define xmalloc(size)
Definition ppm_xfuncs.h:66
enum yac_location location
Definition basic_grid.h:16
struct yac_interp_field tgt_field
Definition interp_grid.c:26
size_t num_src_fields
Definition interp_grid.c:27
struct yac_dist_grid_pair * grid_pair
Definition interp_grid.c:25
struct yac_interp_field src_fields[]
Definition interp_grid.c:28
int * cell_to_vertex
int collection_size
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)
double coordinates_x[]
size_t num_cells[2]
double coordinates_y[]
#define TEST_EXIT_CODE
Definition tests.h:14
#define PUT_ERR(string)
Definition tests.h:10
@ NUM_TEST_FUNCTIONS
Definition toy_scrip.c:136
struct @1 test_functions[]
struct yac_basic_grid ** grids
Definition yac.c:152
void yac_mpi_finalize()
Definition yac_mpi.c:108
void yac_yaxt_init(MPI_Comm comm)
Definition yac_mpi.c:40
void yac_mpi_init()
Definition yac_mpi.c:77
YAC_INT yac_int
Definition yac_types.h:15