YAC 3.13.0
Yet Another Coupler
Loading...
Searching...
No Matches
test_interpolation_parallel6.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 <limits.h>
7#include <math.h>
8
9#include <mpi.h>
10#include <yaxt.h>
11#include "tests.h"
12#include "test_common.h"
15
21#define UNSET (1337.0)
22#define FRAC_MASK_VALUE (-1337.0)
23
24#define SCALE(X) \
25 ((X) * scaling_factors[scaling_factor_idx] + \
26 scaling_summand[scaling_summand_idx])
27
28static int utest_check_interpolation(
29 struct yac_interpolation * interp,
30 double * src_data_raw, size_t num_src_fields, size_t src_field_size,
31 double * ref_tgt_data_raw, size_t tgt_field_size, size_t src_collection_size,
32 struct yac_collection_selection * collection_selection);
33static int utest_check_interpolation_frac(
34 struct yac_interpolation * interp,
35 double * src_data_raw, double * src_frac_mask_raw,
36 size_t num_src_fields, size_t src_field_size,
37 double * ref_tgt_data_raw, size_t tgt_field_size, size_t src_collection_size,
38 struct yac_collection_selection * collection_selection);
39
40int main(void) {
41
42 MPI_Init(NULL, NULL);
43 xt_initialize(MPI_COMM_WORLD);
44
45 int comm_rank, comm_size;
46 MPI_Comm_rank(MPI_COMM_WORLD, &comm_rank);
47 MPI_Comm_size(MPI_COMM_WORLD, &comm_size);
48 MPI_Barrier(MPI_COMM_WORLD);
49
50 double scaling_factors[] = {1.0, 0.1, -1.0, 0.5, 10.0};
51 double scaling_summand[] = {0.0, -1.0, 1.0, 10.0};
52 enum {
53 num_procs = 4,
54 NUM_SCALING_FACTORS = sizeof(scaling_factors) / sizeof(scaling_factors[0]),
55 NUM_SCALING_SUMMAND = sizeof(scaling_summand) / sizeof(scaling_summand[0]),
56 };
57
58 if (comm_size != num_procs) {
59 PUT_ERR("ERROR: wrong number of processes");
60 xt_finalize();
61 MPI_Finalize();
62 return TEST_EXIT_CODE;
63 }
64
65 { // testing yac_interpolation_add_fixed
66 // 16 target points, all odd positions receive -1 and the others 1
67
68 enum {num_fixed_values = 2};
69 double value[num_fixed_values] = {-1.0, 1.0};
70 size_t count[num_fixed_values] = {8, 8};
71 size_t pos[num_fixed_values][8] = {{0, 2, 4, 6, 8, 10, 12, 14},
72 {1, 3, 5, 7, 9, 11, 13, 15}};
73
74 for (int scaling_factor_idx = 0;
75 scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
76 for (int scaling_summand_idx = 0;
77 scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
78
79 enum {collection_size = 1};
82 struct yac_interpolation * interp =
85 scaling_factors[scaling_factor_idx],
86 scaling_summand[scaling_summand_idx]);
87 for (int i = 0; i < num_fixed_values; ++i)
89 interp, value[i], count[i], pos[i]);
90
91 enum {num_src_fields = 1};
92 enum {src_field_size = 0};
93 double * src_data_raw = NULL;
94 enum {tgt_field_size = 16};
95 double ref_tgt_data_raw[collection_size][tgt_field_size] =
96 {{-1,1,-1,1,-1,1,-1,1,-1,1,-1,1,-1,1,-1,1}};
97
98 if (utest_check_interpolation(
99 interp, src_data_raw, num_src_fields, src_field_size,
100 &ref_tgt_data_raw[0][0], tgt_field_size,
102 PUT_ERR("ERROR in yac_interpolation_add_fixed");
103
106 }
107 }
108 }
109
110 { // testing yac_interpolation_add_direct
111 // rank 0 pos {0,2} -> rank 1 pos {0,1}
112 // rank 0 pos {0,2} -> rank 2 pos {0,2}
113 // rank 1 pos {0,2} -> rank 2 pos {1,3}
114 // (rank 0 is only source
115 // rank 1 is source and target
116 // rank 2 is only target
117 // rank 3 is neither source nor target)
118 Xt_int src_indices[num_procs][4] =
119 {{0,1,2,3},{4,5,6,7},{8,9,10,11},{12,13,14,15}};
120 size_t num_src_indices[num_procs] = {4,4,0,0};
121 Xt_int dst_indices[num_procs][4] = {{-1},{0,2},{0,4,2,6},{-1}};
122 size_t num_dst_indices[num_procs] = {0,2,4,0};
123 Xt_idxlist src_idxlist =
124 xt_idxvec_new(src_indices[comm_rank], num_src_indices[comm_rank]);
125 Xt_idxlist dst_idxlist =
126 xt_idxvec_new(dst_indices[comm_rank], num_dst_indices[comm_rank]);
127 Xt_xmap xmap
128 = xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
129 Xt_redist redist = xt_redist_p2p_new(xmap, MPI_DOUBLE);
130 xt_xmap_delete(xmap);
131 xt_idxlist_delete(dst_idxlist);
132 xt_idxlist_delete(src_idxlist);
133
134 // collection selection configurations
135 struct {
136 size_t N;
137 size_t * indices;
138 } collection_selection_configs[] =
139 {{.N = 3, .indices = (size_t[]){0, 2, 3}},
140 {.N = 4, .indices = (size_t[]){0, 1, 2, 3}},
141 {.N = 2, .indices = NULL}};
142 enum {
143 NUM_COLL_SEL_CONFIGS =
144 sizeof(collection_selection_configs)/
145 sizeof(collection_selection_configs[0]) };
146
147 // loop over collection selections
148 for (size_t sel_idx = 0; sel_idx < NUM_COLL_SEL_CONFIGS; ++sel_idx) {
149
152 collection_selection_configs[sel_idx].N,
153 collection_selection_configs[sel_idx].indices);
154
155 for (int scaling_factor_idx = 0;
156 scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
157 for (int scaling_summand_idx = 0;
158 scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
159
160 struct yac_interpolation * interp =
163 scaling_factors[scaling_factor_idx],
164 scaling_summand[scaling_summand_idx]);
165 yac_interpolation_add_direct(interp, redist);
166
167 enum {NUM_SRC_FIELDS = 1, SRC_COLLECTION_SIZE = 5};
168 size_t const src_field_size[num_procs] = {4,4,0,0};
169 double src_data_raw[num_procs][NUM_SRC_FIELDS][4] =
170 {{{0,1,2,3}},{{4,5,6,7}},{{-1}},{{-1}}};
171 size_t const tgt_field_size[num_procs] = {0,4,4,0};
172 double ref_tgt_data_raw[num_procs][4] =
173 {{UNSET},
174 {SCALE(0),SCALE(2),UNSET,UNSET},
175 {SCALE(0),SCALE(4),SCALE(2),SCALE(6)},
176 {UNSET}};
177
178 if (utest_check_interpolation(
179 interp,
180 (src_field_size[comm_rank] > 0)?
181 (&src_data_raw[comm_rank][0][0]):NULL,
182 NUM_SRC_FIELDS, src_field_size[comm_rank],
183 (tgt_field_size[comm_rank] > 0)?
184 (&ref_tgt_data_raw[comm_rank][0]):NULL,
185 tgt_field_size[comm_rank], SRC_COLLECTION_SIZE,
187 PUT_ERR("ERROR in yac_interpolation_add_direct");
188
190 } // NUM_SCALING_SUMMAND
191 } // NUM_SCALING_FACTORS
193 } // NUM_COLL_SEL_CONFIGS
194 xt_redist_delete(redist);
195 }
196
197 { // testing yac_interpolation_add_direct with dynamic fractional masking
198 // rank 0 pos {0,2} -> rank 1 pos {0,1} (pos 2 is masked)
199 // rank 0 pos {0,2} -> rank 2 pos {0,2} (pos 2 is masked)
200 // rank 1 pos {0,2} -> rank 2 pos {1,3}
201 // (rank 0 is only source
202 // rank 1 is source and target
203 // rank 2 is only target
204 // rank 3 is neither source nor target)
205 Xt_int src_indices[num_procs][4] =
206 {{0,1,2,3},{4,5,6,7},{8,9,10,11},{12,13,14,15}};
207 size_t num_src_indices[num_procs] = {4,4,0,0};
208 Xt_int dst_indices[num_procs][4] = {{-1},{0,2},{0,4,2,6},{-1}};
209 size_t num_dst_indices[num_procs] = {0,2,4,0};
210 Xt_idxlist src_idxlist =
211 xt_idxvec_new(src_indices[comm_rank], num_src_indices[comm_rank]);
212 Xt_idxlist dst_idxlist =
213 xt_idxvec_new(dst_indices[comm_rank], num_dst_indices[comm_rank]);
214 Xt_xmap xmap
215 = xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
216 Xt_redist redist = xt_redist_p2p_new(xmap, MPI_DOUBLE);
217 xt_xmap_delete(xmap);
218 xt_idxlist_delete(dst_idxlist);
219 xt_idxlist_delete(src_idxlist);
220
221 for (int scaling_factor_idx = 0;
222 scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
223 for (int scaling_summand_idx = 0;
224 scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
225
226 enum {collection_size = 1};
229 struct yac_interpolation * interp =
232 scaling_factors[scaling_factor_idx],
233 scaling_summand[scaling_summand_idx]);
234 yac_interpolation_add_direct(interp, redist);
235
236 enum {num_src_fields = 1};
237 size_t const src_field_size[num_procs] = {4,4,0,0};
238 double src_data_raw[num_procs][num_src_fields][4] =
239 {{{0,1,2,3}},{{4,5,6,7}},{{-1}},{{-1}}};
240 double src_frac_mask_raw[num_procs][num_src_fields][4] =
241 {{{1.0,1.0,0.0,1.0}},{{0.5,0.5,0.5,0.5}},{{-1}},{{-1}}};
242 size_t const tgt_field_size[num_procs] = {0,4,4,0};
243 double ref_tgt_data_raw[num_procs][collection_size][4] =
244 {{{UNSET}},
246 {{SCALE(0),SCALE(4),FRAC_MASK_VALUE,SCALE(6)}},
247 {{UNSET}}};
248
249 if (utest_check_interpolation_frac(
250 interp,
251 (src_field_size[comm_rank] > 0)?
252 (&src_data_raw[comm_rank][0][0]):NULL,
253 (src_field_size[comm_rank] > 0)?
254 (&src_frac_mask_raw[comm_rank][0][0]):NULL,
255 num_src_fields, src_field_size[comm_rank],
256 (tgt_field_size[comm_rank] > 0)?
257 (&ref_tgt_data_raw[comm_rank][0][0]):NULL,
258 tgt_field_size[comm_rank], collection_size,
260 PUT_ERR("ERROR in yac_interpolation_add_direct");
261
264 }
265 }
266 xt_redist_delete(redist);
267 }
268
269 { // testing yac_interpolation_add_direct_mf
270 // rank 0 field 0 id {0,2} field 1 id {1,3} -> rank 1 pos {0,2,3,1}
271 // rank 0 field 0 id {1} field 1 id {2} -> rank 2 pos {0,2}
272 // rank 1 field 0 id {4} field 1 id {5} -> rank 2 pos {1,3}
273 // (rank 0 is only source
274 // rank 1 is source and target
275 // rank 2 is only target
276 // rank 3 is neither source nor target)
277 enum {num_src_fields = 2};
278 enum {num_src_indices = 4};
279 Xt_int src_indices[num_procs][num_src_indices] =
280 {{0,1,2,3},{4,5,6,7},{8,9,10,11},{12,13,14,15}};
281 Xt_int dst_indices[num_procs][2][2] = {{{-1},{-1}},
282 {{0,2},{1,3}},
283 {{1,4},{2,5}},
284 {{-1},{-1}}};
285 int src_offsets[num_procs] = {0,1,2,3};
286 int dst_offsets[num_procs][num_src_fields][2] = {{{-1},{-1}},
287 {{0,2},{3,1}},
288 {{0,1},{2,3}},
289 {{-1},{-1}}};
290 size_t num_dst_indices[num_procs][num_src_fields] =
291 {{0,0},{2,2},{2,2},{0,0}};
292 Xt_redist redists[num_src_fields];
293 for (int i = 0; i < num_src_fields; ++i) {
294 Xt_idxlist src_idxlist =
295 xt_idxvec_new(src_indices[comm_rank], num_src_indices);
296 Xt_idxlist dst_idxlist =
297 xt_idxvec_new(dst_indices[comm_rank][i], num_dst_indices[comm_rank][i]);
298 Xt_xmap xmap =
299 xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
300 redists[i] =
301 xt_redist_p2p_off_new(
302 xmap, src_offsets, dst_offsets[comm_rank][i], MPI_DOUBLE);
303 xt_xmap_delete(xmap);
304 xt_idxlist_delete(dst_idxlist);
305 xt_idxlist_delete(src_idxlist);
306 }
307
308 for (int scaling_factor_idx = 0;
309 scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
310 for (int scaling_summand_idx = 0;
311 scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
312
313 enum {collection_size = 1};
316 struct yac_interpolation * interp =
319 scaling_factors[scaling_factor_idx],
320 scaling_summand[scaling_summand_idx]);
321 yac_interpolation_add_direct_mf(interp, redists, num_src_fields);
322
323 size_t const src_field_size[num_procs] = {4,4,4,4};
324 double src_data_raw[num_procs][num_src_fields][num_src_indices] =
325 {{{0,1,2,3},{00,10,20,30}},
326 {{4,5,6,7},{40,50,60,70}},
327 {{8,9,10,11},{80,90,100,110}},
328 {{12,13,14,15},{120,130,140,150}}};
329 size_t const tgt_field_size[num_procs] = {0,4,4,0};
330 double ref_tgt_data_raw[num_procs][collection_size][4] =
331 {{{UNSET}},
332 {{SCALE(0),SCALE(30),SCALE(2),SCALE(10)}},
333 {{SCALE(1),SCALE(4),SCALE(20),SCALE(50)}},
334 {{UNSET}}};
335
336 if (utest_check_interpolation(
337 interp,
338 (src_field_size[comm_rank] > 0)?
339 (&src_data_raw[comm_rank][0][0]):NULL,
340 num_src_fields, src_field_size[comm_rank],
341 (tgt_field_size[comm_rank] > 0)?
342 (&ref_tgt_data_raw[comm_rank][0][0]):NULL,
343 tgt_field_size[comm_rank], collection_size,
345 PUT_ERR("ERROR in yac_interpolation_add_direct_mf");
346
349 }
350 }
351 xt_redist_delete(redists[1]);
352 xt_redist_delete(redists[0]);
353 }
354
355 { // testing yac_interpolation_add_direct_mf with dynamic fractional mask
356 // rank 0 field 0 id {0,2} field 1 id {1,3} -> rank 1 pos {0,2,3,1}
357 // rank 0 field 0 id {1} field 1 id {2} -> rank 2 pos {0,2}
358 // rank 1 field 0 id {4} field 1 id {5} -> rank 2 pos {1,3}
359 // (rank 0 is only source
360 // rank 1 is source and target
361 // rank 2 is only target
362 // rank 3 is neither source nor target)
363 enum {num_src_fields = 2};
364 enum {num_src_indices = 4};
365 Xt_int src_indices[num_procs][num_src_indices] =
366 {{0,1,2,3},{4,5,6,7},{8,9,10,11},{12,13,14,15}};
367 Xt_int dst_indices[num_procs][2][2] = {{{-1},{-1}},
368 {{0,2},{1,3}},
369 {{1,4},{2,5}},
370 {{-1},{-1}}};
371 int src_offsets[num_procs] = {0,1,2,3};
372 int dst_offsets[num_procs][num_src_fields][2] = {{{-1},{-1}},
373 {{0,2},{3,1}},
374 {{0,1},{2,3}},
375 {{-1},{-1}}};
376 size_t num_dst_indices[num_procs][num_src_fields] =
377 {{0,0},{2,2},{2,2},{0,0}};
378 Xt_redist redists[num_src_fields];
379 for (int i = 0; i < num_src_fields; ++i) {
380 Xt_idxlist src_idxlist =
381 xt_idxvec_new(src_indices[comm_rank], num_src_indices);
382 Xt_idxlist dst_idxlist =
383 xt_idxvec_new(dst_indices[comm_rank][i], num_dst_indices[comm_rank][i]);
384 Xt_xmap xmap =
385 xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
386 redists[i] =
387 xt_redist_p2p_off_new(
388 xmap, src_offsets, dst_offsets[comm_rank][i], MPI_DOUBLE);
389 xt_xmap_delete(xmap);
390 xt_idxlist_delete(dst_idxlist);
391 xt_idxlist_delete(src_idxlist);
392 }
393
394 for (int scaling_factor_idx = 0;
395 scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
396 for (int scaling_summand_idx = 0;
397 scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
398
399 enum {collection_size = 1};
402 struct yac_interpolation * interp =
405 scaling_factors[scaling_factor_idx],
406 scaling_summand[scaling_summand_idx]);
407 yac_interpolation_add_direct_mf(interp, redists, num_src_fields);
408
409 size_t const src_field_size[num_procs] = {4,4,4,4};
410 double src_data_raw[num_procs][num_src_fields][num_src_indices] =
411 {{{0,1,2,3},{00,10,20,30}},
412 {{4,5,6,7},{40,50,60,70}},
413 {{8,9,10,11},{80,90,100,110}},
414 {{12,13,14,15},{120,130,140,150}}};
415 double src_frac_mask_raw[num_procs][num_src_fields][num_src_indices] =
416 {{{1.0,0.0,1.0,1.0},{1.0,1.0,1.0,1.0}},
417 {{0.0,1.0,1.0,1.0},{1.0,0.0,1.0,1.0}},
418 {{1.0,1.0,1.0,1.0},{1.0,1.0,1.0,1.0}},
419 {{1.0,1.0,1.0,1.0},{1.0,1.0,1.0,1.0}}};
420 size_t const tgt_field_size[num_procs] = {0,4,4,0};
421 double ref_tgt_data_raw[num_procs][collection_size][4] =
422 {{{UNSET}},
423 {{SCALE(0),SCALE(30),SCALE(2),SCALE(10)}},
425 {{UNSET}}};
426
427 if (utest_check_interpolation_frac(
428 interp,
429 (src_field_size[comm_rank] > 0)?
430 (&src_data_raw[comm_rank][0][0]):NULL,
431 (src_field_size[comm_rank] > 0)?
432 (&src_frac_mask_raw[comm_rank][0][0]):NULL,
433 num_src_fields, src_field_size[comm_rank],
434 (tgt_field_size[comm_rank] > 0)?
435 (&ref_tgt_data_raw[comm_rank][0][0]):NULL,
436 tgt_field_size[comm_rank], collection_size,
438 PUT_ERR("ERROR in yac_interpolation_add_direct_mf");
439
442 }
443 }
444 xt_redist_delete(redists[1]);
445 xt_redist_delete(redists[0]);
446 }
447
448 { // test yac_interpolation_add_sum_at_src and
449 // yac_interpolation_add_weight_sum_mvp_at_src
450 // distribution of source and destination (only field_id 0) ids
451 // rank | 0 | 1 | 2 | 3
452 // ------------|----------------|----------------|----------------|---------------
453 // field_id 0 | 0, 1, 2, 3 | 4, 5, 6, 7 | 8, 9, 10, 11 | 12, 13, 14, 15
454 // field_id 1 | 16, 17, 18, 19 | 20, 21, 22, 23 | 24, 25, 26, 27 | 28, 29, 30, 31
455 //
456 // tgt id | src ids | generate on rank | optional weights
457 // -------|--------------|------------------|--------------------
458 // 4 | 0, 17, 2, 19 | 0 | 0.1, 0.2, 0.3, 0.4
459 // 5 | 2, 19, 20, 5 | 1 | 0.5, 0.6, 0.7, 0.8
460 // 6 | 6, 22, 7, 23 | 1 | 0.9, 1.0, 1.1, 1.2
461 // 7 | 23 | 1 | 1.3
462 // 8 | 2, 19, 22, 7 | 0 | 1.4, 1.5, 1.6, 1.7
463 // 9 | 0, 17, 2, 19 | 0 | 0.1, 0.2, 0.3, 0.4
464 // 10 | 5, 6 | 1 | 2.2, 2.3
465 // 11 | 21, 22 | 1 | 2.4, 2.5
466
467 // at first the data for all stencils is gathered on a single process
468 // this is done using the halo_redists
469 enum {num_src_fields = 2};
470 Xt_redist halo_redists[num_src_fields];
471 {
472 enum {num_src_indices = 4};
473 Xt_int src_indices[num_procs][num_src_fields][num_src_indices] =
474 {{{0,1,2,3},{16,17,18,19}},
475 {{4,5,6,7},{20,21,22,23}},
476 {{8,9,10,11},{24,25,26,27}},
477 {{12,13,14,15},{28,29,30,31}}};
478 Xt_int dst_indices[num_procs][num_src_fields][1] =
479 {{{7},{22}},
480 {{2},{19}},
481 {{-1},{-1}},
482 {{-1},{-1}}};
483 size_t num_dst_indices[num_procs][num_src_fields] =
484 {{1,1},{1,1},{0,0},{0,0}};
485 for (int field_id = 0; field_id < num_src_fields; ++field_id) {
486 Xt_idxlist src_idxlist =
487 xt_idxvec_new(src_indices[comm_rank][field_id], num_src_indices);
488 Xt_idxlist dst_idxlist =
489 xt_idxvec_new(
490 dst_indices[comm_rank][field_id],
491 num_dst_indices[comm_rank][field_id]);
492 Xt_xmap xmap =
493 xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
494 halo_redists[field_id] = xt_redist_p2p_new(xmap, MPI_DOUBLE);
495 xt_xmap_delete(xmap);
496 xt_idxlist_delete(dst_idxlist);
497 xt_idxlist_delete(src_idxlist);
498 }
499 }
500
501 // number of stencils to be computed locally (tgt id 4 and 9 have the same
502 // stencil -> we compute it only once)
503 size_t tgt_count[num_procs] = {2,5,0,0};
504 // number of source points per stencil
505 size_t num_src_per_tgt[num_procs][5] =
506 {{4,4},{4,4,1,2,2},{(size_t)-1},{(size_t)-1}};
507 double weights[num_procs][13] =
508 {{0.1,0.2,0.3,0.4,
509 1.4,1.5,1.6,1.7},
510 {0.5,0.6,0.7,0.8,
511 0.9,1.0,1.1,1.2,
512 1.3,
513 2.2,2.3,
514 2.4,2.5},{-1},{-1}};
515 // source field index for each source point used in the stencils (SIZE_MAX
516 // indicates that the respective source point is in the halo buffer)
517 size_t src_field_idx[num_procs][13] =
518 {{0,1,0,1,
519 0,1,SIZE_MAX,SIZE_MAX},
520 {SIZE_MAX,SIZE_MAX,1,0,
521 0,1,0,1,
522 1,
523 0,0,
524 1,1},{(size_t)-1},{(size_t)-1}};
525 // source index; position of the respective source point within the source
526 // field or halo buffer
527 size_t src_idx[num_procs][13] =
528 {{0,1,2,3,
529 2,3,1,0},
530 {0,1,0,1,
531 2,2,3,3,
532 3,
533 1,2,
534 1,2},{(size_t)-1},{(size_t)-1}};
535
536 // the results from the stencil computation get send directly to their
537 // final location in the destination buffer on the respective process using
538 // result_redist
539 // (tgt id 4 and 9 have the same stencil)
540 Xt_redist result_redist;
541 {
542 Xt_int src_indices[num_procs][5] = {{4,8},{5,6,7,10,11},{-1},{-1}};
543 size_t num_src_indices[num_procs] = {2,5,0,0};
544 Xt_int dst_indices[num_procs][4] = {{-1},{4,5,6,7},{8,4,10,11},{-1}};
545 size_t num_dst_indices[num_procs] = {0,4,4,0};
546 Xt_idxlist src_idxlist =
547 xt_idxvec_new(src_indices[comm_rank], num_src_indices[comm_rank]);
548 Xt_idxlist dst_idxlist =
549 xt_idxvec_new(dst_indices[comm_rank], num_dst_indices[comm_rank]);
550 Xt_xmap xmap = xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
551 result_redist = xt_redist_p2p_new(xmap, MPI_DOUBLE);
552 xt_xmap_delete(xmap);
553 xt_idxlist_delete(dst_idxlist);
554 xt_idxlist_delete(src_idxlist);
555 }
556
557 for (int scaling_factor_idx = 0;
558 scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
559 for (int scaling_summand_idx = 0;
560 scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
561 for (int with_weights = 0; with_weights < 2; ++with_weights) {
562
563 enum {collection_size = 1};
566 struct yac_interpolation * interp =
569 scaling_factors[scaling_factor_idx],
570 scaling_summand[scaling_summand_idx]);
571 if (with_weights) {
573 interp, halo_redists, tgt_count[comm_rank],
574 num_src_per_tgt[comm_rank], weights[comm_rank],
575 src_field_idx[comm_rank], src_idx[comm_rank], num_src_fields,
576 result_redist);
577 } else {
579 interp, halo_redists, tgt_count[comm_rank],
580 num_src_per_tgt[comm_rank], src_field_idx[comm_rank],
581 src_idx[comm_rank], num_src_fields, result_redist);
582 }
583
584 enum {src_field_size = 4};
585 double src_data_raw[num_procs][2][src_field_size] =
586 {{{0,1,2,3}, {16,17,18,19}},
587 {{4,5,6,7}, {20,21,22,23}},
588 {{8,9,10,11}, {24,25,26,27}},
589 {{12,13,14,15},{28,29,30,31}}};
590 size_t const tgt_field_size[num_procs] = {0,4,4,0};
591 double ref_tgt_data_raw[2][num_procs][collection_size][4] =
592 {{{{UNSET}},
593 {{SCALE(0+17+2+19),SCALE(2+19+20+5),SCALE(6+22+7+23),SCALE(23)}},
594 {{SCALE(2+19+22+7),SCALE(0+17+2+19),SCALE(5+6),SCALE(21+22)}},
595 {{UNSET}}},
596 {{{UNSET}},
597 {{SCALE(0.1*0+0.2*17+0.3*2+0.4*19),
598 SCALE(0.5*2+0.6*19+0.7*20+0.8*5),
599 SCALE(0.9*6+1.0*22+1.1*7+1.2*23),
600 SCALE(1.3*23)}},
601 {{SCALE(1.4*2+1.5*19+1.6*22+1.7*7),
602 SCALE(0.1*0+0.2*17+0.3*2+0.4*19),
603 SCALE(2.2*5+2.3*6),
604 SCALE(2.4*21+2.5*22)}},
605 {{UNSET}}}};
606
607 if (utest_check_interpolation(
608 interp, &src_data_raw[comm_rank][0][0],
609 num_src_fields, src_field_size,
610 &ref_tgt_data_raw[with_weights][comm_rank][0][0],
611 tgt_field_size[comm_rank], collection_size,
613 PUT_ERR("ERROR in yac_interpolation_add_sum_at_src/"
614 "yac_interpolation_add_weight_sum_mvp_at_src");
615
618 }
619 }
620 }
621 xt_redist_delete(result_redist);
622 xt_redist_delete(halo_redists[1]);
623 xt_redist_delete(halo_redists[0]);
624 }
625
626 { // test yac_interpolation_add_sum_at_src and
627 // yac_interpolation_add_weight_sum_mvp_at_src with dynamic fractional masking
628 // distribution of source and destination (only field_id 0) ids
629 // rank | 0 | 1 | 2 | 3
630 // ------------|--------------------|--------------------|--------------------|-------------------
631 // field_id 0 | 0, 1, 2, 3 | 4, 5, 6, 7 | 8, 9, 10, 11 | 12, 13, 14, 15
632 // frac mask | 1.0, 1.0, 1.0, 1.0 | 1.0, 0.5, 0.0, 0.0 | 1.0, 1.0, 1.0, 1.0 | 1.0, 1.0, 1.0, 1.0
633 // field_id 1 | 16, 17, 18, 19 | 20, 21, 22, 23 | 24, 25, 26, 27 | 28, 29, 30, 31
634 // frac_mask | 1.0, 1.0, 1.0, 1.0 | 0.5, 0.0, 0.0, 0.0 | 1.0, 1.0, 1.0, 1.0 | 1.0, 1.0, 1.0, 1.0
635 //
636 // tgt id | src ids | generate on rank | optional weights
637 // -------|--------------|------------------|--------------------
638 // 4 | 0, 17, 2, 19 | 0 | 0.1, 0.2, 0.3, 0.4
639 // 5 | 2, 19, 20, 5 | 1 | 0.5, 0.6, 0.7, 0.8
640 // 6 | 6, 22, 7, 23 | 1 | 0.9, 1.0, 1.1, 1.2
641 // 7 | 23 | 1 | 1.3
642 // 8 | 2, 19, 22, 7 | 0 | 1.4, 1.5, 1.6, 1.7
643 // 9 | 0, 17, 2, 19 | 0 | 0.1, 0.2, 0.3, 0.4
644 // 10 | 5, 6 | 1 | 2.2, 2.3
645 // 11 | 21, 22 | 1 | 2.4, 2.5
646
647 // at first the data for all stencils is gathered on a single process
648 // this is done using the halo_redists
649 enum {num_src_fields = 2};
650 Xt_redist halo_redists[num_src_fields];
651 {
652 enum {num_src_indices = 4};
653 Xt_int src_indices[num_procs][num_src_fields][num_src_indices] =
654 {{{0,1,2,3},{16,17,18,19}},
655 {{4,5,6,7},{20,21,22,23}},
656 {{8,9,10,11},{24,25,26,27}},
657 {{12,13,14,15},{28,29,30,31}}};
658 Xt_int dst_indices[num_procs][num_src_fields][1] =
659 {{{7},{22}},
660 {{2},{19}},
661 {{-1},{-1}},
662 {{-1},{-1}}};
663 size_t num_dst_indices[num_procs][num_src_fields] =
664 {{1,1},{1,1},{0,0},{0,0}};
665 for (int field_id = 0; field_id < num_src_fields; ++field_id) {
666 Xt_idxlist src_idxlist =
667 xt_idxvec_new(src_indices[comm_rank][field_id], num_src_indices);
668 Xt_idxlist dst_idxlist =
669 xt_idxvec_new(
670 dst_indices[comm_rank][field_id],
671 num_dst_indices[comm_rank][field_id]);
672 Xt_xmap xmap =
673 xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
674 halo_redists[field_id] = xt_redist_p2p_new(xmap, MPI_DOUBLE);
675 xt_xmap_delete(xmap);
676 xt_idxlist_delete(dst_idxlist);
677 xt_idxlist_delete(src_idxlist);
678 }
679 }
680
681 // number of stencils to be computed locally (tgt id 4 and 9 have the same
682 // stencil -> we compute it only once)
683 size_t tgt_count[num_procs] = {2,5,0,0};
684 // number of source points per stencil
685 size_t num_src_per_tgt[num_procs][5] =
686 {{4,4},{4,4,1,2,2},{(size_t)-1},{(size_t)-1}};
687 double weights[num_procs][13] =
688 {{0.1,0.2,0.3,0.4,
689 1.4,1.5,1.6,1.7},
690 {0.5,0.6,0.7,0.8,
691 0.9,1.0,1.1,1.2,
692 1.3,
693 2.2,2.3,
694 2.4,2.5},{-1},{-1}};
695 // source field index for each source point used in the stencils (SIZE_MAX
696 // indicates that the respective source point is in the halo buffer)
697 size_t src_field_idx[num_procs][13] =
698 {{0,1,0,1,
699 0,1,SIZE_MAX,SIZE_MAX},
700 {SIZE_MAX,SIZE_MAX,1,0,
701 0,1,0,1,
702 1,
703 0,0,
704 1,1},{(size_t)-1},{(size_t)-1}};
705 // source index; position of the respective source point within the source
706 // field or halo buffer
707 size_t src_idx[num_procs][13] =
708 {{0,1,2,3,
709 2,3,1,0},
710 {0,1,0,1,
711 2,2,3,3,
712 3,
713 1,2,
714 1,2},{(size_t)-1},{(size_t)-1}};
715
716 // the results from the stencil computation get send directly to their
717 // final location in the destination buffer on the respective process using
718 // result_redist
719 // (tgt id 4 and 9 have the same stencil)
720 Xt_redist result_redist;
721 {
722 Xt_int src_indices[num_procs][5] = {{4,8},{5,6,7,10,11},{-1},{-1}};
723 size_t num_src_indices[num_procs] = {2,5,0,0};
724 Xt_int dst_indices[num_procs][4] = {{-1},{4,5,6,7},{8,4,10,11},{-1}};
725 size_t num_dst_indices[num_procs] = {0,4,4,0};
726 Xt_idxlist src_idxlist =
727 xt_idxvec_new(src_indices[comm_rank], num_src_indices[comm_rank]);
728 Xt_idxlist dst_idxlist =
729 xt_idxvec_new(dst_indices[comm_rank], num_dst_indices[comm_rank]);
730 Xt_xmap xmap = xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
731 result_redist = xt_redist_p2p_new(xmap, MPI_DOUBLE);
732 xt_xmap_delete(xmap);
733 xt_idxlist_delete(dst_idxlist);
734 xt_idxlist_delete(src_idxlist);
735 }
736
737 for (int scaling_factor_idx = 0;
738 scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
739 for (int scaling_summand_idx = 0;
740 scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
741 for (int with_weights = 0; with_weights < 2; ++with_weights) {
742
743 enum {collection_size = 1};
746 struct yac_interpolation * interp =
749 scaling_factors[scaling_factor_idx],
750 scaling_summand[scaling_summand_idx]);
751 if (with_weights) {
753 interp, halo_redists, tgt_count[comm_rank],
754 num_src_per_tgt[comm_rank], weights[comm_rank],
755 src_field_idx[comm_rank], src_idx[comm_rank], num_src_fields,
756 result_redist);
757 } else {
759 interp, halo_redists, tgt_count[comm_rank],
760 num_src_per_tgt[comm_rank], src_field_idx[comm_rank],
761 src_idx[comm_rank], num_src_fields, result_redist);
762 }
763
764 enum {src_field_size = 4};
765 double src_data_raw[num_procs][2][src_field_size] =
766 {{{0,1,2,3}, {16,17,18,19}},
767 {{4,5,6,7}, {20,21,22,23}},
768 {{8,9,10,11}, {24,25,26,27}},
769 {{12,13,14,15},{28,29,30,31}}};
770 double src_frac_mask_raw[num_procs][2][src_field_size] =
771 {{{1.0, 1.0, 1.0, 1.0}, {1.0, 1.0, 1.0, 1.0}},
772 {{1.0, 0.5, 0.0, 0.0}, {0.5, 0.0, 0.0, 0.0}},
773 {{1.0, 1.0, 1.0, 1.0}, {1.0, 1.0, 1.0, 1.0}},
774 {{1.0, 1.0, 1.0, 1.0}, {1.0, 1.0, 1.0, 1.0}}};
775 size_t const tgt_field_size[num_procs] = {0,4,4,0};
776 double ref_tgt_data_raw[2][num_procs][collection_size][4] =
777 {{{{UNSET}},
778 {{SCALE((1.0*0+1.0*17+1.0*2+1.0*19)/(1.0+1.0+1.0+1.0)),
779 SCALE((1.0*2+1.0*19+0.5*20+0.5*5)/(1.0+1.0+0.5+0.5)),
781 {{SCALE((1.0*2+1.0*19+0.0*22+0.0*7)/(1.0+1.0+0.0+0.0)),
782 SCALE((1.0*0+1.0*17+1.0*2+1.0*19)/(1.0+1.0+1.0+1.0)),
783 SCALE((0.5*5+0.0*6)/(0.5+0.0)),
785 {{UNSET}}},
786 {{{UNSET}},
787 {{SCALE((1.0*0.1*0+1.0*0.2*17+1.0*0.3*2+1.0*0.4*19)/
788 (1.0*0.1+1.0*0.2+1.0*0.3+1.0*0.4)),
789 SCALE((1.0*0.5*2+1.0*0.6*19+0.5*0.7*20+0.5*0.8*5)/
790 (1.0*0.5+1.0*0.6+0.5*0.7+0.5*0.8)),
792 {{SCALE((1.0*1.4*2+1.0*1.5*19+0.0*1.6*22+0.0*1.7*7)/
793 (1.0*1.4+1.0*1.5+0.0*1.6+0.0*1.7)),
794 SCALE((1.0*0.1*0+1.0*0.2*17+1.0*0.3*2+1.0*0.4*19)/
795 (1.0*0.1+1.0*0.2+1.0*0.3+1.0*0.4)),
796 SCALE((0.5*2.2*5+0.0*2.3*6)/(0.5*2.2+0.0*2.3)),
798 {{UNSET}}}};
799
800 if (utest_check_interpolation_frac(
801 interp, &src_data_raw[comm_rank][0][0],
802 &src_frac_mask_raw[comm_rank][0][0],
803 num_src_fields, src_field_size,
804 &ref_tgt_data_raw[with_weights][comm_rank][0][0],
805 tgt_field_size[comm_rank], collection_size,
807 PUT_ERR("ERROR in yac_interpolation_add_sum_at_src/"
808 "yac_interpolation_add_weight_sum_mvp_at_src frac");
809
812 }
813 }
814 }
815 xt_redist_delete(result_redist);
816 xt_redist_delete(halo_redists[1]);
817 xt_redist_delete(halo_redists[0]);
818 }
819
820 { // test yac_interpolation_add_sum_at_tgt and
821 // yac_interpolation_add_weight_sum_mvp_at_tgt
822 // distribution of source and destination (only field_id 0) ids
823 // rank | 0 | 1 | 2 | 3
824 // ------------|----------------|----------------|----------------|---------------
825 // field_id 0 | 0, 1, 2, 3 | 4, 5, 6, 7 | 8, 9, 10, 11 | 12, 13, 14, 15
826 // field_id 1 | 16, 17, 18, 19 | 20, 21, 22, 23 | 24, 25, 26, 27 | 28, 29, 30, 31
827 //
828 // tgt id | src ids | optional weights
829 // -------|--------------|--------------------
830 // 4 | 0, 17, 2, 19 | 0.1, 0.2, 0.3, 0.4
831 // 5 | 2, 19, 20, 5 | 0.5, 0.6, 0.7, 0.8
832 // 6 | 6, 22, 7, 23 | 0.9, 1.0, 1.1, 1.2
833 // 7 | 23 | 1.3
834 // 8 | 2, 19, 22, 7 | 1.4, 1.5, 1.6, 1.7
835 // 9 | 0, 17, 2, 19 | 0.1, 0.2, 0.3, 0.4
836 // 10 | 5, 6 | 2.2, 2.3
837 // 11 | 21, 22 | 2.4, 2.5
838
839 // at first the source points required for each target point have to be
840 // transfered to the respective target processes using src_redists
841 enum {num_src_fields = 2};
842 Xt_redist src_redists[num_src_fields];
843 {
844 enum {num_src_indices = 4};
845 Xt_int src_indices[num_procs][num_src_fields][num_src_indices] =
846 {{{0,1,2,3},{16,17,18,19}},
847 {{4,5,6,7},{20,21,22,23}},
848 {{8,9,10,11},{24,25,26,27}},
849 {{12,13,14,15},{28,29,30,31}}};
850
851 Xt_int dst_indices[num_procs][num_src_fields][5] =
852 {{{-1},{-1}},
853 {{0,2},{17,19}},
854 {{0,2,5,6,7},{17,19,21,22}},
855 {{-1},{-1}}};
856 size_t num_dst_indices[num_procs][num_src_fields] =
857 {{0,0},{2,2},{5,4},{0,0}};
858 for (int field_id = 0; field_id < num_src_fields; ++field_id) {
859 Xt_idxlist src_idxlist =
860 xt_idxvec_new(src_indices[comm_rank][field_id], num_src_indices);
861 Xt_idxlist dst_idxlist =
862 xt_idxvec_new(
863 dst_indices[comm_rank][field_id],
864 num_dst_indices[comm_rank][field_id]);
865 Xt_xmap xmap =
866 xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
867 src_redists[field_id] = xt_redist_p2p_new(xmap, MPI_DOUBLE);
868 xt_xmap_delete(xmap);
869 xt_idxlist_delete(dst_idxlist);
870 xt_idxlist_delete(src_idxlist);
871 }
872 }
873
874 // positions of target points to be written
875 size_t tgt_pos[num_procs][4] =
876 {{(size_t)-1},{0,1,2,3},{0,1,2,3},{(size_t)-1}};
877 // number of entries in tgt_pos
878 size_t tgt_count[num_procs] = {0,4,4,0};
879 // number of source points per stencil
880 size_t num_src_per_tgt[num_procs][4] =
881 {{(size_t)-1},{4,4,4,1},{4,4,2,2},{(size_t)-1}};
882 double weights[num_procs][13] =
883 {{-1},
884 {0.1,0.2,0.3,0.4,
885 0.5,0.6,0.7,0.8,
886 0.9,1.0,1.1,1.2,
887 1.3},
888 {1.4,1.5,1.6,1.7,
889 0.1,0.2,0.3,0.4,
890 2.2,2.3,
891 2.4,2.5},
892 {-1}};
893 // source field index for each source point used in the stencils (SIZE_MAX
894 // indicates that the respective source point is in the halo buffer)
895 size_t src_field_idx[num_procs][13] =
896 {{(size_t)-1},
897 {SIZE_MAX,SIZE_MAX,SIZE_MAX,SIZE_MAX,
898 SIZE_MAX,SIZE_MAX,1,0,
899 0,1,0,1,
900 1},
901 {SIZE_MAX,SIZE_MAX,SIZE_MAX,SIZE_MAX,
902 SIZE_MAX,SIZE_MAX,SIZE_MAX,SIZE_MAX,
903 SIZE_MAX,SIZE_MAX,
904 SIZE_MAX,SIZE_MAX},
905 {(size_t)-1}};
906 // source index; position of the respective source point within the source
907 // field or halo buffer
908 size_t src_idx[num_procs][13] =
909 {{(size_t)-1},
910 {0,2,1,3, 1,3,0,1, 2,2,3,3, 3},
911 {1,6,8,4, 0,5,1,6, 2,3, 7,8},
912 {(size_t)-1}};
913
914 for (int scaling_factor_idx = 0;
915 scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
916 for (int scaling_summand_idx = 0;
917 scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
918 for (int with_weights = 0; with_weights < 2; ++with_weights) {
919
920 enum {collection_size = 1};
923 struct yac_interpolation * interp =
926 scaling_factors[scaling_factor_idx],
927 scaling_summand[scaling_summand_idx]);
928 if (with_weights) {
930 interp, src_redists, tgt_pos[comm_rank], tgt_count[comm_rank],
931 num_src_per_tgt[comm_rank], weights[comm_rank],
932 src_field_idx[comm_rank], src_idx[comm_rank], num_src_fields);
933 } else {
935 interp, src_redists, tgt_pos[comm_rank], tgt_count[comm_rank],
936 num_src_per_tgt[comm_rank], src_field_idx[comm_rank],
937 src_idx[comm_rank], num_src_fields);
938 }
939
940 enum {src_field_size = 4};
941 double src_data_raw[num_procs][num_src_fields][src_field_size] =
942 {{{0,1,2,3}, {16,17,18,19}},
943 {{4,5,6,7}, {20,21,22,23}},
944 {{8,9,10,11}, {24,25,26,27}},
945 {{12,13,14,15},{28,29,30,31}}};
946 size_t const tgt_field_size[num_procs] = {0,4,4,0};
947 double ref_tgt_data_raw[2][num_procs][collection_size][4] =
948 {{{{UNSET}},
949 {{SCALE(0+17+2+19),SCALE(2+19+20+5),SCALE(6+22+7+23),SCALE(23)}},
950 {{SCALE(2+19+22+7),SCALE(0+17+2+19),SCALE(5+6),SCALE(21+22)}},
951 {{UNSET}}},
952 {{{UNSET}},
953 {{SCALE(0.1*0+0.2*17+0.3*2+0.4*19),
954 SCALE(0.5*2+0.6*19+0.7*20+0.8*5),
955 SCALE(0.9*6+1.0*22+1.1*7+1.2*23),
956 SCALE(1.3*23)}},
957 {{SCALE(1.4*2+1.5*19+1.6*22+1.7*7),
958 SCALE(0.1*0+0.2*17+0.3*2+0.4*19),
959 SCALE(2.2*5+2.3*6),
960 SCALE(2.4*21+2.5*22)}},
961 {{UNSET}}}};
962
963 if (utest_check_interpolation(
964 interp, &src_data_raw[comm_rank][0][0],
965 num_src_fields, src_field_size,
966 &ref_tgt_data_raw[with_weights][comm_rank][0][0],
967 tgt_field_size[comm_rank], collection_size,
969 PUT_ERR("ERROR in yac_interpolation_add_sum_at_tgt/"
970 "yac_interpolation_add_weight_sum_mvp_at_tgt");
971
974 }
975 }
976 }
977 xt_redist_delete(src_redists[1]);
978 xt_redist_delete(src_redists[0]);
979 }
980
981 { // test yac_interpolation_add_sum_at_tgt and
982 // yac_interpolation_add_weight_sum_mvp_at_tgt with dynamic fractional masking
983 // distribution of source and destination (only field_id 0) ids
984 // rank | 0 | 1 | 2 | 3
985 // ------------|--------------------|--------------------|--------------------|-------------------
986 // field_id 0 | 0, 1, 2, 3 | 4, 5, 6, 7 | 8, 9, 10, 11 | 12, 13, 14, 15
987 // frac mask | 1.0, 1.0, 1.0, 1.0 | 1.0, 0.5, 0.0, 0.0 | 1.0, 1.0, 1.0, 1.0 | 1.0, 1.0, 1.0, 1.0
988 // field_id 1 | 16, 17, 18, 19 | 20, 21, 22, 23 | 24, 25, 26, 27 | 28, 29, 30, 31
989 // frac_mask | 1.0, 1.0, 1.0, 1.0 | 0.5, 0.0, 0.0, 0.0 | 1.0, 1.0, 1.0, 1.0 | 1.0, 1.0, 1.0, 1.0
990 //
991 // tgt id | src ids | optional weights
992 // -------|--------------|--------------------
993 // 4 | 0, 17, 2, 19 | 0.1, 0.2, 0.3, 0.4
994 // 5 | 2, 19, 20, 5 | 0.5, 0.6, 0.7, 0.8
995 // 6 | 6, 22, 7, 23 | 0.9, 1.0, 1.1, 1.2
996 // 7 | 23 | 1.3
997 // 8 | 2, 19, 22, 7 | 1.4, 1.5, 1.6, 1.7
998 // 9 | 0, 17, 2, 19 | 0.1, 0.2, 0.3, 0.4
999 // 10 | 5, 6 | 2.2, 2.3
1000 // 11 | 21, 22 | 2.4, 2.5
1001
1002 // at first the source points required for each target point have to be
1003 // transfered to the respective target processes using src_redists
1004 enum {num_src_fields = 2};
1005 Xt_redist src_redists[num_src_fields];
1006 {
1007 enum {num_src_indices = 4};
1008 Xt_int src_indices[num_procs][num_src_fields][num_src_indices] =
1009 {{{0,1,2,3},{16,17,18,19}},
1010 {{4,5,6,7},{20,21,22,23}},
1011 {{8,9,10,11},{24,25,26,27}},
1012 {{12,13,14,15},{28,29,30,31}}};
1013
1014 Xt_int dst_indices[num_procs][num_src_fields][5] =
1015 {{{-1},{-1}},
1016 {{0,2},{17,19}},
1017 {{0,2,5,6,7},{17,19,21,22}},
1018 {{-1},{-1}}};
1019 size_t num_dst_indices[num_procs][num_src_fields] =
1020 {{0,0},{2,2},{5,4},{0,0}};
1021 for (int field_id = 0; field_id < num_src_fields; ++field_id) {
1022 Xt_idxlist src_idxlist =
1023 xt_idxvec_new(src_indices[comm_rank][field_id], num_src_indices);
1024 Xt_idxlist dst_idxlist =
1025 xt_idxvec_new(
1026 dst_indices[comm_rank][field_id],
1027 num_dst_indices[comm_rank][field_id]);
1028 Xt_xmap xmap =
1029 xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
1030 src_redists[field_id] = xt_redist_p2p_new(xmap, MPI_DOUBLE);
1031 xt_xmap_delete(xmap);
1032 xt_idxlist_delete(dst_idxlist);
1033 xt_idxlist_delete(src_idxlist);
1034 }
1035 }
1036
1037 // positions of target points to be written
1038 size_t tgt_pos[num_procs][4] =
1039 {{(size_t)-1},{0,1,2,3},{0,1,2,3},{(size_t)-1}};
1040 // number of entries in tgt_pos
1041 size_t tgt_count[num_procs] = {0,4,4,0};
1042 // number of source points per stencil
1043 size_t num_src_per_tgt[num_procs][4] =
1044 {{(size_t)-1},{4,4,4,1},{4,4,2,2},{(size_t)-1}};
1045 double weights[num_procs][13] =
1046 {{-1},
1047 {0.1,0.2,0.3,0.4,
1048 0.5,0.6,0.7,0.8,
1049 0.9,1.0,1.1,1.2,
1050 1.3},
1051 {1.4,1.5,1.6,1.7,
1052 0.1,0.2,0.3,0.4,
1053 2.2,2.3,
1054 2.4,2.5},
1055 {-1}};
1056 // source field index for each source point used in the stencils (SIZE_MAX
1057 // indicates that the respective source point is in the halo buffer)
1058 size_t src_field_idx[num_procs][13] =
1059 {{(size_t)-1},
1060 {SIZE_MAX,SIZE_MAX,SIZE_MAX,SIZE_MAX,
1061 SIZE_MAX,SIZE_MAX,1,0,
1062 0,1,0,1,
1063 1},
1064 {SIZE_MAX,SIZE_MAX,SIZE_MAX,SIZE_MAX,
1065 SIZE_MAX,SIZE_MAX,SIZE_MAX,SIZE_MAX,
1066 SIZE_MAX,SIZE_MAX,
1067 SIZE_MAX,SIZE_MAX},
1068 {(size_t)-1}};
1069 // source index; position of the respective source point within the source
1070 // field or halo buffer
1071 size_t src_idx[num_procs][13] =
1072 {{(size_t)-1},
1073 {0,2,1,3, 1,3,0,1, 2,2,3,3, 3},
1074 {1,6,8,4, 0,5,1,6, 2,3, 7,8},
1075 {(size_t)-1}};
1076
1077 for (int scaling_factor_idx = 0;
1078 scaling_factor_idx < NUM_SCALING_FACTORS; ++scaling_factor_idx) {
1079 for (int scaling_summand_idx = 0;
1080 scaling_summand_idx < NUM_SCALING_SUMMAND; ++scaling_summand_idx) {
1081 for (int with_weights = 0; with_weights < 2; ++with_weights) {
1082
1083 enum {collection_size = 1};
1086 struct yac_interpolation * interp =
1089 scaling_factors[scaling_factor_idx],
1090 scaling_summand[scaling_summand_idx]);
1091 if (with_weights) {
1093 interp, src_redists, tgt_pos[comm_rank], tgt_count[comm_rank],
1094 num_src_per_tgt[comm_rank], weights[comm_rank],
1095 src_field_idx[comm_rank], src_idx[comm_rank], num_src_fields);
1096 } else {
1098 interp, src_redists, tgt_pos[comm_rank], tgt_count[comm_rank],
1099 num_src_per_tgt[comm_rank], src_field_idx[comm_rank],
1100 src_idx[comm_rank], num_src_fields);
1101 }
1102
1103 enum {src_field_size = 4};
1104 double src_data_raw[num_procs][num_src_fields][src_field_size] =
1105 {{{0,1,2,3}, {16,17,18,19}},
1106 {{4,5,6,7}, {20,21,22,23}},
1107 {{8,9,10,11}, {24,25,26,27}},
1108 {{12,13,14,15},{28,29,30,31}}};
1109 double src_frac_mask_raw[num_procs][2][src_field_size] =
1110 {{{1.0, 1.0, 1.0, 1.0}, {1.0, 1.0, 1.0, 1.0}},
1111 {{1.0, 0.5, 0.0, 0.0}, {0.5, 0.0, 0.0, 0.0}},
1112 {{1.0, 1.0, 1.0, 1.0}, {1.0, 1.0, 1.0, 1.0}},
1113 {{1.0, 1.0, 1.0, 1.0}, {1.0, 1.0, 1.0, 1.0}}};
1114 size_t const tgt_field_size[num_procs] = {0,4,4,0};
1115 double ref_tgt_data_raw[2][num_procs][collection_size][4] =
1116 {{{{UNSET}},
1117 {{SCALE((1.0*0+1.0*17+1.0*2+1.0*19)/(1.0+1.0+1.0+1.0)),
1118 SCALE((1.0*2+1.0*19+0.5*20+0.5*5)/(1.0+1.0+0.5+0.5)),
1120 {{SCALE((1.0*2+1.0*19+0.0*22+0.0*7)/(1.0+1.0+0.0+0.0)),
1121 SCALE((1.0*0+1.0*17+1.0*2+1.0*19)/(1.0+1.0+1.0+1.0)),
1122 SCALE((0.5*5+0.0*6)/(0.5+0.0)),
1124 {{UNSET}}},
1125 {{{UNSET}},
1126 {{SCALE((1.0*0.1*0+1.0*0.2*17+1.0*0.3*2+1.0*0.4*19)/
1127 (1.0*0.1+1.0*0.2+1.0*0.3+1.0*0.4)),
1128 SCALE((1.0*0.5*2+1.0*0.6*19+0.5*0.7*20+0.5*0.8*5)/
1129 (1.0*0.5+1.0*0.6+0.5*0.7+0.5*0.8)),
1131 {{SCALE((1.0*1.4*2+1.0*1.5*19+0.0*1.6*22+0.0*1.7*7)/
1132 (1.0*1.4+1.0*1.5+0.0*1.6+0.0*1.7)),
1133 SCALE((1.0*0.1*0+1.0*0.2*17+1.0*0.3*2+1.0*0.4*19)/
1134 (1.0*0.1+1.0*0.2+1.0*0.3+1.0*0.4)),
1135 SCALE((0.5*2.2*5+0.0*2.3*6)/(0.5*2.2+0.0*2.3)),
1137 {{UNSET}}}};
1138
1139 if (utest_check_interpolation_frac(
1140 interp,
1141 &src_data_raw[comm_rank][0][0],
1142 &src_frac_mask_raw[comm_rank][0][0],
1143 num_src_fields, src_field_size,
1144 &ref_tgt_data_raw[with_weights][comm_rank][0][0],
1145 tgt_field_size[comm_rank], collection_size,
1147 PUT_ERR("ERROR in yac_interpolation_add_sum_at_tgt/"
1148 "yac_interpolation_add_weight_sum_mvp_at_tgt frac");
1149
1152 }
1153 }
1154 }
1155 xt_redist_delete(src_redists[1]);
1156 xt_redist_delete(src_redists[0]);
1157 }
1158
1159 { // testing yac_interpolation_execute_get_test for a stack comprised of
1160 // a fixed and a direct interpolation
1161 // (fixed will return INACTIVE instantly, while direct will require the
1162 // associated put to be called. Results of both interpolation operations
1163 // were combined incorrectly at some point.)
1164
1165 enum {COLLECTION_SIZE = 1};
1166 double const scaling_factor = 1.0;
1167 double const scaling_summand = 0.0;
1170 struct yac_interpolation * interp =
1173 scaling_factor, scaling_summand);
1175
1176 // fixed:
1177 // pos "0" at rank "0" is set to 0.0
1178 double fixed_value = 4.0;
1179 size_t count[num_procs] = {1, 0, 0, 0};
1180 size_t pos[num_procs][1] = {{3},{SIZE_MAX},{SIZE_MAX},{SIZE_MAX}};
1182 interp, fixed_value, count[comm_rank], pos[comm_rank]);
1183
1184 // direct:
1185 // rank 1 pos {0} -> rank 0 pos {1}
1186 // rank 2 pos {0} -> rank 0 pos {2}
1187 // rank 3 pos {0} -> rank 0 pos {3}
1188 Xt_int src_indices[num_procs][1] = {{-1},{1},{2},{3}};
1189 size_t num_src_indices[num_procs] = {0,1,1,1};
1190 Xt_int dst_indices[num_procs][3] = {{1,2,3},{-1},{-1},{-1}};
1191 size_t num_dst_indices[num_procs] = {3,0,0,0};
1192 Xt_idxlist src_idxlist =
1193 xt_idxvec_new(src_indices[comm_rank], num_src_indices[comm_rank]);
1194 Xt_idxlist dst_idxlist =
1195 xt_idxvec_new(dst_indices[comm_rank], num_dst_indices[comm_rank]);
1196 Xt_xmap xmap
1197 = xt_xmap_all2all_new(src_idxlist, dst_idxlist, MPI_COMM_WORLD);
1198 Xt_redist redist = xt_redist_p2p_new(xmap, MPI_DOUBLE);
1199 xt_xmap_delete(xmap);
1200 xt_idxlist_delete(dst_idxlist);
1201 xt_idxlist_delete(src_idxlist);
1202 yac_interpolation_add_direct(interp, redist);
1203 xt_redist_delete(redist);
1204
1205 // setup source and target data
1206 enum {NUM_SRC_FIELDS = 1};
1207 double *** src_data = xmalloc(COLLECTION_SIZE * sizeof(src_data));
1208 src_data[0] = xmalloc(COLLECTION_SIZE * NUM_SRC_FIELDS * sizeof(*src_data));
1209 src_data[0][0] = xmalloc(1 * sizeof(**src_data));
1210 src_data[0][0][0] = (double)comm_rank;
1211 double tgt_data_raw[4] = {UNSET, UNSET, UNSET, UNSET};
1212 double ** tgt_data = xmalloc(COLLECTION_SIZE * sizeof(tgt_data));
1213 tgt_data[0] = tgt_data_raw;
1214
1215 // there should be not active interpolation exchange
1218 PUT_ERR("ERROR in yac_interpolation_execute_put/get_test");
1219
1220 // rank "0" is the only target
1221 if (comm_rank == 0) {
1222
1223 // initiate asynchronous get operation
1224 yac_interpolation_execute_get_async(interp, tgt_data);
1225
1226 // at this point no put should have been executed, fixed interpolation
1227 // operation will already be completed, but direct not
1228 // -> get_test should return ACTIVE
1230 PUT_ERR("ERROR in yac_interpolation_execute_get_test");
1231 }
1232
1233 // ensure that no put was called before get_test was tested
1234 MPI_Barrier(MPI_COMM_WORLD);
1235
1236 // rank "1", "2", and "3" are source
1237 if (comm_rank != 0) {
1238 yac_interpolation_execute_put(interp, src_data);
1239 }
1240
1241 // ensure that all operations have been completed
1243
1244 // there should be not active interpolation exchange
1247 PUT_ERR("ERROR in yac_interpolation_execute_put/get_test");
1248
1249 // check results
1250 double ref_tgt_data_raw[num_procs][4] =
1251 {{1.0, 2.0, 3.0, 4.0}, {UNSET, UNSET, UNSET, UNSET},
1253 for (size_t i = 0; i < sizeof(tgt_data_raw)/sizeof(tgt_data_raw[0]); ++i) {
1254 if (tgt_data_raw[i] != ref_tgt_data_raw[comm_rank][i]) {
1255 PUT_ERR("ERROR in interpolation results");
1256 }
1257 }
1258
1259 // cleanup
1260 free(tgt_data);
1261 free(src_data[0][0]);
1262 free(src_data[0]);
1263 free(src_data);
1265 }
1266
1267 xt_finalize();
1268 MPI_Finalize();
1269
1270 return TEST_EXIT_CODE;
1271}
1272
1273static void init_tgt_data(
1274 double ** tgt_data, size_t tgt_field_size, size_t collection_size) {
1275
1276 if (tgt_data == NULL) return;
1277
1278 for (size_t i = 0; i < collection_size; ++i)
1279 for (size_t j = 0; j < tgt_field_size; ++j)
1280 tgt_data[i][j] = UNSET;
1281}
1282
1284 double ** tgt_data, double * ref_tgt_data_raw,
1285 size_t tgt_field_size, size_t collection_size) {
1286
1287 int diff_count = 0;
1288 if (tgt_data != NULL)
1289 for (size_t i = 0; i < collection_size; ++i)
1290 for (size_t j = 0; j < tgt_field_size; ++j)
1291 if (fabs(tgt_data[i][j] - ref_tgt_data_raw[j]) > 1e-9)
1292 ++diff_count;
1293
1294 return diff_count;
1295}
1296
1297static int utest_check_interpolation_execute(
1298 struct yac_interpolation * interp, double *** src_data, double ** tgt_data,
1299 double * ref_tgt_data_raw, size_t tgt_field_size, size_t collection_size) {
1300
1301 int diff_count = 0;
1302
1303 // there should be not active interpolation exchange
1306 PUT_ERR("ERROR in yac_interpolation_execute_put/get_test");
1307
1308 // synchronous exchange
1309 init_tgt_data(tgt_data, tgt_field_size, collection_size);
1310 yac_interpolation_execute(interp, src_data, tgt_data);
1311 diff_count +=
1312 check_tgt_data(tgt_data, ref_tgt_data_raw, tgt_field_size, collection_size);
1313
1314 MPI_Barrier(MPI_COMM_WORLD);
1315
1316 // there should be no active interpolation exchange
1319 PUT_ERR("ERROR in yac_interpolation_execute_put/get_test");
1320
1321 // first put then get
1322 if (src_data != NULL) {
1323 yac_interpolation_execute_put(interp, src_data);
1324
1325 // this should return at some point
1326 if (tgt_data == NULL)
1327 while (!yac_interpolation_execute_put_test(interp));
1328 }
1329 if (tgt_data != NULL) yac_interpolation_execute_get(interp, tgt_data);
1330 diff_count +=
1331 check_tgt_data(tgt_data, ref_tgt_data_raw, tgt_field_size, collection_size);
1332
1333 // first async get then put
1334 if (tgt_data != NULL)
1335 yac_interpolation_execute_get_async(interp, tgt_data);
1336 if (src_data != NULL)
1337 yac_interpolation_execute_put(interp, src_data);
1338 // this should return at some point
1339 while (!yac_interpolation_execute_put_test(interp));
1340 diff_count +=
1341 check_tgt_data(tgt_data, ref_tgt_data_raw, tgt_field_size, collection_size);
1342
1343 // first async get then put
1344 if (tgt_data != NULL)
1345 yac_interpolation_execute_get_async(interp, tgt_data);
1346 if (src_data != NULL)
1347 yac_interpolation_execute_put(interp, src_data);
1348 // this should return at some point
1350 check_tgt_data(tgt_data, ref_tgt_data_raw, tgt_field_size, collection_size);
1351
1352 return diff_count;
1353}
1354
1355static int utest_check_interpolation_execute_frac(
1356 struct yac_interpolation * interp, double *** src_data,
1357 double *** src_frac_mask, double ** tgt_data, double * ref_tgt_data_raw,
1358 size_t tgt_field_size, size_t collection_size) {
1359
1360 int diff_count = 0;
1361
1362 // there should be not active interpolation exchange
1365 PUT_ERR("ERROR in yac_interpolation_execute_put/get_test");
1366
1367 // synchronous exchange
1368 init_tgt_data(tgt_data, tgt_field_size, collection_size);
1369 yac_interpolation_execute_frac(interp, src_data, src_frac_mask, tgt_data);
1370 diff_count +=
1371 check_tgt_data(tgt_data, ref_tgt_data_raw, tgt_field_size, collection_size);
1372
1373 MPI_Barrier(MPI_COMM_WORLD);
1374
1375 // there should be not active interpolation exchange
1378 PUT_ERR("ERROR in yac_interpolation_execute_put/get_test");
1379
1380 // first put then get
1381 if (src_data != NULL) {
1382 yac_interpolation_execute_put_frac(interp, src_data, src_frac_mask);
1383 // this should return at some point
1384 if (tgt_data == NULL)
1385 while (!yac_interpolation_execute_put_test(interp));
1386 }
1387 if (tgt_data != NULL) yac_interpolation_execute_get(interp, tgt_data);
1388 diff_count +=
1389 check_tgt_data(tgt_data, ref_tgt_data_raw, tgt_field_size, collection_size);
1390
1391 return diff_count;
1392}
1393
1394static int utest_check_interpolation(
1395 struct yac_interpolation * interp,
1396 double * src_data_raw, size_t num_src_fields, size_t src_field_size,
1397 double * ref_tgt_data_raw, size_t tgt_field_size, size_t src_collection_size,
1399
1400 // generate source input data pointers
1401 double *** src_data = NULL;
1402 if (src_data_raw != NULL) {
1403 src_data = xmalloc(src_collection_size * sizeof(src_data));
1404 double ** src_field_data =
1405 xmalloc(src_collection_size * num_src_fields * sizeof(*src_field_data));
1406 for (size_t i = 0; i < src_collection_size;
1407 ++i, src_field_data += num_src_fields) {
1408 src_data[i] = src_field_data;
1409 for (size_t j = 0, src_field_offset = 0; j < num_src_fields;
1410 ++j, src_field_offset += src_field_size)
1411 src_data[i][j] = src_data_raw + src_field_offset;
1412 }
1413 }
1414
1415 // generate target input data pointers
1416 size_t tgt_collection_size =
1418 double ** tgt_data = NULL;
1419 if (ref_tgt_data_raw != NULL) {
1420 double * tgt_data_raw =
1421 xmalloc(tgt_collection_size * tgt_field_size * sizeof(*tgt_data_raw));
1422 for (size_t i = 0; i < tgt_collection_size * tgt_field_size; ++i)
1423 tgt_data_raw[i] = UNSET;
1424 tgt_data = xmalloc(tgt_collection_size * sizeof(*tgt_data));
1425 for (size_t i = 0; i < tgt_collection_size; ++i,
1426 tgt_data_raw += tgt_field_size)
1427 tgt_data[i] = tgt_data_raw;
1428 }
1429
1430 // check results
1431 int diff_count =
1432 utest_check_interpolation_execute(
1433 interp, src_data, tgt_data, ref_tgt_data_raw,
1434 tgt_field_size, tgt_collection_size);
1435
1436 struct yac_interpolation * interp_copy = yac_interpolation_copy(interp);
1437 diff_count +=
1438 utest_check_interpolation_execute(
1439 interp, src_data, tgt_data, ref_tgt_data_raw,
1440 tgt_field_size, tgt_collection_size);
1441 yac_interpolation_delete(interp_copy);
1442
1443 diff_count +=
1444 abs(
1448
1449 // clean up
1450 if (ref_tgt_data_raw != NULL) {
1451 free(tgt_data[0]);
1452 free(tgt_data);
1453 }
1454 if (src_data_raw != NULL) {
1455 free(src_data[0]);
1456 free(src_data);
1457 }
1458
1459 return diff_count;
1460}
1461
1462static int utest_check_interpolation_frac(
1463 struct yac_interpolation * interp,
1464 double * src_data_raw, double * src_frac_mask_raw,
1465 size_t num_src_fields, size_t src_field_size,
1466 double * ref_tgt_data_raw, size_t tgt_field_size, size_t src_collection_size,
1468
1469 // generate source input data pointers
1470 double *** src_data = NULL;
1471 if (src_data_raw != NULL) {
1472 src_data = xmalloc(src_collection_size * sizeof(src_data));
1473 double ** src_field_data =
1474 xmalloc(src_collection_size * num_src_fields * sizeof(*src_field_data));
1475 for (size_t i = 0; i < src_collection_size;
1476 ++i, src_field_data += num_src_fields) {
1477 src_data[i] = src_field_data;
1478 for (size_t j = 0, src_field_offset = 0; j < num_src_fields;
1479 ++j, src_field_offset += src_field_size)
1480 src_data[i][j] = src_data_raw + src_field_offset;
1481 }
1482 }
1483
1484 // generate source fractional mask data pointers
1485 size_t tgt_collection_size =
1487 double *** src_frac_data = NULL;
1488 if (src_frac_mask_raw != NULL) {
1489 src_frac_data = xmalloc(tgt_collection_size * sizeof(src_frac_data));
1490 double ** src_frac_mask_data =
1491 xmalloc(
1492 tgt_collection_size * num_src_fields * sizeof(*src_frac_mask_data));
1493 for (size_t i = 0; i < tgt_collection_size;
1494 ++i, src_frac_mask_data += num_src_fields) {
1495 src_frac_data[i] = src_frac_mask_data;
1496 for (size_t j = 0, src_field_offset = 0; j < num_src_fields;
1497 ++j, src_field_offset += src_field_size)
1498 src_frac_data[i][j] = src_frac_mask_raw + src_field_offset;
1499 }
1500 }
1501
1502 // apply fractional mask to source data
1503 if ((src_data != NULL) && (src_frac_data != NULL))
1504 for (size_t i = 0; i < src_collection_size; ++i)
1505 for (size_t j = 0; j < num_src_fields; ++j)
1506 for (size_t k = 0; k < src_field_size; ++k)
1507 src_data[i][j][k] *= src_frac_data[i][j][k];
1508
1509 // generate target input data pointers
1510 double ** tgt_data = NULL;
1511 if (ref_tgt_data_raw != NULL) {
1512 double * tgt_data_raw =
1513 xmalloc(tgt_collection_size * tgt_field_size * sizeof(*tgt_data_raw));
1514 for (size_t i = 0; i < tgt_collection_size * tgt_field_size; ++i)
1515 tgt_data_raw[i] = UNSET;
1516 tgt_data = xmalloc(tgt_collection_size * sizeof(*tgt_data));
1517 for (size_t i = 0; i < tgt_collection_size;
1518 ++i, tgt_data_raw += tgt_field_size)
1519 tgt_data[i] = tgt_data_raw;
1520 }
1521
1522 // check results
1523 int diff_count =
1524 utest_check_interpolation_execute_frac(
1525 interp, src_data, src_frac_data, tgt_data, ref_tgt_data_raw,
1526 tgt_field_size, tgt_collection_size);
1527
1528 struct yac_interpolation * interp_copy = yac_interpolation_copy(interp);
1529 diff_count +=
1530 utest_check_interpolation_execute_frac(
1531 interp, src_data, src_frac_data, tgt_data, ref_tgt_data_raw,
1532 tgt_field_size, tgt_collection_size);
1533 yac_interpolation_delete(interp_copy);
1534
1535 // clean up
1536 if (ref_tgt_data_raw != NULL) {
1537 free(tgt_data[0]);
1538 free(tgt_data);
1539 }
1540 if (src_frac_mask_raw != NULL) {
1541 free(src_frac_data[0]);
1542 free(src_frac_data);
1543 }
1544 if (src_data_raw != NULL) {
1545 free(src_data[0]);
1546 free(src_data);
1547 }
1548
1549 return diff_count;
1550}
size_t yac_collection_selection_get_collection_size(struct yac_collection_selection const *collection_selection)
Get the size of the collection selection.
void yac_collection_selection_delete(struct yac_collection_selection *collection_selection)
Delete a collection selection object.
int yac_collection_selection_compare(struct yac_collection_selection const *a, struct yac_collection_selection const *b)
Compare two collection selections.
struct yac_collection_selection * yac_collection_selection_new(size_t collection_size, size_t const *selection_indices)
Create a new collection selection.
struct @34::@35 value
void yac_interpolation_add_sum_at_src(struct yac_interpolation *interp, Xt_redist *halo_redists, size_t tgt_count, size_t *num_src_per_tgt, size_t *src_field_idx, size_t *src_idx, size_t num_src_fields, Xt_redist result_redist)
Add a sum operator where accumulation occurs on source processes.
struct yac_interpolation * yac_interpolation_copy(struct yac_interpolation *interp)
Create a deep copy of an interpolation object.
void yac_interpolation_add_weight_sum_mvp_at_tgt(struct yac_interpolation *interp, Xt_redist *src_redists, size_t *tgt_pos, size_t tgt_count, size_t *num_src_per_tgt, double *weights, size_t *src_field_idx, size_t *src_idx, size_t num_src_fields)
Add a weighted sum operator (distributed matrix-vector product), which computes the product at the ta...
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_add_direct_mf(struct yac_interpolation *interp, Xt_redist *redists, size_t num_src_fields)
Add a direct redistribution operator for multiple source fields.
void yac_interpolation_delete(struct yac_interpolation *interp)
Free an interpolation object and release all resources.
void yac_interpolation_add_sum_at_tgt(struct yac_interpolation *interp, Xt_redist *src_redists, size_t *tgt_pos, size_t tgt_count, size_t *num_src_per_tgt, size_t *src_field_idx, size_t *src_idx, size_t num_src_fields)
Add a sum operator where accumulation occurs on target processes.
int yac_interpolation_execute_put_test(struct yac_interpolation *interp)
Test whether the asynchronous put phase has completed.
void yac_interpolation_add_fixed(struct yac_interpolation *interp, double value, size_t count, size_t *pos)
Add a fixed-value operator to an interpolation.
void yac_interpolation_execute_wait(struct yac_interpolation *interp)
Wait for completion of pending asynchronous interpolation operations.
void yac_interpolation_execute_frac(struct yac_interpolation *interp, double ***src_fields, double ***src_frac_masks, double **tgt_field)
Execute interpolation with fractional masks and write results to the target field.
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_get_async(struct yac_interpolation *interp, double **tgt_field)
Complete interpolation asynchronously and write results to the target field (get phase).
void yac_interpolation_add_direct(struct yac_interpolation *interp, Xt_redist redist)
Add a direct redistribution operator.
struct yac_interpolation * yac_interpolation_new(struct yac_collection_selection const *collection_selection, double frac_mask_fallback_value, double scale_factor, double scale_summand)
Create a new interpolation object.
void yac_interpolation_execute_put_frac(struct yac_interpolation *interp, double ***src_fields, double ***src_frac_masks)
Provide source field data with fractional masks and start asynchronous execution of interpolation (pu...
void yac_interpolation_execute_put(struct yac_interpolation *interp, double ***src_fields)
Provide source field data and start asynchronous execution of interpolation (put phase).
void yac_interpolation_add_weight_sum_mvp_at_src(struct yac_interpolation *interp, Xt_redist *halo_redists, size_t tgt_count, size_t *num_src_per_tgt, double *weights, size_t *src_field_idx, size_t *src_idx, size_t num_src_fields, Xt_redist result_redist)
Add a weighted sum operator (distributed matrix-vector product), which computes the productes at the ...
int yac_interpolation_execute_get_test(struct yac_interpolation *interp)
Test whether the asynchronous get phase has completed.
double const YAC_FRAC_MASK_NO_VALUE
struct yac_collection_selection const * yac_interpolation_get_collection_selection(struct yac_interpolation const *interp)
Return the collection selection associated with an interpolation.
Defines internal basic interpolation definitions.
#define xmalloc(size)
Definition ppm_xfuncs.h:66
struct yac_collection_selection * collection_selection
Selection of field collections to which this interpolation applies.
@ COLLECTION_SIZE
int collection_size
#define N
static int check_tgt_data(double **tgt_data, double *ref_tgt_data_raw, size_t tgt_field_size, size_t collection_size)
#define FRAC_MASK_VALUE
static void init_tgt_data(double **tgt_data, size_t tgt_field_size, size_t collection_size)
#define SCALE(X)
#define TEST_EXIT_CODE
Definition tests.h:14
#define PUT_ERR(string)
Definition tests.h:10
int * field_id