YAC 3.13.0
Yet Another Coupler
Loading...
Searching...
No Matches
test_interpolation_parallel1_c.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 <string.h>
8
9#include <mpi.h>
10#include <yaxt.h>
11
12#include "tests.h"
13#include "test_common.h"
14#include "dist_grid_utils.h"
17
23#define YAC_RAD (0.01745329251994329576923690768489) // M_PI / 180
24
25#define MAX_COLLECTION_SIZE (10)
26
27char const grid_name_src[] = "src_grid";
28char const grid_name_tgt[] = "tgt_grid";
29
30static void utest_target_main(
31 MPI_Comm target_comm, enum yac_interp_weights_reorder_type reorder_type);
32
33static void utest_source_main(
34 MPI_Comm source_comm, enum yac_interp_weights_reorder_type reorder_type);
35
36int main(int argc, char *argv[]) {
37
38 if (argc != 2) {
39 PUT_ERR("wrong number of arguments\n");
40 return TEST_EXIT_CODE;
41 }
42
43 enum yac_interp_weights_reorder_type reorder_type =
44 (strcmp(argv[1], "src") == 0)?YAC_MAPPING_ON_SRC:YAC_MAPPING_ON_TGT;
45
46 if ((reorder_type != YAC_MAPPING_ON_SRC) && strcmp(argv[1], "tgt")) {
47 PUT_ERR("invalid argument (has to be either \"src\" or \"tgt\")\n");
48 return TEST_EXIT_CODE;
49 }
50
51 MPI_Init(NULL, NULL);
52
53 xt_initialize(MPI_COMM_WORLD);
54
55 int comm_rank, comm_size;
56 MPI_Comm_rank(MPI_COMM_WORLD, &comm_rank);
57 MPI_Comm_size(MPI_COMM_WORLD, &comm_size);
58 MPI_Barrier(MPI_COMM_WORLD);
59
60 if (comm_size != 4) {
61 PUT_ERR("ERROR: wrong number of processes");
62 xt_finalize();
63 MPI_Finalize();
64 return TEST_EXIT_CODE;
65 }
66
67 // split processes into source and target
68
69 int tgt_flag = comm_rank < 2;
70
71 MPI_Comm split_comm;
72 MPI_Comm_split(MPI_COMM_WORLD, tgt_flag, 0, &split_comm);
73
74 if (tgt_flag) utest_target_main(split_comm, reorder_type);
75 else utest_source_main(split_comm, reorder_type);
76
77 MPI_Comm_free(&split_comm);
78 xt_finalize();
79 MPI_Finalize();
80
81 return TEST_EXIT_CODE;
82}
83
84/*
85 * The source grid is distributed among 2 processes
86 *
87 * The global source grid has 5x4 cells:
88 *
89 * 24--44--25--45--26--46--27--47--28--48--29
90 * | | | | | |
91 * 34 15 36 16 38 17 40 18 42 19 43
92 * | | | | | |
93 * 18--33--19--35--20--37--21--39--22--41--23
94 * | | | | | |
95 * 23 10 25 11 27 12 29 13 31 14 32
96 * | | | | | |
97 * 12--22--13--24--14--26--15--28--16--30--17
98 * | | | | | |
99 * 12 05 14 06 16 07 18 08 20 09 21
100 * | | | | | |
101 * 06--11--07--13--08--15--09--17--10--19--11
102 * | | | | | |
103 * 01 00 03 01 05 02 07 03 09 04 10
104 * | | | | | |
105 * 00--01--01--02--02--04--03--06--04--08--05
106 */
107static void utest_source_main(
108 MPI_Comm source_comm, enum yac_interp_weights_reorder_type reorder_type) {
109
110 int my_source_rank;
111 MPI_Comm_rank(source_comm, &my_source_rank);
112
113 double coordinates_x[] = {0.0, 1.0, 2.0, 3.0, 4.0, 5.0};
114 double coordinates_y[] = {0.0, 1.0, 2.0, 3.0, 4.0};
115 size_t const num_cells[2] = {5,4};
116 size_t local_start[2][2] = {{0,0},{3,0}};
117 size_t local_count[2][2] = {{3,4},{2,4}};
118 int with_halo = 1;
119 for (size_t i = 0; i <= num_cells[0]; ++i) coordinates_x[i] *= YAC_RAD;
120 for (size_t i = 0; i <= num_cells[1]; ++i) coordinates_y[i] *= YAC_RAD;
121
125 local_start[my_source_rank], local_count[my_source_rank], with_halo);
126 struct yac_basic_grid * src_grid =
128 struct yac_basic_grid * tgt_grid =
130
131 struct yac_dist_grid_pair * grid_pair =
132 yac_dist_grid_pair_new(src_grid, tgt_grid, MPI_COMM_WORLD);
133
134 struct yac_interp_field src_fields[] =
135 {{.location = YAC_LOC_CORNER, .coordinates_idx = SIZE_MAX, .masks_idx = SIZE_MAX}};
136 size_t num_src_fields = sizeof(src_fields) / sizeof(src_fields[0]);
138 {.location = YAC_LOC_CORNER, .coordinates_idx = SIZE_MAX, .masks_idx = SIZE_MAX};
139
140 struct yac_interp_grid * interp_grid =
143
144 struct interp_method * method_stack[] =
147 NULL};
148
149 struct yac_interp_weights * weights =
150 yac_interp_method_do_search(method_stack, interp_grid);
151
152 yac_interp_method_delete(method_stack);
153
154 // -------------------
155 // set up source data
156 // -------------------
157
158 // src_data dimensions [collection_idx]
159 // [pointset_idx]
160 // [local_idx]
161 double *** src_data = xmalloc(MAX_COLLECTION_SIZE * sizeof(*src_data));
162 for (size_t collection_idx = 0; collection_idx < MAX_COLLECTION_SIZE;
163 ++collection_idx) {
164 src_data[collection_idx] = xmalloc(1 * sizeof(**src_data));
165 src_data[collection_idx][0] =
166 xmalloc(grid_data.num_vertices * sizeof(***src_data));
167 for (size_t i = 0; i < grid_data.num_vertices; ++i)
168 src_data[collection_idx][0][i] =
169 (grid_data.core_vertex_mask[i])?
170 ((double)(grid_data.vertex_ids[i]) + (double)(collection_idx * 30)):
171 (-1.0);
172 }
173
175 ++collection_size) {
176
177 //---------------------------
178 // set up field interpolation
179 //---------------------------
180
181 struct yac_interpolation * interpolation =
183 weights, reorder_type, collection_size,
184 YAC_FRAC_MASK_NO_VALUE, 1.0, 0.0, NULL, 1, 1);
185 struct yac_interpolation * interpolation_copy =
186 yac_interpolation_copy(interpolation);
187
188 //---------------------
189 // do the interpolation
190 //---------------------
191
192 yac_interpolation_execute_put(interpolation, src_data);
193 yac_interpolation_execute_put(interpolation_copy, src_data);
194
195 yac_interpolation_delete(interpolation_copy);
196 yac_interpolation_delete(interpolation);
197 }
198
199 //--------
200 // cleanup
201 //--------
202
203
204 for (size_t collection_idx = 0; collection_idx < MAX_COLLECTION_SIZE;
205 ++collection_idx) {
206 free(src_data[collection_idx][0]);
207 free(src_data[collection_idx]);
208 }
209 free(src_data);
210
211 yac_basic_grid_delete(tgt_grid);
212 yac_basic_grid_delete(src_grid);
214 yac_interp_grid_delete(interp_grid);
215 yac_dist_grid_pair_delete(grid_pair);
216}
217
218/*
219 * The target grid is distributed among 2 processes
220 *
221 * The global target grid has 6x3 cells:
222 *
223 * 21--39--22--40--23--41--24--42--25--43--26--44--27
224 * | | | | | | |
225 * 27 12 29 13 31 14 33 15 35 16 37 17 38
226 * | | | | | | |
227 * 14--26--15--28--16--30--17--32--18--34--19--36--20
228 * | | | | | | |
229 * 14 06 16 07 18 08 20 09 22 10 24 11 25
230 * | | | | | | |
231 * 07--13--08--15--09--17--10--19--11--21--12--23--13
232 * | | | | | | |
233 * 01 00 03 01 05 02 07 03 09 04 11 05 12
234 * | | | | | | |
235 * 00--00--01--02--02--04--03--06--04--08--05--10--06
236 */
237static void utest_target_main(
238 MPI_Comm target_comm, enum yac_interp_weights_reorder_type reorder_type) {
239
240 int my_target_rank;
241 MPI_Comm_rank(target_comm, &my_target_rank);
242
243 double coordinates_x[] = {0.5,1.5,2.5,3.5,4.5,5.5,6.5};
244 double coordinates_y[] = {0.5,1.5,2.5,3.5};
245 size_t const num_cells[2] = {6,3};
246 size_t local_start[2][2] = {{0,0},{3,0}};
247 size_t local_count[2][2] = {{3,3},{3,3}};
248 int with_halo = 0;
249 for (size_t i = 0; i <= num_cells[0]; ++i) coordinates_x[i] *= YAC_RAD;
250 for (size_t i = 0; i <= num_cells[1]; ++i) coordinates_y[i] *= YAC_RAD;
251
255 local_start[my_target_rank], local_count[my_target_rank], with_halo);
256 struct yac_basic_grid * tgt_grid =
258 struct yac_basic_grid * src_grid =
260
261 struct yac_dist_grid_pair * grid_pair =
262 yac_dist_grid_pair_new(tgt_grid, src_grid, MPI_COMM_WORLD);
263
264 struct yac_interp_field src_fields[] =
265 {{.location = YAC_LOC_CORNER, .coordinates_idx = SIZE_MAX, .masks_idx = SIZE_MAX}};
266 size_t num_src_fields = sizeof(src_fields) / sizeof(src_fields[0]);
268 {.location = YAC_LOC_CORNER, .coordinates_idx = SIZE_MAX, .masks_idx = SIZE_MAX};
269
270 struct yac_interp_grid * interp_grid =
273
274 struct interp_method * method_stack[] =
277 NULL};
278
279 struct yac_interp_weights * weights =
280 yac_interp_method_do_search(method_stack, interp_grid);
281
282 yac_interp_method_delete(method_stack);
283
284 //---------------------
285 // do the interpolation
286 //---------------------
287
288 double target_field[MAX_COLLECTION_SIZE][16];
289 double * target_data[MAX_COLLECTION_SIZE];
290 for (unsigned i = 0; i < MAX_COLLECTION_SIZE; ++i)
291 target_data[i] = target_field[i];
292
294 ++collection_size) {
295
296 //---------------------------
297 // set up field interpolation
298 //---------------------------
299
300 struct yac_interpolation * interpolation[2];
301 interpolation[0] =
303 weights, reorder_type, collection_size,
304 YAC_FRAC_MASK_NO_VALUE, 1.0, 0.0, NULL, 1, 1);
305 interpolation[1] = yac_interpolation_copy(interpolation[0]);
306
307 for (size_t k = 0; k < 2; ++k) {
308
309 for (unsigned i = 0; i < collection_size; ++i)
310 for (unsigned j = 0; j < 16; ++j)
311 target_field[i][j] = -1;
312
313 yac_interpolation_execute_get(interpolation[k], target_data);
314
315 //----------------------------
316 // check interpolation results
317 //----------------------------
318
319 double ref_target_data[2][16] = {{3.5,4.5,5.5,6.5,
320 9.5,10.5,11.5,12.5,
321 15.5,16.5,17.5,18.5,
322 21.5,22.5,23.5,24.5},
323 {6.5,7.5,1337,1337,
324 12.5,13.5,1337,1337,
325 18.5,19.5,1337,1337,
326 24.5,25.5,1337,1337}};
327
328 for (unsigned i = 0; i < collection_size; ++i) {
329 for (unsigned j = 0; j < 16; ++j) {
330 if ((double_are_equal(ref_target_data[my_target_rank][j], -1.0)) ||
331 (double_are_equal(ref_target_data[my_target_rank][j], 1337.0))) {
332 if (double_are_unequal(target_data[i][j],
333 ref_target_data[my_target_rank][j]))
334 PUT_ERR("error in interpolated data on target side\n");
335 } else {
337 target_data[i][j],
338 ref_target_data[my_target_rank][j] + (double)(i * 30)))
339 PUT_ERR("error in interpolated data on target side\n");
340 }
341 }
342 }
343 }
344
345 for (int i = 0; i < 2; ++i)
346 yac_interpolation_delete(interpolation[i]);
347 }
348
349 //---------
350 // clean up
351 //---------
352
353 yac_basic_grid_delete(src_grid);
354 yac_basic_grid_delete(tgt_grid);
356 yac_interp_grid_delete(interp_grid);
357 yac_dist_grid_pair_delete(grid_pair);
358}
359
struct yac_basic_grid * yac_basic_grid_new(char const *name, struct yac_basic_grid_data grid_data)
Definition basic_grid.c:53
struct yac_basic_grid * yac_basic_grid_empty_new(char const *name)
Definition basic_grid.c:66
void yac_basic_grid_delete(struct yac_basic_grid *grid)
Definition basic_grid.c:73
void yac_dist_grid_pair_delete(struct yac_dist_grid_pair *grid_pair)
Definition dist_grid.c:2315
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:2063
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)
Definition interp_grid.c:30
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_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
struct yac_interpolation * yac_interpolation_copy(struct yac_interpolation *interp)
Create a deep copy of an interpolation object.
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
@ YAC_LOC_CORNER
Definition location.h:15
#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:25
size_t num_src_fields
Definition interp_grid.c:26
struct yac_dist_grid_pair * grid_pair
Definition interp_grid.c:24
struct yac_interp_field src_fields[]
Definition interp_grid.c:27
int double_are_unequal(double a, double b)
int double_are_equal(double a, double b)
int collection_size
static MPI_Comm split_comm
char const grid_name_tgt[]
#define MAX_COLLECTION_SIZE
char const grid_name_src[]
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