YAC 3.12.0
Yet Another Coupler
Loading...
Searching...
No Matches
test_interp_method_rbf_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 <stdio.h>
7#include <math.h>
8
9#include "tests.h"
12#include "dist_grid_utils.h"
13#include "yac_mpi.h"
14#include "geometry.h"
15
16#include <mpi.h>
17#include <yaxt.h>
18#include <netcdf.h>
19
25char const src_grid_name[] = "src_grid";
26char const tgt_grid_name[] = "tgt_grid";
27
28static void utest_target_main(MPI_Comm target_comm);
29static void utest_source_main(MPI_Comm source_comm);
30
31int main (void) {
32
33 MPI_Init(NULL, NULL);
34
35 xt_initialize(MPI_COMM_WORLD);
36
37 int comm_rank, comm_size;
38 MPI_Comm_rank(MPI_COMM_WORLD, &comm_rank);
39 MPI_Comm_size(MPI_COMM_WORLD, &comm_size);
40 MPI_Barrier(MPI_COMM_WORLD);
41
42 if (comm_size != 2) {
43 PUT_ERR("ERROR: wrong number of processes");
44 xt_finalize();
45 MPI_Finalize();
46 return TEST_EXIT_CODE;
47 }
48
49 int tgt_flag = comm_rank < 1;
50
51 MPI_Comm split_comm;
52 MPI_Comm_split(MPI_COMM_WORLD, tgt_flag, 0, &split_comm);
53
54 if (tgt_flag) utest_target_main(split_comm);
55 else utest_source_main(split_comm);
56
57 MPI_Comm_free(&split_comm);
58 xt_finalize();
59 MPI_Finalize();
60
61 return TEST_EXIT_CODE;
62}
63
64static void utest_source_main(MPI_Comm source_comm) {
65
66 // 1 and 5 nearest neighbours per target point with fixed value fallback
67
68 // corner and cell ids for a 7 x 7 grid (x == target point position)
69
70 // 56-----57-----58-----59-----60-----61-----62-----63
71 // | x| x| x| x| x| x| x|
72 // | 42 | 43 | 44 | 45 | 46 | 47 | 48 |
73 // | | | | | | | |
74 // 48-----49-----50-----51-----52-----53-----54-----55
75 // | x| x| x| x| x| x| x|
76 // | 35 | 36 | 37 | 38 | 39 | 40 | 41 |
77 // | | | | | | | |
78 // 40-----41-----42-----43-----44-----45-----46-----47
79 // | x| x| x| x| x| x| x|
80 // | 28 | 29 | 30 | 31 | 32 | 33 | 34 |
81 // | | | | | | | |
82 // 32-----33-----34-----35-----36-----37-----38-----39
83 // | x| x| x| x| x| x| x|
84 // | 21 | 22 | 23 | 24 | 25 | 26 | 27 |
85 // | | | | | | | |
86 // 24-----25-----26-----27-----28-----29-----30-----31
87 // | x| x| x| x| x| x| x|
88 // | 14 | 15 | 16 | 17 | 18 | 19 | 20 |
89 // | | | | | | | |
90 // 16-----17-----18-----19-----20-----21-----22-----23
91 // | x| x| x| x| x| x| x|
92 // | 07 | 08 | 09 | 10 | 11 | 12 | 13 |
93 // | | | | | | | |
94 // 08-----09-----10-----11-----12-----13-----14-----15
95 // | x| x| x| x| x| x| x|
96 // | 00 | 01 | 02 | 03 | 04 | 05 | 06 |
97 // | | | | | | | |
98 // 00-----01-----02-----03-----04-----05-----06-----07
99 //
100 // the grid is distributed among the processes as follows:
101 // (index == process)
102 //
103 // 3---3---3---3---3---3---3---3
104 // | 3 | 3 | 3 | 3 | 3 | 3 | 3 |
105 // 3---3---3---3---3---3---3---3
106 // | 3 | 3 | 3 | 3 | 3 | 3 | 3 |
107 // 3---3---3---1---2---2---3---3
108 // | 1 | 1 | 1 | 2 | 2 | 2 | 2 |
109 // 1---1---1---2---2---2---2---2
110 // | 1 | 1 | 1 | 2 | 2 | 2 | 2 |
111 // 1---1---1---1---2---2---2---2
112 // | 1 | 1 | 1 | 2 | 2 | 2 | 2 |
113 // 1---1---1---0---0---0---2---2
114 // | 0 | 0 | 0 | 0 | 0 | 0 | 0 |
115 // 0---0---0---0---0---0---0---0
116 // | 0 | 0 | 0 | 0 | 0 | 0 | 0 |
117 // 0---0---0---0---0---0---0---0
118 //
119 // the source mask looks as follows (# == masked point)
120 //
121 // +---+---+---+---+---+---+---#
122 // | | | | | | | |
123 // +---+---+---+---+---+---#---+
124 // | | | | | | | |
125 // +---+---+---+---+---#---+---+
126 // | | | | | | | |
127 // +---+---+---+---#---+---+---+
128 // | | | | | | | |
129 // +---+---+---#---+---+---+---+
130 // | | | | | | | |
131 // +---+---#---+---+---+---+---+
132 // | | | | | | | |
133 // +---#---+---+---+---+---+---+
134 // | | | | | | | |
135 // #---+---+---+---+---+---+---+
136
137 int my_source_rank;
138 MPI_Comm_rank(source_comm, &my_source_rank);
139
140 double coordinates_x[] = {-2.0, -1.0, 0.0, 1.0, 2.0};
141 double coordinates_y[] = {-2.0, -1.0, 0.0, 1.0, 2.0};
142 size_t const num_cells[2] = {4,4};
143 size_t local_start[2] = {0,0};
144 size_t local_count[2] = {4,4};
145 int with_halo = 0;
146 for (size_t i = 0; i <= num_cells[0]; ++i) coordinates_x[i] *= YAC_RAD;
147 for (size_t i = 0; i <= num_cells[1]; ++i) coordinates_y[i] *= YAC_RAD;
148
152 local_start, local_count, with_halo);
153
154 struct yac_basic_grid * grids[2] =
158
159 struct yac_dist_grid_pair * grid_pair =
160 yac_dist_grid_pair_new(grids[0], grids[1], MPI_COMM_WORLD);
161
162 struct yac_interp_field src_fields[] =
163 {{.location = YAC_LOC_CORNER, .coordinates_idx = SIZE_MAX, .masks_idx = SIZE_MAX}};
164 size_t num_src_fields = sizeof(src_fields) / sizeof(src_fields[0]);
166 {.location = YAC_LOC_CELL, .coordinates_idx = 0, .masks_idx = SIZE_MAX};
167
168 struct yac_interp_grid * interp_grid =
171
172 struct interp_method * method_stack[2] =
174 (struct yac_nnn_config){
175 .type = YAC_INTERP_NNN_RBF, .n = 9,
176 .max_search_distance = YAC_INTERP_RBF_MAX_SEARCH_DISTANCE_DEFAULT,
177 .data.rbf_scale = YAC_INTERP_RBF_SCALE_DEFAULT}),
178 NULL};
179
180 struct yac_interp_weights * weights =
181 yac_interp_method_do_search(method_stack, interp_grid);
182
183 yac_interp_method_delete(method_stack);
184
185 enum yac_interp_weights_reorder_type reorder_type[2] =
187
188 for (size_t i = 0; i < 2; ++i) {
189 for (size_t collection_size = 1; collection_size < 4;
190 collection_size += 2) {
191
192 struct yac_interpolation * interpolation =
194 weights, reorder_type[i], collection_size,
195 YAC_FRAC_MASK_NO_VALUE, 1.0, 0.0, NULL, 1, 1);
196
197 // check generated interpolation
198 {
199 double *** src_data = xmalloc(collection_size * sizeof(*src_data));
200 double ref_src_data[25] = {0,1,2,3,4,
201 1,2,3,4,5,
202 2,3,4,5,6,
203 3,4,5,6,7,
204 4,5,6,7,8};
205
206 for (size_t collection_idx = 0; collection_idx < collection_size;
207 ++collection_idx) {
208 // only one field
209 src_data[collection_idx] = xmalloc(1 * sizeof(**src_data));
210 src_data[collection_idx][0] =
211 xmalloc(grid_data.num_vertices * sizeof(***src_data));
212 for (size_t i = 0; i < grid_data.num_vertices; ++i)
213 src_data[collection_idx][0][i] =
214 ref_src_data[i] + (double)(collection_idx * 9);
215 }
216
217 yac_interpolation_execute(interpolation, src_data, NULL);
218
219 for (size_t collection_idx = 0; collection_idx < collection_size;
220 ++collection_idx) {
221 free(src_data[collection_idx][0]);
222 free(src_data[collection_idx]);
223 }
224 free(src_data);
225 }
226
227 yac_interpolation_delete(interpolation);
228 }
229 }
230
232 yac_interp_grid_delete(interp_grid);
233 yac_dist_grid_pair_delete(grid_pair);
236}
237
238static void utest_target_main(MPI_Comm target_comm) {
239
240 UNUSED(target_comm);
241
242 double coordinates_x[] = {-1.5,-0.5,0.5,1.5};
243 double coordinates_y[] = {-1.5,-0.5,0.5,1.5};
244 double cell_coordinates_x[] = {-1.0,0.0,1.0};
245 double cell_coordinates_y[] = {-1.0,0.0,1.0};
246 yac_coordinate_pointer cell_coordinates = xmalloc(9 * sizeof(*cell_coordinates));
247 size_t const num_cells[2] = {3,3};
248 size_t local_start[2] = {0,0};
249 size_t local_count[2] = {3,3};
250 int with_halo = 0;
251 for (size_t i = 0; i <= num_cells[0]; ++i) coordinates_x[i] *= YAC_RAD;
252 for (size_t i = 0; i <= num_cells[1]; ++i) coordinates_y[i] *= YAC_RAD;
253 for (size_t i = 0, k = 0; i < num_cells[1]; ++i)
254 for (size_t j = 0; j < num_cells[0]; ++j, ++k)
256 cell_coordinates_x[j], cell_coordinates_y[i], cell_coordinates[k]);
257
261 local_start, local_count, with_halo);
262
263 struct yac_basic_grid * grids[2] =
267 grids[0], YAC_LOC_CELL, cell_coordinates);
268
269 struct yac_dist_grid_pair * grid_pair =
270 yac_dist_grid_pair_new(grids[0], grids[1], MPI_COMM_WORLD);
271
272 struct yac_interp_field src_fields[] =
273 {{.location = YAC_LOC_CORNER, .coordinates_idx = SIZE_MAX, .masks_idx = SIZE_MAX}};
274 size_t num_src_fields = sizeof(src_fields) / sizeof(src_fields[0]);
276 {.location = YAC_LOC_CELL, .coordinates_idx = 0, .masks_idx = SIZE_MAX};
277
278 struct yac_interp_grid * interp_grid =
281
282 struct interp_method * method_stack[2] =
284 (struct yac_nnn_config){
285 .type = YAC_INTERP_NNN_RBF, .n = 9,
286 .max_search_distance = YAC_INTERP_RBF_MAX_SEARCH_DISTANCE_DEFAULT,
287 .data.rbf_scale = YAC_INTERP_RBF_SCALE_DEFAULT}),
288 NULL};
289
290 struct yac_interp_weights * weights =
291 yac_interp_method_do_search(method_stack, interp_grid);
292
293 yac_interp_method_delete(method_stack);
294
295 enum yac_interp_weights_reorder_type reorder_type[2] =
297
298 for (size_t i = 0; i < 2; ++i) {
299 for (size_t collection_size = 1; collection_size < 4;
300 collection_size += 2) {
301
302 struct yac_interpolation * interpolation =
304 weights, reorder_type[i], collection_size,
305 YAC_FRAC_MASK_NO_VALUE, 1.0, 0.0, NULL, 1, 1);
306
307 // check generated interpolation
308 {
309 double ref_tgt_field[9] = {2,3,4, 3,4,5, 4,5,6};
310
311 double ** tgt_data = xmalloc(collection_size * sizeof(*tgt_data));
312 for (size_t collection_idx = 0; collection_idx < collection_size;
313 ++collection_idx) {
314 tgt_data[collection_idx] =
315 xmalloc(grid_data.num_cells * sizeof(**tgt_data));
316 for (size_t j = 0; j < 9; ++j) tgt_data[collection_idx][j] = -1.0;
317 }
318
319 yac_interpolation_execute(interpolation, NULL, tgt_data);
320
321 for (size_t collection_idx = 0; collection_idx < collection_size;
322 ++collection_idx)
323 for (size_t j = 0, offset = collection_idx * 9;
324 j < grid_data.num_cells; ++j)
325 if (fabs((ref_tgt_field[j] + (double)offset) -
326 tgt_data[collection_idx][j]) > 1e-3)
327 PUT_ERR("wrong interpolation result");
328
329 for (size_t collection_idx = 0; collection_idx < collection_size;
330 ++collection_idx)
331 free(tgt_data[collection_idx]);
332 free(tgt_data);
333 }
334
335 yac_interpolation_delete(interpolation);
336 }
337 }
338
340 yac_interp_grid_delete(interp_grid);
341 yac_dist_grid_pair_delete(grid_pair);
344}
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_coordinates_nocpy(struct yac_basic_grid *grid, enum yac_location location, yac_coordinate_pointer coordinates)
Definition basic_grid.c:208
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
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_nnn_new(struct yac_nnn_config config)
@ YAC_INTERP_NNN_RBF
radial basis functions
#define YAC_INTERP_RBF_MAX_SEARCH_DISTANCE_DEFAULT
#define YAC_INTERP_RBF_SCALE_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)
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 collection_size
static MPI_Comm split_comm
char const src_grid_name[]
char const tgt_grid_name[]
double coordinates_x[]
size_t num_cells[2]
double cell_coordinates_y[]
double coordinates_y[]
double cell_coordinates_x[]
#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
double(* yac_coordinate_pointer)[3]
Definition yac_types.h:19