YAC 3.13.0
Yet Another Coupler
Loading...
Searching...
No Matches
test_interpolation_parallel1.F90
Go to the documentation of this file.
1! Copyright (c) 2024 The YAC Authors
2!
3! SPDX-License-Identifier: BSD-3-Clause
4
8
9#include "test_macros.inc"
10
11PROGRAM main
12
13 USE utest
14 USE yac_core
15 USE mpi
16 use, INTRINSIC :: iso_c_binding
17
18 IMPLICIT NONE
19
20 INTERFACE
21
22 ! fortran interface to C routine from dist_grid_utils
23 FUNCTION yac_generate_basic_grid_reg2d_c( &
24 name, coordinates_x, coordinates_y, num_cells, &
25 local_start, local_count, with_halo) &
26 bind( c, name='yac_generate_basic_grid_reg2d' )
27
28 use, INTRINSIC :: iso_c_binding
29
30 CHARACTER(KIND=c_char) :: name(*)
31 REAL(kind=c_double) :: coordinates_x(*)
32 REAL(kind=c_double) :: coordinates_y(*)
33 INTEGER(kind=c_size_t) :: num_cells(2)
34 INTEGER(kind=c_size_t) :: local_start(2)
35 INTEGER(kind=c_size_t) :: local_count(2)
36 INTEGER(kind=c_int), value :: with_halo
37
38 TYPE(c_ptr) :: yac_generate_basic_grid_reg2d_c
39 END FUNCTION yac_generate_basic_grid_reg2d_c
40 END INTERFACE
41
42 CHARACTER(LEN=*), PARAMETER :: grid_name_src = "src_grid"
43 CHARACTER(LEN=*), PARAMETER :: grid_name_tgt = "tgt_grid"
44 INTEGER(kind=c_size_t), PARAMETER :: collection_size = 10
45
46 INTEGER, PARAMETER :: max_opt_arg_len = 1024
47 CHARACTER(max_opt_arg_len) :: reorder_type_arg
48 INTEGER(KIND=c_int) :: reorder_type
49 INTEGER :: arg_len
50
51 INTEGER :: comm_rank, comm_size, ierror
52
53 LOGICAL :: tgt_flag
54
55 INTEGER :: split_comm
56
57 REAL(KIND=c_double), PARAMETER :: &
58 yac_rad = 0.01745329251994329576923690768489_c_double ! M_PI / 180
59
60 CALL start_test('interpolation_parallel1')
61
62 CALL yac_mpi_init_c()
63 CALL yac_yaxt_init_c(mpi_comm_world)
64
65 CALL test(command_argument_count() == 1)
66 CALL get_command_argument(1, reorder_type_arg, arg_len)
67 CALL test((reorder_type_arg == "src") .OR. (reorder_type_arg == "tgt"))
68 reorder_type = &
69 merge(yac_mapping_on_src, yac_mapping_on_tgt, reorder_type_arg == "src")
70
71 CALL mpi_comm_rank(mpi_comm_world, comm_rank, ierror)
72 CALL mpi_comm_size(mpi_comm_world, comm_size, ierror)
73 CALL test(comm_size == 4)
74 CALL mpi_barrier(mpi_comm_world, ierror)
75
76 ! split processes into source and target
77
78 tgt_flag = comm_rank < 2
79
80 CALL mpi_comm_split( &
81 mpi_comm_world, merge(1,0,tgt_flag), 0, split_comm, ierror)
82
83 IF (tgt_flag) THEN
84 CALL target_main(mpi_comm_world, split_comm, reorder_type)
85 ELSE
86 CALL source_main(mpi_comm_world, split_comm, reorder_type)
87 END IF
88
89 CALL mpi_comm_free(split_comm, ierror)
90
92
93 CALL stop_test
94 CALL exit_tests
95
96CONTAINS
97
98 !
99 ! The source grid is distributed among 2 processes
100 !
101 ! The global source grid has 5x4 cells:
102 !
103 ! 24--44--25--45--26--46--27--47--28--48--29
104 ! | | | | | |
105 ! 34 15 36 16 38 17 40 18 42 19 43
106 ! | | | | | |
107 ! 18--33--19--35--20--37--21--39--22--41--23
108 ! | | | | | |
109 ! 23 10 25 11 27 12 29 13 31 14 32
110 ! | | | | | |
111 ! 12--22--13--24--14--26--15--28--16--30--17
112 ! | | | | | |
113 ! 12 05 14 06 16 07 18 08 20 09 21
114 ! | | | | | |
115 ! 06--11--07--13--08--15--09--17--10--19--11
116 ! | | | | | |
117 ! 01 00 03 01 05 02 07 03 09 04 10
118 ! | | | | | |
119 ! 00--01--01--02--02--04--03--06--04--08--05
120 !
121 SUBROUTINE source_main(global_comm, source_comm, reorder_type)
122
123 INTEGER, INTENT(IN) :: global_comm
124 INTEGER, INTENT(IN) :: source_comm
125 INTEGER(KIND=c_int) :: reorder_type
126
127 INTEGER :: my_source_rank, ierror
128 INTEGER(KIND=c_size_t) :: i
129
130 REAL(KIND=c_double), PARAMETER :: &
131 coordinates_x(6) = (/0.0_c_double, 1.0_c_double, 2.0_c_double, &
132 3.0_c_double, 4.0_c_double, 5.0_c_double/) * yac_rad
133 REAL(KIND=c_double), PARAMETER :: &
134 coordinates_y(5) = (/0.0_c_double, 1.0_c_double, 2.0_c_double, &
135 3.0_c_double, 4.0_c_double/) * yac_rad
136 INTEGER(KIND=c_size_t), PARAMETER :: num_cells(2) = (/5,4/)
137 INTEGER(KIND=c_size_t), PARAMETER :: &
138 local_start(2,2) = reshape((/0,0,3,0/),(/2,2/))
139 INTEGER(KIND=c_size_t), PARAMETER :: &
140 local_count(2,2) = reshape((/3,4,2,4/),(/2,2/))
141 INTEGER(KIND=c_int), PARAMETER :: with_halo = 1
142
143 INTEGER, PARAMETER :: &
144 src_data_int(25,2) = &
145 reshape( &
146 (/0,1,2,3,-1, &
147 6,7,8,9,-1, &
148 12,13,14,15,-1, &
149 18,19,20,21,-1, &
150 24,25,26,27,-1, &
151 -1,3,4,5, &
152 -1,9,10,11, &
153 -1,15,16,17, &
154 -1,21,22,23, &
155 -1,27,28,29, &
156 -1,-1,-1,-1,-1/), (/25,2/))
157 REAL(KIND=c_double), ALLOCATABLE, TARGET :: src_data(:,:)
158 INTEGER(KIND=c_size_t) :: data_size
159 TYPE(c_ptr), TARGET :: src_field_(collection_size)
160 TYPE(c_ptr), TARGET :: src_field_collection(collection_size)
161 TYPE(c_ptr) :: src_field
162
163 TYPE(c_ptr) :: src_grid, tgt_grid
164 TYPE(c_ptr) :: grid_pair
165 TYPE(c_ptr) :: interp_grid
166 TYPE(c_ptr) :: interp_stack_config
167 TYPE(c_ptr) :: interp_method_stack
168 TYPE(c_ptr) :: interp_weights
169 TYPE(c_ptr) :: interpolation
170
171 CALL mpi_comm_rank(source_comm, my_source_rank, ierror)
172
173 src_grid = &
174 yac_generate_basic_grid_reg2d_c( &
175 trim(grid_name_src) // c_null_char, &
176 coordinates_x, coordinates_y, num_cells, &
177 local_start(:, my_source_rank + 1), &
178 local_count(:, my_source_rank + 1), with_halo)
179 tgt_grid = &
180 yac_basic_grid_empty_new_c(trim(grid_name_tgt) // c_null_char)
181
182 grid_pair = yac_dist_grid_pair_new_c(src_grid, tgt_grid, global_comm)
183
184 interp_grid = &
186 grid_pair, trim(grid_name_src) // c_null_char, &
187 trim(grid_name_tgt) // c_null_char, 1_c_size_t, &
188 (/yac_loc_corner/), (/-1_c_size_t/), (/-1_c_size_t/), &
189 yac_loc_corner, -1_c_size_t, -1_c_size_t)
190
191 interp_stack_config = yac_interp_stack_config_new_c()
193 interp_stack_config, yac_interp_avg_arithmetic, 0_c_int)
195 interp_stack_config, 1337.0_c_double)
196
197 interp_method_stack = &
198 yac_interp_stack_config_generate_c(interp_stack_config)
199
200 interp_weights = &
201 yac_interp_method_do_search_c(interp_method_stack, interp_grid)
202
203 interpolation = &
205 interp_weights, reorder_type, collection_size, &
207 1.0_c_double, 0.0_c_double, c_null_char, 1_c_int, 1_c_int)
208
209 ! -------------------
210 ! set up source data
211 ! -------------------
212
213 ! src_data dimensions [collection_idx]
214 ! [pointset_idx]
215 ! [local_idx]
217 ALLOCATE(src_data(data_size, collection_size))
218 DO i = 1, collection_size
219 src_data(1:data_size,i) = &
220 REAL( &
221 src_data_int(1:data_size,my_source_rank+1) + INT(i * 30), c_double)
222 src_field_(i) = c_loc(src_data(:,i))
223 src_field_collection(i) = c_loc(src_field_(i))
224 END DO
225 src_field = c_loc(src_field_collection)
226
227 CALL yac_interpolation_execute_put_c(interpolation, src_field)
228
229 DEALLOCATE(src_data)
230 CALL yac_interpolation_delete_c(interpolation)
231 CALL yac_interp_weights_delete_c(interp_weights)
232 CALL yac_interp_method_delete_c(interp_method_stack);
233 CALL yac_interp_stack_config_delete_c(interp_stack_config)
234 CALL yac_interp_grid_delete_c(interp_grid)
235 CALL yac_dist_grid_pair_delete_c(grid_pair)
236 CALL yac_basic_grid_delete_c(tgt_grid)
237 CALL yac_basic_grid_delete_c(src_grid)
238
239 END SUBROUTINE
240
241 !
242 ! The target grid is distributed among 2 processes
243 !
244 ! The global target grid has 6x3 cells:
245 !
246 ! 21--39--22--40--23--41--24--42--25--43--26--44--27
247 ! | | | | | | |
248 ! 27 12 29 13 31 14 33 15 35 16 37 17 38
249 ! | | | | | | |
250 ! 14--26--15--28--16--30--17--32--18--34--19--36--20
251 ! | | | | | | |
252 ! 14 06 16 07 18 08 20 09 22 10 24 11 25
253 ! | | | | | | |
254 ! 07--13--08--15--09--17--10--19--11--21--12--23--13
255 ! | | | | | | |
256 ! 01 00 03 01 05 02 07 03 09 04 11 05 12
257 ! | | | | | | |
258 ! 00--00--01--02--02--04--03--06--04--08--05--10--06
259 !
260 SUBROUTINE target_main(global_comm, target_comm, reorder_type)
261
262
263 INTEGER, INTENT(IN) :: global_comm
264 INTEGER, INTENT(IN) :: target_comm
265 INTEGER(KIND=c_int) :: reorder_type
266
267 INTEGER :: my_target_rank, ierror
268 INTEGER(KIND=c_size_t) :: i, j
269
270 REAL(KIND=c_double), PARAMETER :: &
271 coordinates_x(7) = (/0.5_c_double, 1.5_c_double, 2.5_c_double, &
272 3.5_c_double, 4.5_c_double, 5.5_c_double, &
273 6.5_c_double/) * yac_rad
274 REAL(KIND=c_double), PARAMETER :: &
275 coordinates_y(4) = (/0.5_c_double, 1.5_c_double, &
276 2.5_c_double, 3.5_c_double/) * yac_rad
277 INTEGER(KIND=c_size_t), PARAMETER :: num_cells(2) = (/6,3/)
278 INTEGER(KIND=c_size_t), PARAMETER :: &
279 local_start(2,2) = reshape((/0,0,3,0/),(/2,2/))
280 INTEGER(KIND=c_size_t), PARAMETER :: &
281 local_count(2,2) = reshape((/3,3,3,3/),(/2,2/))
282 INTEGER(KIND=c_int), PARAMETER :: with_halo = 0
283
284 REAL(KIND=c_double), TARGET :: tgt_data(16, collection_size)
285 REAL(KIND=c_double), PARAMETER :: &
286 ref_tgt_data(16, 2) = &
287 reshape( &
288 (/3.5_c_double,4.5_c_double,5.5_c_double,6.5_c_double, &
289 9.5_c_double,10.5_c_double,11.5_c_double,12.5_c_double, &
290 15.5_c_double,16.5_c_double,17.5_c_double,18.5_c_double, &
291 21.5_c_double,22.5_c_double,23.5_c_double,24.5_c_double, &
292 6.5_c_double,7.5_c_double,1337.0_c_double,1337.0_c_double, &
293 12.5_c_double,13.5_c_double,1337.0_c_double,1337.0_c_double, &
294 18.5_c_double,19.5_c_double,1337.0_c_double,1337.0_c_double, &
295 24.5_c_double,25.5_c_double,1337.0_c_double,1337.0_c_double/), &
296 (/16,2/))
297 TYPE(c_ptr), TARGET :: tgt_field_collection(collection_size)
298 TYPE(c_ptr) :: tgt_field
299 REAL(kind=c_double) :: ref_tgt_value
300
301 TYPE(c_ptr) :: src_grid, tgt_grid
302 TYPE(c_ptr) :: grid_pair
303 TYPE(c_ptr) :: interp_grid
304 TYPE(c_ptr) :: interp_stack_config
305 TYPE(c_ptr) :: interp_method_stack
306 TYPE(c_ptr) :: interp_weights
307 TYPE(c_ptr) :: interpolation
308
309 CALL mpi_comm_rank(target_comm, my_target_rank, ierror)
310
311 tgt_grid = &
312 yac_generate_basic_grid_reg2d_c( &
313 trim(grid_name_tgt) // c_null_char, &
314 coordinates_x, coordinates_y, num_cells, &
315 local_start(:, my_target_rank + 1), &
316 local_count(:, my_target_rank + 1), with_halo)
317 src_grid = &
318 yac_basic_grid_empty_new_c(trim(grid_name_src) // c_null_char)
319
320 grid_pair = yac_dist_grid_pair_new_c(tgt_grid, src_grid, global_comm)
321
322 interp_grid = &
324 grid_pair, trim(grid_name_src) // c_null_char, &
325 trim(grid_name_tgt) // c_null_char, 1_c_size_t, &
326 (/yac_loc_corner/), (/-1_c_size_t/), (/-1_c_size_t/), &
327 yac_loc_corner, -1_c_size_t, -1_c_size_t)
328
329 interp_stack_config = yac_interp_stack_config_new_c()
331 interp_stack_config, yac_interp_avg_arithmetic, 0_c_int)
333 interp_stack_config, 1337.0_c_double)
334
335 interp_method_stack = &
336 yac_interp_stack_config_generate_c(interp_stack_config)
337
338 interp_weights = &
339 yac_interp_method_do_search_c(interp_method_stack, interp_grid)
340
341 interpolation = &
343 interp_weights, reorder_type, collection_size, &
345 1.0_c_double, 0.0_c_double, c_null_char, 1_c_int, 1_c_int)
346
347 tgt_data = -1.0_c_double
348 DO i = 1, collection_size
349 tgt_field_collection(i) = c_loc(tgt_data(:,i))
350 END DO
351 tgt_field = c_loc(tgt_field_collection)
352
353 CALL yac_interpolation_execute_get_c(interpolation, tgt_field)
354
355 DO i = 1, collection_size
356 DO j = 1, 16_c_size_t
357 IF (ref_tgt_data(j, my_target_rank+1) == 1337.0_c_double) THEN
358 ref_tgt_value = 1337.0_c_double
359 ELSE
360 ref_tgt_value = &
361 ref_tgt_data(j, my_target_rank+1) + &
362 REAL(i * 30_c_size_t, c_double)
363 END IF
364 CALL test(tgt_data(j,i) == ref_tgt_value)
365 END DO
366 END DO
367
368 CALL yac_interpolation_delete_c(interpolation)
369 CALL yac_interp_weights_delete_c(interp_weights)
370 CALL yac_interp_method_delete_c(interp_method_stack);
371 CALL yac_interp_stack_config_delete_c(interp_stack_config)
372 CALL yac_interp_grid_delete_c(interp_grid)
373 CALL yac_dist_grid_pair_delete_c(grid_pair)
374 CALL yac_basic_grid_delete_c(src_grid)
375 CALL yac_basic_grid_delete_c(tgt_grid)
376
377 END SUBROUTINE
378
379END PROGRAM
static void merge(char *base_a, size_t num_a, int a_ascending, char *base_b, size_t num_b, int b_ascending, size_t size, int(*compar)(const void *, const void *), char *target)
Definition mergesort.c:56
Definition utest.F90:5
subroutine, public start_test(name)
Definition utest.F90:20
subroutine, public stop_test()
Definition utest.F90:27
subroutine, public exit_tests()
Definition utest.F90:81
@ yac_loc_corner
Definition yac_core.F90:34
@ yac_interp_avg_arithmetic
Definition yac_core.F90:53
@ yac_mapping_on_src
Definition yac_core.F90:162
@ yac_mapping_on_tgt
Definition yac_core.F90:163