YAC 3.12.0
Yet Another Coupler
Loading...
Searching...
No Matches
test_interp_method_creep_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
7#include "tests.h"
8#include "test_common.h"
9#include "geometry.h"
15#include "dist_grid_utils.h"
16#include "yac_mpi.h"
17
18#include <mpi.h>
19#include <yaxt.h>
20#include <netcdf.h>
21
27static char const * grid_names[2] = {"src_grid", "tgt_grid"};
28
29int main(void) {
30
31 MPI_Init(NULL, NULL);
32
33 xt_initialize(MPI_COMM_WORLD);
34
35 int comm_rank, comm_size;
36 MPI_Comm_rank(MPI_COMM_WORLD, &comm_rank);
37 MPI_Comm_size(MPI_COMM_WORLD, &comm_size);
38 MPI_Barrier(MPI_COMM_WORLD);
39
40 if (comm_size != 3) {
41 PUT_ERR("ERROR: wrong number of processes");
42 xt_finalize();
43 MPI_Finalize();
44 return TEST_EXIT_CODE;
45 }
46
47 MPI_Comm split_comm;
48 MPI_Comm_split(MPI_COMM_WORLD, comm_rank < 2, 0, &split_comm);
49
50 int split_comm_rank, split_comm_size;
51 MPI_Comm_rank(split_comm, &split_comm_rank);
52 MPI_Comm_size(split_comm, &split_comm_size);
53
54 {// stack that starts with linear point interpolation, followed by creep
55 // fill and fixed as backup with two source processes and a single target
56 // process
57
58 // the global source grid is a 3x2 grid:
59 // 08--14--09--15--10--16--11
60 // | | | |
61 // 08 03 10 04 12 05 13
62 // | | | |
63 // 04--07--05--09--06--11--07
64 // | | | |
65 // 01 00 03 01 05 02 06
66 // | | | |
67 // 00--00--01--02--02--04--03
68 //
69 //---------------
70 // setup
71 //---------------
72
73 int is_tgt = split_comm_size == 1;
74 double coordinates_x[] = {0.0,1.0,2.0,3.0};
75 double coordinates_y[] = {0.0,1.0,2.0,3.0,4.0,5.0,6.0,7.0,8.0};
76 double tgt_cell_coordinates_x[] = {0.5,1.5,2.5};
77 double tgt_cell_coordinates_y[] = {0.5,1.5,2.5,3.5,4.5,5.5,6.5,7.5};
78 double tgt_cell_coordinates[24][3];
79 int tgt_cell_mask[24] = {1, 1, 1,
80 1, 1, 1,
81 1, 1, 1,
82 1, 0, 1,
83 1, 0, 0,
84 1, 1, 1,
85 0, 0, 0,
86 1, 1, 1};
87 size_t const num_cells[2][2] = {{3,2}, {3,8}};
88 size_t local_start[2][2][2] = {{{0,0},{1,0}}, {{0,0}}};
89 size_t local_count[2][2][2] = {{{2,2},{2,2}}, {{3,8}}};
90 int with_halo = 0;
91 for (size_t i = 0; i <= num_cells[is_tgt][0]; ++i)
93 for (size_t i = 0; i <= num_cells[is_tgt][1]; ++i)
95
99 local_start[is_tgt][split_comm_rank],
100 local_count[is_tgt][split_comm_rank], with_halo);
101
102 struct yac_basic_grid * grids[2] =
105
106 if (is_tgt) {
107 for (size_t i = 0, k = 0; i < num_cells[is_tgt][1]; ++i)
108 for (size_t j = 0; j < num_cells[is_tgt][0]; ++j, ++k)
110 tgt_cell_coordinates_x[j], tgt_cell_coordinates_y[i],
111 tgt_cell_coordinates[k]);
113 grids[0], YAC_LOC_CELL, tgt_cell_coordinates, grid_data.num_cells);
115 grids[0], YAC_LOC_CELL, tgt_cell_mask, grid_data.num_cells, NULL);
116 }
117
118 struct yac_dist_grid_pair * grid_pair =
119 yac_dist_grid_pair_new(grids[is_tgt], grids[is_tgt^1], MPI_COMM_WORLD);
120 struct yac_interp_field src_fields[] =
121 {{.location = YAC_LOC_CORNER, .coordinates_idx = SIZE_MAX, .masks_idx = SIZE_MAX}};
122 size_t num_src_fields = sizeof(src_fields) / sizeof(src_fields[0]);
124 {.location = YAC_LOC_CELL, .coordinates_idx = 0, .masks_idx = 0};
125
126 struct yac_interp_grid * interp_grid =
129
130 for (int creep_distance = -1; creep_distance <= 7; ++creep_distance) {
131
132 struct interp_method * method_stack[] =
134 yac_interp_method_creep_new(creep_distance),
135 yac_interp_method_fixed_new(-1.0), NULL};
136
137 struct yac_interp_weights * weights =
138 yac_interp_method_do_search(method_stack, interp_grid);
139
140 enum yac_interp_weights_reorder_type reorder_type[2] =
142
143 for (size_t i = 0; i < 2; ++i) {
144
145 struct yac_interpolation * interpolation =
147 weights, reorder_type[i], 1,
148 YAC_FRAC_MASK_NO_VALUE, 1.0, 0.0, NULL, 1, 1);
149
150 // check generated interpolation
151 {
152 double * src_field = NULL;
153 double ** src_fields = &src_field;
154 double * tgt_field = NULL;
155 double * ref_tgt_field = NULL;
156 double ref_global_tgt_field[9][24] =
157 {{ 2.5, 3.5, 4.5, // creep_distance = -1
158 6.5, 7.5, 8.5,
159 6.5, 7.5, 8.5,
160 6.5, -2.0, 8.5,
161 6.5, -2.0, -2.0,
162 6.5, 6.5, 6.5,
163 -2.0, -2.0, -2.0,
164 -1.0, -1.0, -1.0},
165 { 2.5, 3.5, 4.5, // creep_distance = 0
166 6.5, 7.5, 8.5,
167 -1.0, -1.0, -1.0,
168 -1.0, -2.0, -1.0,
169 -1.0, -2.0, -2.0,
170 -1.0, -1.0, -1.0,
171 -2.0, -2.0, -2.0,
172 -1.0, -1.0, -1.0},
173 { 2.5, 3.5, 4.5, // creep_distance = 1
174 6.5, 7.5, 8.5,
175 6.5, 7.5, 8.5,
176 -1.0, -2.0, -1.0,
177 -1.0, -2.0, -2.0,
178 -1.0, -1.0, -1.0,
179 -2.0, -2.0, -2.0,
180 -1.0, -1.0, -1.0},
181 { 2.5, 3.5, 4.5, // creep_distance = 2
182 6.5, 7.5, 8.5,
183 6.5, 7.5, 8.5,
184 6.5, -2.0, 8.5,
185 -1.0, -2.0, -2.0,
186 -1.0, -1.0, -1.0,
187 -2.0, -2.0, -2.0,
188 -1.0, -1.0, -1.0},
189 { 2.5, 3.5, 4.5, // creep_distance = 3
190 6.5, 7.5, 8.5,
191 6.5, 7.5, 8.5,
192 6.5, -2.0, 8.5,
193 6.5, -2.0, -2.0,
194 -1.0, -1.0, -1.0,
195 -2.0, -2.0, -2.0,
196 -1.0, -1.0, -1.0},
197 { 2.5, 3.5, 4.5, // creep_distance = 4
198 6.5, 7.5, 8.5,
199 6.5, 7.5, 8.5,
200 6.5, -2.0, 8.5,
201 6.5, -2.0, -2.0,
202 6.5, -1.0, -1.0,
203 -2.0, -2.0, -2.0,
204 -1.0, -1.0, -1.0},
205 { 2.5, 3.5, 4.5, // creep_distance = 5
206 6.5, 7.5, 8.5,
207 6.5, 7.5, 8.5,
208 6.5, -2.0, 8.5,
209 6.5, -2.0, -2.0,
210 6.5, 6.5, -1.0,
211 -2.0, -2.0, -2.0,
212 -1.0, -1.0, -1.0},
213 { 2.5, 3.5, 4.5, // creep_distance = 6
214 6.5, 7.5, 8.5,
215 6.5, 7.5, 8.5,
216 6.5, -2.0, 8.5,
217 6.5, -2.0, -2.0,
218 6.5, 6.5, 6.5,
219 -2.0, -2.0, -2.0,
220 -1.0, -1.0, -1.0},
221 { 2.5, 3.5, 4.5, // creep_distance = 7
222 6.5, 7.5, 8.5,
223 6.5, 7.5, 8.5,
224 6.5, -2.0, 8.5,
225 6.5, -2.0, -2.0,
226 6.5, 6.5, 6.5,
227 -2.0, -2.0, -2.0,
228 -1.0, -1.0, -1.0}};
229
230 if (is_tgt) {
231 tgt_field = xmalloc(grid_data.num_cells * sizeof(*tgt_field));
232 ref_tgt_field = xmalloc(grid_data.num_cells * sizeof(*ref_tgt_field));
233 for (size_t i = 0; i < grid_data.num_cells; ++i) {
234 tgt_field[i] = -2.0;
235 ref_tgt_field[i] =
236 ref_global_tgt_field[creep_distance+1][grid_data.cell_ids[i]];
237 }
238 } else {
239 src_field = xmalloc(grid_data.num_vertices * sizeof(*src_field));
240 for (size_t i = 0; i < grid_data.num_vertices; ++i)
241 src_field[i] = (double)(grid_data.vertex_ids[i]);
242 }
243
244 yac_interpolation_execute(interpolation, &src_fields, &tgt_field);
245
246 if (is_tgt)
247 for (size_t i = 0; i < grid_data.num_cells; ++i)
248 if (fabs(ref_tgt_field[i] - tgt_field[i]) > 1e-9)
249 PUT_ERR("wrong interpolation result");
250
251 free(src_field);
252 free(tgt_field);
253 free(ref_tgt_field);
254 }
255
256 yac_interpolation_delete(interpolation);
257 }
258
260 yac_interp_method_delete(method_stack);
261 }
262 yac_interp_grid_delete(interp_grid);
263 yac_dist_grid_pair_delete(grid_pair);
266 }
267
268 {// stack that starts with linear point interpolation, followed by creep
269 // fill and fixed as backup with two source processes and a single target
270 // process
271
272 // the global source grid is a 7x7 grid
273 // the global target grid is a 6x6 grid
274 //
275 //---------------
276 // setup
277 //---------------
278
279 int is_tgt = split_comm_size == 1;
280 double coordinates_x[2][8] = {{0.0,1.0,2.0,3.0,4.0,5.0,6.0,7.0},
281 {0.5,1.5,2.5,3.5,4.5,5.5,6.5}};
282 double coordinates_y[2][8] = {{0.0,1.0,2.0,3.0,4.0,5.0,6.0,7.0},
283 {0.5,1.5,2.5,3.5,4.5,5.5,6.5}};
284 int vertex_mask[2][8*8] =
285 {{1,1,1,1,1,1,1,1,
286 1,1,1,1,1,1,1,1,
287 1,1,0,0,0,0,1,1,
288 1,1,0,0,0,0,1,1,
289 1,1,0,0,0,0,1,1,
290 1,1,0,0,0,0,1,1,
291 1,1,1,1,1,1,1,1,
292 1,1,1,1,1,1,1,1},
293 {1,1,1,1,1,1,1,
294 1,1,1,1,1,1,1,
295 1,1,1,0,1,1,1,
296 1,1,0,1,0,1,1,
297 1,1,1,0,1,1,1,
298 1,1,1,1,1,1,1,
299 1,1,1,1,1,1,1}};
300 size_t const num_cells[2][2] = {{7,7}, {6,6}};
301 size_t local_start[2][2][2] = {{{0,0},{3,0}}, {{0,0}}};
302 size_t local_count[2][2][2] = {{{3,7},{4,7}}, {{6,6}}};
303 int with_halo = 0;
304 for (size_t i = 0; i <= num_cells[is_tgt][0]; ++i)
305 coordinates_x[is_tgt][i] *= YAC_RAD;
306 for (size_t i = 0; i <= num_cells[is_tgt][1]; ++i)
307 coordinates_y[is_tgt][i] *= YAC_RAD;
308
311 coordinates_x[is_tgt], coordinates_y[is_tgt], num_cells[is_tgt],
312 local_start[is_tgt][split_comm_rank],
313 local_count[is_tgt][split_comm_rank], with_halo);
314
315 struct yac_basic_grid * grids[2] =
318
319 int temp_vertex_mask[8*8];
320 for (size_t i = 0; i < grid_data.num_vertices; ++i)
321 temp_vertex_mask[i] = vertex_mask[is_tgt][grid_data.vertex_ids[i]];
323 grids[0], YAC_LOC_CORNER, temp_vertex_mask, grid_data.num_vertices, NULL);
324
325 struct yac_dist_grid_pair * grid_pair =
326 yac_dist_grid_pair_new(grids[is_tgt], grids[is_tgt^1], MPI_COMM_WORLD);
327 struct yac_interp_field src_fields[] =
328 {{.location = YAC_LOC_CORNER, .coordinates_idx = SIZE_MAX, .masks_idx = 0}};
329 size_t num_src_fields = sizeof(src_fields) / sizeof(src_fields[0]);
331 {.location = YAC_LOC_CORNER, .coordinates_idx = SIZE_MAX, .masks_idx = 0};
332
333 struct yac_interp_grid * interp_grid =
336
337 for (int creep_distance = -1; creep_distance <= 3; ++creep_distance) {
338
339 struct interp_method * method_stack[] =
341 yac_interp_method_creep_new(creep_distance),
342 yac_interp_method_fixed_new(-1.0), NULL};
343
344 struct yac_interp_weights * weights =
345 yac_interp_method_do_search(method_stack, interp_grid);
346
347 enum yac_interp_weights_reorder_type reorder_type[2] =
349
350 for (size_t i = 0; i < 2; ++i) {
351
352 struct yac_interpolation * interpolation =
354 weights, reorder_type[i], 1,
355 YAC_FRAC_MASK_NO_VALUE, 1.0, 0.0, NULL, 1, 1);
356
357 // check generated interpolation
358 {
359 double * src_field = NULL;
360 double ** src_fields = &src_field;
361 double * tgt_field = NULL;
362 double * ref_tgt_field = NULL;
363 double ref_global_tgt_field[5][7*7] =
364 {{ 4.5, 5.5, 6.5, 7.5, 8.5, 9.5, 10.5, // creep_distance = -1
365 12.5, 9.0, 6.5, 7.5, 8.5, 14.0, 18.5,
366 20.5, 20.5, 13.5, -2.0, 17.5, 26.5, 26.5,
367 28.5, 28.5, -2.0, -1.0, -2.0, 34.5, 34.5,
368 36.5, 36.5, 45.5, -2.0, 49.5, 42.5, 42.5,
369 44.5, 49.0, 54.5, 55.5, 56.5, 54.0, 50.5,
370 52.5, 53.5, 54.5, 55.5, 56.5, 57.5, 58.5},
371 { 4.5, 5.5, 6.5, 7.5, 8.5, 9.5, 10.5, // creep_distance = 0
372 12.5, -1.0, -1.0, -1.0, -1.0, -1.0, 18.5,
373 20.5, -1.0, -1.0, -2.0, -1.0, -1.0, 26.5,
374 28.5, -1.0, -2.0, -1.0, -2.0, -1.0, 34.5,
375 36.5, -1.0, -1.0, -2.0, -1.0, -1.0, 42.5,
376 44.5, -1.0, -1.0, -1.0, -1.0, -1.0, 50.5,
377 52.5, 53.5, 54.5, 55.5, 56.5, 57.5, 58.5},
378 { 4.5, 5.5, 6.5, 7.5, 8.5, 9.5, 10.5, // creep_distance = 1
379 12.5, 9.0, 6.5, 7.5, 8.5, 14.0, 18.5,
380 20.5, 20.5, -1.0, -2.0, -1.0, 26.5, 26.5,
381 28.5, 28.5, -2.0, -1.0, -2.0, 34.5, 34.5,
382 36.5, 36.5, -1.0, -2.0, -1.0, 42.5, 42.5,
383 44.5, 49.0, 54.5, 55.5, 56.5, 54.0, 50.5,
384 52.5, 53.5, 54.5, 55.5, 56.5, 57.5, 58.5},
385 { 4.5, 5.5, 6.5, 7.5, 8.5, 9.5, 10.5, // creep_distance = 2
386 12.5, 9.0, 6.5, 7.5, 8.5, 14.0, 18.5,
387 20.5, 20.5, 13.5, -2.0, 17.5, 26.5, 26.5,
388 28.5, 28.5, -2.0, -1.0, -2.0, 34.5, 34.5,
389 36.5, 36.5, 45.5, -2.0, 49.5, 42.5, 42.5,
390 44.5, 49.0, 54.5, 55.5, 56.5, 54.0, 50.5,
391 52.5, 53.5, 54.5, 55.5, 56.5, 57.5, 58.5},
392 { 4.5, 5.5, 6.5, 7.5, 8.5, 9.5, 10.5, // creep_distance = 3
393 12.5, 9.0, 6.5, 7.5, 8.5, 14.0, 18.5,
394 20.5, 20.5, 13.5, -2.0, 17.5, 26.5, 26.5,
395 28.5, 28.5, -2.0, -1.0, -2.0, 34.5, 34.5,
396 36.5, 36.5, 45.5, -2.0, 49.5, 42.5, 42.5,
397 44.5, 49.0, 54.5, 55.5, 56.5, 54.0, 50.5,
398 52.5, 53.5, 54.5, 55.5, 56.5, 57.5, 58.5}};
399
400 if (is_tgt) {
401 tgt_field = xmalloc(grid_data.num_vertices * sizeof(*tgt_field));
402 ref_tgt_field = xmalloc(grid_data.num_vertices * sizeof(*ref_tgt_field));
403 for (size_t i = 0; i < grid_data.num_vertices; ++i) {
404 tgt_field[i] = -2.0;
405 ref_tgt_field[i] =
406 ref_global_tgt_field[creep_distance+1][grid_data.vertex_ids[i]];
407 }
408 } else {
409 src_field = xmalloc(grid_data.num_vertices * sizeof(*src_field));
410 for (size_t i = 0; i < grid_data.num_vertices; ++i)
411 src_field[i] = (double)(grid_data.vertex_ids[i]);
412 }
413
414 yac_interpolation_execute(interpolation, &src_fields, &tgt_field);
415
416 if (is_tgt)
417 for (size_t i = 0; i < grid_data.num_vertices; ++i)
418 if (fabs(ref_tgt_field[i] - tgt_field[i]) > 1e-9)
419 PUT_ERR("wrong interpolation result");
420
421 free(src_field);
422 free(tgt_field);
423 free(ref_tgt_field);
424 }
425
426 yac_interpolation_delete(interpolation);
427 }
428
430 yac_interp_method_delete(method_stack);
431 }
432 yac_interp_grid_delete(interp_grid);
433 yac_dist_grid_pair_delete(grid_pair);
436 }
437
438 {// stack that starts with linear point interpolation, followed by creep
439 // fill and fixed as backup with two source processes and a single target
440 // process
441
442 // the global source grid is a 5x1 grid
443 // the global target grid is a 5x3 grid
444 //
445 //---------------
446 // setup
447 //---------------
448
449 int is_tgt = split_comm_size == 1;
450 double coordinates_x[6] = {0.0,1.0,2.0,3.0,4.0,5.0};
451 double coordinates_y[4] = {0.0,1.0,2.0,3.0};
452 int cell_mask[2][5*3] =
453 {{1,0,1,0,1,},
454 {1,1,1,1,1,
455 1,1,0,1,1,
456 1,1,1,1,1}};
457 size_t const num_cells[2][2] = {{5,1}, {5,3}};
458 size_t local_start[2][2][2] = {{{0,0},{3,0}}, {{0,0}}};
459 size_t local_count[2][2][2] = {{{3,1},{2,1}}, {{5,3}}};
460 int with_halo = 0;
461 for (size_t i = 0; i <= num_cells[is_tgt][0]; ++i)
463 for (size_t i = 0; i <= num_cells[is_tgt][1]; ++i)
465
469 local_start[is_tgt][split_comm_rank],
470 local_count[is_tgt][split_comm_rank], with_halo);
471
472 struct yac_basic_grid * grids[2] =
475
476 int temp_cell_mask[5*3];
477 for (size_t i = 0; i < grid_data.num_cells; ++i)
478 temp_cell_mask[i] = cell_mask[is_tgt][grid_data.cell_ids[i]];
480 grids[0], YAC_LOC_CELL, temp_cell_mask, grid_data.num_cells, NULL);
481
482 struct yac_dist_grid_pair * grid_pair =
483 yac_dist_grid_pair_new(grids[is_tgt], grids[is_tgt^1], MPI_COMM_WORLD);
484 struct yac_interp_field src_fields[] =
485 {{.location = YAC_LOC_CELL, .coordinates_idx = SIZE_MAX, .masks_idx = 0}};
486 size_t num_src_fields = sizeof(src_fields) / sizeof(src_fields[0]);
488 {.location = YAC_LOC_CELL, .coordinates_idx = SIZE_MAX, .masks_idx = 0};
489
490 struct yac_interp_grid * interp_grid =
493
494 for (int creep_distance = -1; creep_distance <= 5; ++creep_distance) {
495
496 struct interp_method * method_stack[] =
498 yac_interp_method_creep_new(creep_distance),
499 yac_interp_method_fixed_new(-1.0), NULL};
500
501 struct yac_interp_weights * weights =
502 yac_interp_method_do_search(method_stack, interp_grid);
503
504 enum yac_interp_weights_reorder_type reorder_type[2] =
506
507 for (size_t i = 0; i < 2; ++i) {
508
509 struct yac_interpolation * interpolation =
511 weights, reorder_type[i], 1,
512 YAC_FRAC_MASK_NO_VALUE, 1.0, 0.0, NULL, 1, 1);
513
514 // check generated interpolation
515 {
516 double * src_field = NULL;
517 double ** src_fields = &src_field;
518 double * tgt_field = NULL;
519 double * ref_tgt_field = NULL;
520 double ref_global_tgt_field[7][7*7] =
521 {{ 0.0, 1.0 , 2.0 , 3.0 , 4.0, // creep_distance = -1
522 0.0, 0.5 , -2.0 , 3.5 , 4.0,
523 0.0, 0.25, 2.0 , 3.75, 4.0},
524 { 0.0, -1.0 , 2.0 , -1.0 , 4.0, // creep_distance = 0
525 -1.0, -1.0 , -2.0 , -1.0 , -1.0,
526 -1.0, -1.0 , -1.0 , -1.0 , -1.0},
527 { 0.0, 1.0 , 2.0 , 3.0 , 4.0, // creep_distance = 1
528 0.0, -1.0 , -2.0 , -1.0 , 4.0,
529 -1.0, -1.0 , -1.0 , -1.0 , -1.0},
530 { 0.0, 1.0 , 2.0 , 3.0 , 4.0, // creep_distance = 2
531 0.0, 0.5 , -2.0 , 3.5 , 4.0,
532 0.0, -1.0 , -1.0 , -1.0 , 4.0},
533 { 0.0, 1.0 , 2.0 , 3.0 , 4.0, // creep_distance = 3
534 0.0, 0.5 , -2.0 , 3.5 , 4.0,
535 0.0, 0.25, -1.0 , 3.75, 4.0},
536 { 0.0, 1.0 , 2.0 , 3.0 , 4.0, // creep_distance = 4
537 0.0, 0.5 , -2.0 , 3.5 , 4.0,
538 0.0, 0.25, 2.0 , 3.75, 4.0},
539 { 0.0, 1.0 , 2.0 , 3.0 , 4.0, // creep_distance = 5
540 0.0, 0.5 , -2.0 , 3.5 , 4.0,
541 0.0, 0.25, 2.0 , 3.75, 4.0}};
542
543 if (is_tgt) {
544 tgt_field = xmalloc(grid_data.num_cells* sizeof(*tgt_field));
545 ref_tgt_field = xmalloc(grid_data.num_cells* sizeof(*ref_tgt_field));
546 for (size_t i = 0; i < grid_data.num_cells; ++i) {
547 tgt_field[i] = -2.0;
548 ref_tgt_field[i] =
549 ref_global_tgt_field[creep_distance+1][grid_data.cell_ids[i]];
550 }
551 } else {
552 src_field = xmalloc(grid_data.num_cells* sizeof(*src_field));
553 for (size_t i = 0; i < grid_data.num_cells; ++i)
554 src_field[i] = (double)(grid_data.cell_ids[i]);
555 }
556
557 yac_interpolation_execute(interpolation, &src_fields, &tgt_field);
558
559 if (is_tgt)
560 for (size_t i = 0; i < grid_data.num_cells; ++i)
561 if (fabs(ref_tgt_field[i] - tgt_field[i]) > 1e-9)
562 PUT_ERR("wrong interpolation result");
563
564 free(src_field);
565 free(tgt_field);
566 free(ref_tgt_field);
567 }
568
569 yac_interpolation_delete(interpolation);
570 }
571
573 yac_interp_method_delete(method_stack);
574 }
575 yac_interp_grid_delete(interp_grid);
576 yac_dist_grid_pair_delete(grid_pair);
579 }
580
581 MPI_Comm_free(&split_comm);
582
583 xt_finalize();
584
585 MPI_Finalize();
586
587 return TEST_EXIT_CODE;
588}
struct yac_basic_grid * yac_basic_grid_new(char const *name, struct yac_basic_grid_data grid_data)
Definition basic_grid.c:50
size_t yac_basic_grid_add_mask(struct yac_basic_grid *grid, enum yac_location location, int const *mask, size_t count, char const *mask_name)
Definition basic_grid.c:284
size_t yac_basic_grid_add_coordinates(struct yac_basic_grid *grid, enum yac_location location, yac_coordinate_pointer coordinates, size_t count)
Definition basic_grid.c:232
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
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
static void LLtoXYZ_deg(double lon, double lat, double p_out[])
Definition geometry.h:269
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)
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_conserv_new(int order, int enforced_conserv, int partial_coverage, enum yac_interp_method_conserv_normalisation normalisation)
@ YAC_INTERP_CONSERV_DESTAREA
struct interp_method * yac_interp_method_creep_new(int creep_distance)
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_mask
static char const * grid_names[2]
static MPI_Comm split_comm
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
struct yac_basic_grid ** grids
Definition yac.c:152