1 | MODULE pmc_parent |
---|
2 | |
---|
3 | !------------------------------------------------------------------------------! |
---|
4 | ! This file is part of PALM. |
---|
5 | ! |
---|
6 | ! PALM is free software: you can redistribute it and/or modify it under the |
---|
7 | ! terms of the GNU General Public License as published by the Free Software |
---|
8 | ! Foundation, either version 3 of the License, or (at your option) any later |
---|
9 | ! version. |
---|
10 | ! |
---|
11 | ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY |
---|
12 | ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR |
---|
13 | ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. |
---|
14 | ! |
---|
15 | ! You should have received a copy of the GNU General Public License along with |
---|
16 | ! PALM. If not, see <http://www.gnu.org/licenses/>. |
---|
17 | ! |
---|
18 | ! Copyright 1997-2017 Leibniz Universitaet Hannover |
---|
19 | !------------------------------------------------------------------------------! |
---|
20 | ! |
---|
21 | ! Current revisions: |
---|
22 | ! ------------------ |
---|
23 | ! |
---|
24 | ! |
---|
25 | ! Former revisions: |
---|
26 | ! ----------------- |
---|
27 | ! $Id: pmc_parent_mod.f90 2101 2017-01-05 16:42:31Z Giersch $ |
---|
28 | ! |
---|
29 | ! 2000 2016-08-20 18:09:15Z knoop |
---|
30 | ! Forced header and separation lines into 80 columns |
---|
31 | ! |
---|
32 | ! 1938 2016-06-13 15:26:05Z hellstea |
---|
33 | ! Minor clean up. |
---|
34 | ! |
---|
35 | ! 1901 2016-05-04 15:39:38Z raasch |
---|
36 | ! Module renamed. Code clean up. The words server/client changed to parent/child. |
---|
37 | ! |
---|
38 | ! 1900 2016-05-04 15:27:53Z raasch |
---|
39 | ! re-formatted to match PALM style |
---|
40 | ! |
---|
41 | ! 1850 2016-04-08 13:29:27Z maronga |
---|
42 | ! Module renamed |
---|
43 | ! |
---|
44 | ! |
---|
45 | ! 1833 2016-04-07 14:23:03Z raasch |
---|
46 | ! gfortran requires pointer attributes for some array declarations, |
---|
47 | ! long line wrapped |
---|
48 | ! |
---|
49 | ! 1808 2016-04-05 19:44:00Z raasch |
---|
50 | ! MPI module used by default on all machines |
---|
51 | ! |
---|
52 | ! 1797 2016-03-21 16:50:28Z raasch |
---|
53 | ! introduction of different datatransfer modes |
---|
54 | ! |
---|
55 | ! 1791 2016-03-11 10:41:25Z raasch |
---|
56 | ! Debug write-statements commented out |
---|
57 | ! |
---|
58 | ! 1786 2016-03-08 05:49:27Z raasch |
---|
59 | ! change in child-parent data transfer: parent now gets data from child |
---|
60 | ! instead that child put's it to the parent |
---|
61 | ! |
---|
62 | ! 1779 2016-03-03 08:01:28Z raasch |
---|
63 | ! kind=dp replaced by wp, |
---|
64 | ! error messages removed or changed to PALM style, dim_order removed |
---|
65 | ! array management changed from linked list to sequential loop |
---|
66 | ! |
---|
67 | ! 1766 2016-02-29 08:37:15Z raasch |
---|
68 | ! modifications to allow for using PALM's pointer version |
---|
69 | ! +new routine pmc_s_set_active_data_array |
---|
70 | ! |
---|
71 | ! 1764 2016-02-28 12:45:19Z raasch |
---|
72 | ! cpp-statement added (nesting can only be used in parallel mode) |
---|
73 | ! |
---|
74 | ! 1762 2016-02-25 12:31:13Z hellstea |
---|
75 | ! Initial revision by K. Ketelsen |
---|
76 | ! |
---|
77 | ! Description: |
---|
78 | ! ------------ |
---|
79 | ! |
---|
80 | ! Parent part of Palm Model Coupler |
---|
81 | !-------------------------------------------------------------------------------! |
---|
82 | |
---|
83 | #if defined( __parallel ) |
---|
84 | USE, INTRINSIC :: ISO_C_BINDING |
---|
85 | |
---|
86 | #if defined( __mpifh ) |
---|
87 | INCLUDE "mpif.h" |
---|
88 | #else |
---|
89 | USE MPI |
---|
90 | #endif |
---|
91 | USE kinds |
---|
92 | USE pmc_general, & |
---|
93 | ONLY: arraydef, childdef, da_namedef, da_namelen, pedef, & |
---|
94 | pmc_g_setname, pmc_max_array, pmc_max_models, pmc_sort |
---|
95 | |
---|
96 | USE pmc_handle_communicator, & |
---|
97 | ONLY: m_model_comm,m_model_rank,m_model_npes, m_to_child_comm, & |
---|
98 | m_world_rank, pmc_parent_for_child |
---|
99 | |
---|
100 | USE pmc_mpi_wrapper, & |
---|
101 | ONLY: pmc_alloc_mem, pmc_bcast, pmc_time |
---|
102 | |
---|
103 | IMPLICIT NONE |
---|
104 | |
---|
105 | PRIVATE |
---|
106 | SAVE |
---|
107 | |
---|
108 | TYPE childindexdef |
---|
109 | INTEGER :: nrpoints !< |
---|
110 | INTEGER, DIMENSION(:,:), ALLOCATABLE :: index_list_2d !< |
---|
111 | END TYPE childindexdef |
---|
112 | |
---|
113 | TYPE(childdef), DIMENSION(pmc_max_models) :: children !< |
---|
114 | TYPE(childindexdef), DIMENSION(pmc_max_models) :: indchildren !< |
---|
115 | |
---|
116 | INTEGER :: next_array_in_list = 0 !< |
---|
117 | |
---|
118 | |
---|
119 | PUBLIC pmc_parent_for_child |
---|
120 | |
---|
121 | |
---|
122 | INTERFACE pmc_parentinit |
---|
123 | MODULE PROCEDURE pmc_parentinit |
---|
124 | END INTERFACE pmc_parentinit |
---|
125 | |
---|
126 | INTERFACE pmc_s_set_2d_index_list |
---|
127 | MODULE PROCEDURE pmc_s_set_2d_index_list |
---|
128 | END INTERFACE pmc_s_set_2d_index_list |
---|
129 | |
---|
130 | INTERFACE pmc_s_clear_next_array_list |
---|
131 | MODULE PROCEDURE pmc_s_clear_next_array_list |
---|
132 | END INTERFACE pmc_s_clear_next_array_list |
---|
133 | |
---|
134 | INTERFACE pmc_s_getnextarray |
---|
135 | MODULE PROCEDURE pmc_s_getnextarray |
---|
136 | END INTERFACE pmc_s_getnextarray |
---|
137 | |
---|
138 | INTERFACE pmc_s_set_dataarray |
---|
139 | MODULE PROCEDURE pmc_s_set_dataarray_2d |
---|
140 | MODULE PROCEDURE pmc_s_set_dataarray_3d |
---|
141 | END INTERFACE pmc_s_set_dataarray |
---|
142 | |
---|
143 | INTERFACE pmc_s_setind_and_allocmem |
---|
144 | MODULE PROCEDURE pmc_s_setind_and_allocmem |
---|
145 | END INTERFACE pmc_s_setind_and_allocmem |
---|
146 | |
---|
147 | INTERFACE pmc_s_fillbuffer |
---|
148 | MODULE PROCEDURE pmc_s_fillbuffer |
---|
149 | END INTERFACE pmc_s_fillbuffer |
---|
150 | |
---|
151 | INTERFACE pmc_s_getdata_from_buffer |
---|
152 | MODULE PROCEDURE pmc_s_getdata_from_buffer |
---|
153 | END INTERFACE pmc_s_getdata_from_buffer |
---|
154 | |
---|
155 | INTERFACE pmc_s_set_active_data_array |
---|
156 | MODULE PROCEDURE pmc_s_set_active_data_array |
---|
157 | END INTERFACE pmc_s_set_active_data_array |
---|
158 | |
---|
159 | PUBLIC pmc_parentinit, pmc_s_clear_next_array_list, pmc_s_fillbuffer, & |
---|
160 | pmc_s_getdata_from_buffer, pmc_s_getnextarray, & |
---|
161 | pmc_s_setind_and_allocmem, pmc_s_set_active_data_array, & |
---|
162 | pmc_s_set_dataarray, pmc_s_set_2d_index_list |
---|
163 | |
---|
164 | CONTAINS |
---|
165 | |
---|
166 | |
---|
167 | SUBROUTINE pmc_parentinit |
---|
168 | |
---|
169 | IMPLICIT NONE |
---|
170 | |
---|
171 | INTEGER :: childid !< |
---|
172 | INTEGER :: i !< |
---|
173 | INTEGER :: j !< |
---|
174 | INTEGER :: istat !< |
---|
175 | |
---|
176 | |
---|
177 | DO i = 1, SIZE( pmc_parent_for_child )-1 |
---|
178 | |
---|
179 | childid = pmc_parent_for_child( i ) |
---|
180 | |
---|
181 | children(childid)%model_comm = m_model_comm |
---|
182 | children(childid)%inter_comm = m_to_child_comm(childid) |
---|
183 | |
---|
184 | ! |
---|
185 | !-- Get rank and size |
---|
186 | CALL MPI_COMM_RANK( children(childid)%model_comm, & |
---|
187 | children(childid)%model_rank, istat ) |
---|
188 | CALL MPI_COMM_SIZE( children(childid)%model_comm, & |
---|
189 | children(childid)%model_npes, istat ) |
---|
190 | CALL MPI_COMM_REMOTE_SIZE( children(childid)%inter_comm, & |
---|
191 | children(childid)%inter_npes, istat ) |
---|
192 | |
---|
193 | ! |
---|
194 | !-- Intra communicater is used for MPI_GET |
---|
195 | CALL MPI_INTERCOMM_MERGE( children(childid)%inter_comm, .FALSE., & |
---|
196 | children(childid)%intra_comm, istat ) |
---|
197 | CALL MPI_COMM_RANK( children(childid)%intra_comm, & |
---|
198 | children(childid)%intra_rank, istat ) |
---|
199 | |
---|
200 | ALLOCATE( children(childid)%pes(children(childid)%inter_npes)) |
---|
201 | |
---|
202 | ! |
---|
203 | !-- Allocate array of TYPE arraydef for all child PEs to store information |
---|
204 | !-- of the transfer array |
---|
205 | DO j = 1, children(childid)%inter_npes |
---|
206 | ALLOCATE( children(childid)%pes(j)%array_list(pmc_max_array) ) |
---|
207 | ENDDO |
---|
208 | |
---|
209 | CALL get_da_names_from_child (childid) |
---|
210 | |
---|
211 | ENDDO |
---|
212 | |
---|
213 | END SUBROUTINE pmc_parentinit |
---|
214 | |
---|
215 | |
---|
216 | |
---|
217 | SUBROUTINE pmc_s_set_2d_index_list( childid, index_list ) |
---|
218 | |
---|
219 | IMPLICIT NONE |
---|
220 | |
---|
221 | INTEGER, INTENT(IN) :: childid !< |
---|
222 | INTEGER, DIMENSION(:,:), INTENT(INOUT) :: index_list !< |
---|
223 | |
---|
224 | INTEGER :: ian !< |
---|
225 | INTEGER :: ic !< |
---|
226 | INTEGER :: ie !< |
---|
227 | INTEGER :: ip !< |
---|
228 | INTEGER :: is !< |
---|
229 | INTEGER :: istat !< |
---|
230 | INTEGER :: n !< |
---|
231 | |
---|
232 | |
---|
233 | IF ( m_model_rank == 0 ) THEN |
---|
234 | |
---|
235 | ! |
---|
236 | !-- Sort to ascending parent PE order |
---|
237 | CALL pmc_sort( index_list, 6 ) |
---|
238 | |
---|
239 | is = 1 |
---|
240 | DO ip = 0, m_model_npes-1 |
---|
241 | |
---|
242 | ! |
---|
243 | !-- Split into parent PEs |
---|
244 | ie = is - 1 |
---|
245 | |
---|
246 | ! |
---|
247 | !-- There may be no entry for this PE |
---|
248 | IF ( is <= SIZE( index_list,2 ) .AND. ie >= 0 ) THEN |
---|
249 | |
---|
250 | DO WHILE ( index_list(6,ie+1 ) == ip ) |
---|
251 | ie = ie + 1 |
---|
252 | IF ( ie == SIZE( index_list,2 ) ) EXIT |
---|
253 | ENDDO |
---|
254 | |
---|
255 | ian = ie - is + 1 |
---|
256 | |
---|
257 | ELSE |
---|
258 | is = -1 |
---|
259 | ie = -2 |
---|
260 | ian = 0 |
---|
261 | ENDIF |
---|
262 | |
---|
263 | ! |
---|
264 | !-- Send data to other parent PEs |
---|
265 | IF ( ip == 0 ) THEN |
---|
266 | indchildren(childid)%nrpoints = ian |
---|
267 | IF ( ian > 0) THEN |
---|
268 | ALLOCATE( indchildren(childid)%index_list_2d(6,ian) ) |
---|
269 | indchildren(childid)%index_list_2d(:,1:ian) = & |
---|
270 | index_list(:,is:ie) |
---|
271 | ENDIF |
---|
272 | ELSE |
---|
273 | CALL MPI_SEND( ian, 1, MPI_INTEGER, ip, 1000, m_model_comm, & |
---|
274 | istat ) |
---|
275 | IF ( ian > 0) THEN |
---|
276 | CALL MPI_SEND( index_list(1,is), 6*ian, MPI_INTEGER, ip, & |
---|
277 | 1001, m_model_comm, istat ) |
---|
278 | ENDIF |
---|
279 | ENDIF |
---|
280 | is = ie + 1 |
---|
281 | |
---|
282 | ENDDO |
---|
283 | |
---|
284 | ELSE |
---|
285 | |
---|
286 | CALL MPI_RECV( indchildren(childid)%nrpoints, 1, MPI_INTEGER, 0, 1000, & |
---|
287 | m_model_comm, MPI_STATUS_IGNORE, istat ) |
---|
288 | ian = indchildren(childid)%nrpoints |
---|
289 | |
---|
290 | IF ( ian > 0 ) THEN |
---|
291 | ALLOCATE( indchildren(childid)%index_list_2d(6,ian) ) |
---|
292 | CALL MPI_RECV( indchildren(childid)%index_list_2d, 6*ian, & |
---|
293 | MPI_INTEGER, 0, 1001, m_model_comm, & |
---|
294 | MPI_STATUS_IGNORE, istat) |
---|
295 | ENDIF |
---|
296 | |
---|
297 | ENDIF |
---|
298 | |
---|
299 | CALL set_pe_index_list( childid, children(childid), & |
---|
300 | indchildren(childid)%index_list_2d, & |
---|
301 | indchildren(childid)%nrpoints ) |
---|
302 | |
---|
303 | END SUBROUTINE pmc_s_set_2d_index_list |
---|
304 | |
---|
305 | |
---|
306 | |
---|
307 | SUBROUTINE pmc_s_clear_next_array_list |
---|
308 | |
---|
309 | IMPLICIT NONE |
---|
310 | |
---|
311 | next_array_in_list = 0 |
---|
312 | |
---|
313 | END SUBROUTINE pmc_s_clear_next_array_list |
---|
314 | |
---|
315 | |
---|
316 | |
---|
317 | LOGICAL FUNCTION pmc_s_getnextarray( childid, myname ) |
---|
318 | |
---|
319 | ! |
---|
320 | !-- List handling is still required to get minimal interaction with |
---|
321 | !-- pmc_interface |
---|
322 | !-- TODO: what does "still" mean? Is there a chance to change this! |
---|
323 | CHARACTER(LEN=*), INTENT(OUT) :: myname !< |
---|
324 | INTEGER(iwp), INTENT(IN) :: childid !< |
---|
325 | |
---|
326 | TYPE(arraydef), POINTER :: ar |
---|
327 | TYPE(pedef), POINTER :: ape |
---|
328 | |
---|
329 | next_array_in_list = next_array_in_list + 1 |
---|
330 | |
---|
331 | ! |
---|
332 | !-- Array names are the same on all children PEs, so take first PE to get the name |
---|
333 | ape => children(childid)%pes(1) |
---|
334 | |
---|
335 | IF ( next_array_in_list > ape%nr_arrays ) THEN |
---|
336 | |
---|
337 | ! |
---|
338 | !-- All arrays are done |
---|
339 | pmc_s_getnextarray = .FALSE. |
---|
340 | RETURN |
---|
341 | ENDIF |
---|
342 | |
---|
343 | ar => ape%array_list(next_array_in_list) |
---|
344 | myname = ar%name |
---|
345 | |
---|
346 | ! |
---|
347 | !-- Return true if legal array |
---|
348 | !-- TODO: what does this comment mean? Can there be non-legal arrays?? |
---|
349 | pmc_s_getnextarray = .TRUE. |
---|
350 | |
---|
351 | END FUNCTION pmc_s_getnextarray |
---|
352 | |
---|
353 | |
---|
354 | |
---|
355 | SUBROUTINE pmc_s_set_dataarray_2d( childid, array, array_2 ) |
---|
356 | |
---|
357 | IMPLICIT NONE |
---|
358 | |
---|
359 | INTEGER,INTENT(IN) :: childid !< |
---|
360 | |
---|
361 | REAL(wp), INTENT(IN), DIMENSION(:,:), POINTER :: array !< |
---|
362 | REAL(wp), INTENT(IN), DIMENSION(:,:), POINTER, OPTIONAL :: array_2 !< |
---|
363 | |
---|
364 | INTEGER :: nrdims !< |
---|
365 | INTEGER, DIMENSION(4) :: dims !< |
---|
366 | TYPE(C_PTR) :: array_adr !< |
---|
367 | TYPE(C_PTR) :: second_adr !< |
---|
368 | |
---|
369 | |
---|
370 | dims = 1 |
---|
371 | nrdims = 2 |
---|
372 | dims(1) = SIZE( array,1 ) |
---|
373 | dims(2) = SIZE( array,2 ) |
---|
374 | array_adr = C_LOC( array ) |
---|
375 | |
---|
376 | IF ( PRESENT( array_2 ) ) THEN |
---|
377 | second_adr = C_LOC(array_2) |
---|
378 | CALL pmc_s_setarray( childid, nrdims, dims, array_adr, & |
---|
379 | second_adr = second_adr) |
---|
380 | ELSE |
---|
381 | CALL pmc_s_setarray( childid, nrdims, dims, array_adr ) |
---|
382 | ENDIF |
---|
383 | |
---|
384 | END SUBROUTINE pmc_s_set_dataarray_2d |
---|
385 | |
---|
386 | |
---|
387 | |
---|
388 | SUBROUTINE pmc_s_set_dataarray_3d( childid, array, nz_cl, nz, array_2 ) |
---|
389 | |
---|
390 | IMPLICIT NONE |
---|
391 | |
---|
392 | INTEGER, INTENT(IN) :: childid !< |
---|
393 | INTEGER, INTENT(IN) :: nz !< |
---|
394 | INTEGER, INTENT(IN) :: nz_cl !< |
---|
395 | |
---|
396 | REAL(wp), INTENT(IN), DIMENSION(:,:,:), POINTER :: array !< |
---|
397 | REAL(wp), INTENT(IN), DIMENSION(:,:,:), POINTER, OPTIONAL :: array_2 !< |
---|
398 | |
---|
399 | INTEGER :: nrdims !< |
---|
400 | INTEGER, DIMENSION(4) :: dims !< |
---|
401 | TYPE(C_PTR) :: array_adr !< |
---|
402 | TYPE(C_PTR) :: second_adr !< |
---|
403 | |
---|
404 | ! |
---|
405 | !-- TODO: the next assignment seems to be obsolete. Please check! |
---|
406 | dims = 1 |
---|
407 | dims = 0 |
---|
408 | nrdims = 3 |
---|
409 | dims(1) = SIZE( array,1 ) |
---|
410 | dims(2) = SIZE( array,2 ) |
---|
411 | dims(3) = SIZE( array,3 ) |
---|
412 | dims(4) = nz_cl+dims(1)-nz ! works for first dimension 1:nz and 0:nz+1 |
---|
413 | |
---|
414 | array_adr = C_LOC(array) |
---|
415 | |
---|
416 | ! |
---|
417 | !-- In PALM's pointer version, two indices have to be stored internally. |
---|
418 | !-- The active address of the data array is set in swap_timelevel. |
---|
419 | IF ( PRESENT( array_2 ) ) THEN |
---|
420 | second_adr = C_LOC( array_2 ) |
---|
421 | CALL pmc_s_setarray( childid, nrdims, dims, array_adr, & |
---|
422 | second_adr = second_adr) |
---|
423 | ELSE |
---|
424 | CALL pmc_s_setarray( childid, nrdims, dims, array_adr ) |
---|
425 | ENDIF |
---|
426 | |
---|
427 | END SUBROUTINE pmc_s_set_dataarray_3d |
---|
428 | |
---|
429 | |
---|
430 | |
---|
431 | SUBROUTINE pmc_s_setind_and_allocmem( childid ) |
---|
432 | |
---|
433 | USE control_parameters, & |
---|
434 | ONLY: message_string |
---|
435 | |
---|
436 | IMPLICIT NONE |
---|
437 | |
---|
438 | ! |
---|
439 | !-- Naming convention for appendices: _pc -> parent to child transfer |
---|
440 | !-- _cp -> child to parent transfer |
---|
441 | !-- send -> parent to child transfer |
---|
442 | !-- recv -> child to parent transfer |
---|
443 | INTEGER, INTENT(IN) :: childid !< |
---|
444 | |
---|
445 | INTEGER :: arlen !< |
---|
446 | INTEGER :: i !< |
---|
447 | INTEGER :: ierr !< |
---|
448 | INTEGER :: istat !< |
---|
449 | INTEGER :: j !< |
---|
450 | INTEGER :: myindex !< |
---|
451 | INTEGER :: rcount !< count MPI requests |
---|
452 | INTEGER :: tag !< |
---|
453 | |
---|
454 | INTEGER(idp) :: bufsize !< size of MPI data window |
---|
455 | INTEGER(KIND=MPI_ADDRESS_KIND) :: winsize !< |
---|
456 | |
---|
457 | INTEGER, DIMENSION(1024) :: req !< |
---|
458 | |
---|
459 | TYPE(C_PTR) :: base_ptr !< |
---|
460 | TYPE(pedef), POINTER :: ape !< |
---|
461 | TYPE(arraydef), POINTER :: ar !< |
---|
462 | |
---|
463 | REAL(wp),DIMENSION(:), POINTER, SAVE :: base_array_pc !< base array for parent to child transfer |
---|
464 | REAL(wp),DIMENSION(:), POINTER, SAVE :: base_array_cp !< base array for child to parent transfer |
---|
465 | |
---|
466 | ! |
---|
467 | !-- Parent to child direction |
---|
468 | myindex = 1 |
---|
469 | rcount = 0 |
---|
470 | bufsize = 8 |
---|
471 | |
---|
472 | ! |
---|
473 | !-- First stride: compute size and set index |
---|
474 | DO i = 1, children(childid)%inter_npes |
---|
475 | |
---|
476 | ape => children(childid)%pes(i) |
---|
477 | tag = 200 |
---|
478 | |
---|
479 | DO j = 1, ape%nr_arrays |
---|
480 | |
---|
481 | ar => ape%array_list(j) |
---|
482 | IF ( ar%nrdims == 2 ) THEN |
---|
483 | arlen = ape%nrele |
---|
484 | ELSEIF ( ar%nrdims == 3 ) THEN |
---|
485 | arlen = ape%nrele * ar%a_dim(4) |
---|
486 | ELSE |
---|
487 | arlen = -1 |
---|
488 | ENDIF |
---|
489 | ar%sendindex = myindex |
---|
490 | |
---|
491 | tag = tag + 1 |
---|
492 | rcount = rcount + 1 |
---|
493 | CALL MPI_ISEND( myindex, 1, MPI_INTEGER, i-1, tag, & |
---|
494 | children(childid)%inter_comm, req(rcount), ierr ) |
---|
495 | |
---|
496 | ! |
---|
497 | !-- Maximum of 1024 outstanding requests |
---|
498 | !-- TODO: what does this limit mean? |
---|
499 | IF ( rcount == 1024 ) THEN |
---|
500 | CALL MPI_WAITALL( rcount, req, MPI_STATUSES_IGNORE, ierr ) |
---|
501 | rcount = 0 |
---|
502 | ENDIF |
---|
503 | |
---|
504 | myindex = myindex + arlen |
---|
505 | bufsize = bufsize + arlen |
---|
506 | ar%sendsize = arlen |
---|
507 | |
---|
508 | ENDDO |
---|
509 | |
---|
510 | IF ( rcount > 0 ) THEN |
---|
511 | CALL MPI_WAITALL( rcount, req, MPI_STATUSES_IGNORE, ierr ) |
---|
512 | ENDIF |
---|
513 | |
---|
514 | ENDDO |
---|
515 | |
---|
516 | ! |
---|
517 | !-- Create RMA (One Sided Communication) window for data buffer parent to |
---|
518 | !-- child transfer. |
---|
519 | !-- The buffer of MPI_GET (counterpart of transfer) can be PE-local, i.e. |
---|
520 | !-- it can but must not be part of the MPI RMA window. Only one RMA window is |
---|
521 | !-- required to prepare the data for |
---|
522 | !-- parent -> child transfer on the parent side |
---|
523 | !-- and for |
---|
524 | !-- child -> parent transfer on the child side |
---|
525 | CALL pmc_alloc_mem( base_array_pc, bufsize ) |
---|
526 | children(childid)%totalbuffersize = bufsize * wp |
---|
527 | |
---|
528 | winsize = bufsize * wp |
---|
529 | CALL MPI_WIN_CREATE( base_array_pc, winsize, wp, MPI_INFO_NULL, & |
---|
530 | children(childid)%intra_comm, & |
---|
531 | children(childid)%win_parent_child, ierr ) |
---|
532 | |
---|
533 | ! |
---|
534 | !-- Open window to set data |
---|
535 | CALL MPI_WIN_FENCE( 0, children(childid)%win_parent_child, ierr ) |
---|
536 | |
---|
537 | ! |
---|
538 | !-- Second stride: set buffer pointer |
---|
539 | DO i = 1, children(childid)%inter_npes |
---|
540 | |
---|
541 | ape => children(childid)%pes(i) |
---|
542 | |
---|
543 | DO j = 1, ape%nr_arrays |
---|
544 | |
---|
545 | ar => ape%array_list(j) |
---|
546 | ar%sendbuf = C_LOC( base_array_pc(ar%sendindex) ) |
---|
547 | |
---|
548 | IF ( ar%sendindex + ar%sendsize > bufsize ) THEN |
---|
549 | WRITE( message_string, '(a,i4,4i7,1x,a)' ) & |
---|
550 | 'parent buffer too small ',i, & |
---|
551 | ar%sendindex,ar%sendsize,ar%sendindex+ar%sendsize, & |
---|
552 | bufsize,trim(ar%name) |
---|
553 | CALL message( 'pmc_s_setind_and_allocmem', 'PA0429', 3, 2, 0, 6, 0 ) |
---|
554 | ENDIF |
---|
555 | ENDDO |
---|
556 | ENDDO |
---|
557 | |
---|
558 | ! |
---|
559 | !-- Child to parent direction |
---|
560 | bufsize = 8 |
---|
561 | |
---|
562 | ! |
---|
563 | !-- First stride: compute size and set index |
---|
564 | DO i = 1, children(childid)%inter_npes |
---|
565 | |
---|
566 | ape => children(childid)%pes(i) |
---|
567 | tag = 300 |
---|
568 | |
---|
569 | DO j = 1, ape%nr_arrays |
---|
570 | |
---|
571 | ar => ape%array_list(j) |
---|
572 | |
---|
573 | ! |
---|
574 | !-- Receive index from child |
---|
575 | tag = tag + 1 |
---|
576 | CALL MPI_RECV( myindex, 1, MPI_INTEGER, i-1, tag, & |
---|
577 | children(childid)%inter_comm, MPI_STATUS_IGNORE, ierr ) |
---|
578 | |
---|
579 | IF ( ar%nrdims == 3 ) THEN |
---|
580 | bufsize = MAX( bufsize, ape%nrele * ar%a_dim(4) ) |
---|
581 | ELSE |
---|
582 | bufsize = MAX( bufsize, ape%nrele ) |
---|
583 | ENDIF |
---|
584 | ar%recvindex = myindex |
---|
585 | |
---|
586 | ENDDO |
---|
587 | |
---|
588 | ENDDO |
---|
589 | |
---|
590 | ! |
---|
591 | !-- Create RMA (one sided communication) data buffer. |
---|
592 | !-- The buffer for MPI_GET can be PE local, i.e. it can but must not be part of |
---|
593 | !-- the MPI RMA window |
---|
594 | CALL pmc_alloc_mem( base_array_cp, bufsize, base_ptr ) |
---|
595 | children(childid)%totalbuffersize = bufsize * wp |
---|
596 | |
---|
597 | CALL MPI_BARRIER( children(childid)%intra_comm, ierr ) |
---|
598 | |
---|
599 | ! |
---|
600 | !-- Second stride: set buffer pointer |
---|
601 | DO i = 1, children(childid)%inter_npes |
---|
602 | |
---|
603 | ape => children(childid)%pes(i) |
---|
604 | |
---|
605 | DO j = 1, ape%nr_arrays |
---|
606 | ar => ape%array_list(j) |
---|
607 | ar%recvbuf = base_ptr |
---|
608 | ENDDO |
---|
609 | |
---|
610 | ENDDO |
---|
611 | |
---|
612 | END SUBROUTINE pmc_s_setind_and_allocmem |
---|
613 | |
---|
614 | |
---|
615 | |
---|
616 | SUBROUTINE pmc_s_fillbuffer( childid, waittime ) |
---|
617 | |
---|
618 | IMPLICIT NONE |
---|
619 | |
---|
620 | INTEGER, INTENT(IN) :: childid !< |
---|
621 | |
---|
622 | REAL(wp), INTENT(OUT), OPTIONAL :: waittime !< |
---|
623 | |
---|
624 | INTEGER :: ierr !< |
---|
625 | INTEGER :: ij !< |
---|
626 | INTEGER :: ip !< |
---|
627 | INTEGER :: istat !< |
---|
628 | INTEGER :: j !< |
---|
629 | INTEGER :: myindex !< |
---|
630 | |
---|
631 | INTEGER, DIMENSION(1) :: buf_shape |
---|
632 | |
---|
633 | REAL(wp) :: t1 !< |
---|
634 | REAL(wp) :: t2 !< |
---|
635 | REAL(wp), POINTER, DIMENSION(:) :: buf !< |
---|
636 | REAL(wp), POINTER, DIMENSION(:,:) :: data_2d !< |
---|
637 | REAL(wp), POINTER, DIMENSION(:,:,:) :: data_3d !< |
---|
638 | |
---|
639 | TYPE(pedef), POINTER :: ape !< |
---|
640 | TYPE(arraydef), POINTER :: ar !< |
---|
641 | |
---|
642 | ! |
---|
643 | !-- Synchronization of the model is done in pmci_synchronize. |
---|
644 | !-- Therefor the RMA window can be filled without |
---|
645 | !-- sychronization at this point and a barrier is not necessary. |
---|
646 | !-- Please note that waittime has to be set in pmc_s_fillbuffer AND |
---|
647 | !-- pmc_c_getbuffer |
---|
648 | IF ( PRESENT( waittime) ) THEN |
---|
649 | t1 = pmc_time() |
---|
650 | CALL MPI_BARRIER( children(childid)%intra_comm, ierr ) |
---|
651 | t2 = pmc_time() |
---|
652 | waittime = t2- t1 |
---|
653 | ENDIF |
---|
654 | |
---|
655 | DO ip = 1, children(childid)%inter_npes |
---|
656 | |
---|
657 | ape => children(childid)%pes(ip) |
---|
658 | |
---|
659 | DO j = 1, ape%nr_arrays |
---|
660 | |
---|
661 | ar => ape%array_list(j) |
---|
662 | myindex = 1 |
---|
663 | |
---|
664 | IF ( ar%nrdims == 2 ) THEN |
---|
665 | |
---|
666 | buf_shape(1) = ape%nrele |
---|
667 | CALL C_F_POINTER( ar%sendbuf, buf, buf_shape ) |
---|
668 | CALL C_F_POINTER( ar%data, data_2d, ar%a_dim(1:2) ) |
---|
669 | DO ij = 1, ape%nrele |
---|
670 | buf(myindex) = data_2d(ape%locind(ij)%j,ape%locind(ij)%i) |
---|
671 | myindex = myindex + 1 |
---|
672 | ENDDO |
---|
673 | |
---|
674 | ELSEIF ( ar%nrdims == 3 ) THEN |
---|
675 | |
---|
676 | buf_shape(1) = ape%nrele*ar%a_dim(4) |
---|
677 | CALL C_F_POINTER( ar%sendbuf, buf, buf_shape ) |
---|
678 | CALL C_F_POINTER( ar%data, data_3d, ar%a_dim(1:3) ) |
---|
679 | DO ij = 1, ape%nrele |
---|
680 | buf(myindex:myindex+ar%a_dim(4)-1) = & |
---|
681 | data_3d(1:ar%a_dim(4),ape%locind(ij)%j,ape%locind(ij)%i) |
---|
682 | myindex = myindex + ar%a_dim(4) |
---|
683 | ENDDO |
---|
684 | |
---|
685 | ENDIF |
---|
686 | |
---|
687 | ENDDO |
---|
688 | |
---|
689 | ENDDO |
---|
690 | |
---|
691 | ! |
---|
692 | !-- Buffer is filled |
---|
693 | CALL MPI_BARRIER( children(childid)%intra_comm, ierr ) |
---|
694 | |
---|
695 | END SUBROUTINE pmc_s_fillbuffer |
---|
696 | |
---|
697 | |
---|
698 | |
---|
699 | SUBROUTINE pmc_s_getdata_from_buffer( childid, waittime ) |
---|
700 | |
---|
701 | IMPLICIT NONE |
---|
702 | |
---|
703 | INTEGER, INTENT(IN) :: childid !< |
---|
704 | REAL(wp), INTENT(OUT), OPTIONAL :: waittime !< |
---|
705 | |
---|
706 | INTEGER :: ierr !< |
---|
707 | INTEGER :: ij !< |
---|
708 | INTEGER :: ip !< |
---|
709 | INTEGER :: istat !< |
---|
710 | INTEGER :: j !< |
---|
711 | INTEGER :: myindex !< |
---|
712 | INTEGER :: nr !< |
---|
713 | INTEGER :: target_pe !< |
---|
714 | INTEGER(kind=MPI_ADDRESS_KIND) :: target_disp !< |
---|
715 | |
---|
716 | INTEGER, DIMENSION(1) :: buf_shape !< |
---|
717 | |
---|
718 | REAL(wp) :: t1 !< |
---|
719 | REAL(wp) :: t2 !< |
---|
720 | REAL(wp), POINTER, DIMENSION(:) :: buf !< |
---|
721 | REAL(wp), POINTER, DIMENSION(:,:) :: data_2d !< |
---|
722 | REAL(wp), POINTER, DIMENSION(:,:,:) :: data_3d !< |
---|
723 | |
---|
724 | TYPE(pedef), POINTER :: ape !< |
---|
725 | TYPE(arraydef), POINTER :: ar !< |
---|
726 | |
---|
727 | |
---|
728 | t1 = pmc_time() |
---|
729 | |
---|
730 | ! |
---|
731 | !-- Wait for child to fill buffer |
---|
732 | CALL MPI_BARRIER( children(childid)%intra_comm, ierr ) |
---|
733 | t2 = pmc_time() - t1 |
---|
734 | IF ( PRESENT( waittime ) ) waittime = t2 |
---|
735 | |
---|
736 | ! |
---|
737 | !-- TODO: check next statement |
---|
738 | !-- Fence might do it, test later |
---|
739 | !-- CALL MPI_WIN_FENCE( 0, children(childid)%win_parent_child, ierr) |
---|
740 | CALL MPI_BARRIER( children(childid)%intra_comm, ierr ) |
---|
741 | |
---|
742 | DO ip = 1, children(childid)%inter_npes |
---|
743 | |
---|
744 | ape => children(childid)%pes(ip) |
---|
745 | |
---|
746 | DO j = 1, ape%nr_arrays |
---|
747 | |
---|
748 | ar => ape%array_list(j) |
---|
749 | |
---|
750 | IF ( ar%recvindex < 0 ) CYCLE |
---|
751 | |
---|
752 | IF ( ar%nrdims == 2 ) THEN |
---|
753 | nr = ape%nrele |
---|
754 | ELSEIF ( ar%nrdims == 3 ) THEN |
---|
755 | nr = ape%nrele * ar%a_dim(4) |
---|
756 | ENDIF |
---|
757 | |
---|
758 | buf_shape(1) = nr |
---|
759 | CALL C_F_POINTER( ar%recvbuf, buf, buf_shape ) |
---|
760 | |
---|
761 | ! |
---|
762 | !-- MPI passive target RMA |
---|
763 | IF ( nr > 0 ) THEN |
---|
764 | target_disp = ar%recvindex - 1 |
---|
765 | |
---|
766 | ! |
---|
767 | !-- Child PEs are located behind parent PEs |
---|
768 | target_pe = ip - 1 + m_model_npes |
---|
769 | CALL MPI_WIN_LOCK( MPI_LOCK_SHARED, target_pe, 0, & |
---|
770 | children(childid)%win_parent_child, ierr ) |
---|
771 | CALL MPI_GET( buf, nr, MPI_REAL, target_pe, target_disp, nr, & |
---|
772 | MPI_REAL, children(childid)%win_parent_child, ierr ) |
---|
773 | CALL MPI_WIN_UNLOCK( target_pe, & |
---|
774 | children(childid)%win_parent_child, ierr ) |
---|
775 | ENDIF |
---|
776 | |
---|
777 | myindex = 1 |
---|
778 | IF ( ar%nrdims == 2 ) THEN |
---|
779 | |
---|
780 | CALL C_F_POINTER( ar%data, data_2d, ar%a_dim(1:2) ) |
---|
781 | DO ij = 1, ape%nrele |
---|
782 | data_2d(ape%locind(ij)%j,ape%locind(ij)%i) = buf(myindex) |
---|
783 | myindex = myindex + 1 |
---|
784 | ENDDO |
---|
785 | |
---|
786 | ELSEIF ( ar%nrdims == 3 ) THEN |
---|
787 | |
---|
788 | CALL C_F_POINTER( ar%data, data_3d, ar%a_dim(1:3)) |
---|
789 | DO ij = 1, ape%nrele |
---|
790 | data_3d(1:ar%a_dim(4),ape%locind(ij)%j,ape%locind(ij)%i) = & |
---|
791 | buf(myindex:myindex+ar%a_dim(4)-1) |
---|
792 | myindex = myindex + ar%a_dim(4) |
---|
793 | ENDDO |
---|
794 | |
---|
795 | ENDIF |
---|
796 | |
---|
797 | ENDDO |
---|
798 | |
---|
799 | ENDDO |
---|
800 | |
---|
801 | END SUBROUTINE pmc_s_getdata_from_buffer |
---|
802 | |
---|
803 | |
---|
804 | |
---|
805 | SUBROUTINE get_da_names_from_child( childid ) |
---|
806 | |
---|
807 | ! |
---|
808 | !-- Get data array description and name from child |
---|
809 | IMPLICIT NONE |
---|
810 | |
---|
811 | INTEGER, INTENT(IN) :: childid !< |
---|
812 | |
---|
813 | TYPE(da_namedef) :: myname !< |
---|
814 | |
---|
815 | DO |
---|
816 | CALL pmc_bcast( myname%couple_index, 0, comm=m_to_child_comm(childid) ) |
---|
817 | IF ( myname%couple_index == -1 ) EXIT |
---|
818 | CALL pmc_bcast( myname%parentdesc, 0, comm=m_to_child_comm(childid) ) |
---|
819 | CALL pmc_bcast( myname%nameonparent, 0, comm=m_to_child_comm(childid) ) |
---|
820 | CALL pmc_bcast( myname%childdesc, 0, comm=m_to_child_comm(childid) ) |
---|
821 | CALL pmc_bcast( myname%nameonchild, 0, comm=m_to_child_comm(childid) ) |
---|
822 | |
---|
823 | CALL pmc_g_setname( children(childid), myname%couple_index, & |
---|
824 | myname%nameonparent ) |
---|
825 | ENDDO |
---|
826 | |
---|
827 | END SUBROUTINE get_da_names_from_child |
---|
828 | |
---|
829 | |
---|
830 | |
---|
831 | SUBROUTINE pmc_s_setarray(childid, nrdims, dims, array_adr, second_adr ) |
---|
832 | |
---|
833 | ! |
---|
834 | !-- Set array for child inter PE 0 |
---|
835 | IMPLICIT NONE |
---|
836 | |
---|
837 | INTEGER, INTENT(IN) :: childid !< |
---|
838 | INTEGER, INTENT(IN) :: nrdims !< |
---|
839 | INTEGER, INTENT(IN), DIMENSION(:) :: dims !< |
---|
840 | |
---|
841 | TYPE(C_PTR), INTENT(IN) :: array_adr !< |
---|
842 | TYPE(C_PTR), INTENT(IN), OPTIONAL :: second_adr !< |
---|
843 | |
---|
844 | INTEGER :: i !< local counter |
---|
845 | |
---|
846 | TYPE(pedef), POINTER :: ape !< |
---|
847 | TYPE(arraydef), POINTER :: ar !< |
---|
848 | |
---|
849 | |
---|
850 | DO i = 1, children(childid)%inter_npes |
---|
851 | |
---|
852 | ape => children(childid)%pes(i) |
---|
853 | ar => ape%array_list(next_array_in_list) |
---|
854 | ar%nrdims = nrdims |
---|
855 | ar%a_dim = dims |
---|
856 | ar%data = array_adr |
---|
857 | |
---|
858 | IF ( PRESENT( second_adr ) ) THEN |
---|
859 | ar%po_data(1) = array_adr |
---|
860 | ar%po_data(2) = second_adr |
---|
861 | ELSE |
---|
862 | ar%po_data(1) = C_NULL_PTR |
---|
863 | ar%po_data(2) = C_NULL_PTR |
---|
864 | ENDIF |
---|
865 | |
---|
866 | ENDDO |
---|
867 | |
---|
868 | END SUBROUTINE pmc_s_setarray |
---|
869 | |
---|
870 | |
---|
871 | |
---|
872 | SUBROUTINE pmc_s_set_active_data_array( childid, iactive ) |
---|
873 | |
---|
874 | IMPLICIT NONE |
---|
875 | |
---|
876 | INTEGER, INTENT(IN) :: childid !< |
---|
877 | INTEGER, INTENT(IN) :: iactive !< |
---|
878 | |
---|
879 | INTEGER :: i !< |
---|
880 | INTEGER :: ip !< |
---|
881 | INTEGER :: j !< |
---|
882 | |
---|
883 | TYPE(pedef), POINTER :: ape !< |
---|
884 | TYPE(arraydef), POINTER :: ar !< |
---|
885 | |
---|
886 | DO ip = 1, children(childid)%inter_npes |
---|
887 | |
---|
888 | ape => children(childid)%pes(ip) |
---|
889 | |
---|
890 | DO j = 1, ape%nr_arrays |
---|
891 | |
---|
892 | ar => ape%array_list(j) |
---|
893 | IF ( iactive == 1 .OR. iactive == 2 ) THEN |
---|
894 | ar%data = ar%po_data(iactive) |
---|
895 | ENDIF |
---|
896 | |
---|
897 | ENDDO |
---|
898 | |
---|
899 | ENDDO |
---|
900 | |
---|
901 | END SUBROUTINE pmc_s_set_active_data_array |
---|
902 | |
---|
903 | |
---|
904 | |
---|
905 | SUBROUTINE set_pe_index_list( childid, mychild, index_list, nrp ) |
---|
906 | |
---|
907 | IMPLICIT NONE |
---|
908 | |
---|
909 | INTEGER, INTENT(IN) :: childid !< |
---|
910 | INTEGER, INTENT(IN), DIMENSION(:,:) :: index_list !< |
---|
911 | INTEGER, INTENT(IN) :: nrp !< |
---|
912 | |
---|
913 | TYPE(childdef), INTENT(INOUT) :: mychild !< |
---|
914 | |
---|
915 | INTEGER :: i !< |
---|
916 | INTEGER :: ierr !< |
---|
917 | INTEGER :: ind !< |
---|
918 | INTEGER :: indwin !< |
---|
919 | INTEGER :: indwin2 !< |
---|
920 | INTEGER :: i2 !< |
---|
921 | INTEGER :: j !< |
---|
922 | INTEGER :: rempe !< |
---|
923 | INTEGER(KIND=MPI_ADDRESS_KIND) :: winsize !< |
---|
924 | |
---|
925 | INTEGER, DIMENSION(mychild%inter_npes) :: remind !< |
---|
926 | |
---|
927 | INTEGER, DIMENSION(:), POINTER :: remindw !< |
---|
928 | INTEGER, DIMENSION(:), POINTER :: rldef !< |
---|
929 | |
---|
930 | TYPE(pedef), POINTER :: ape !< |
---|
931 | |
---|
932 | ! |
---|
933 | !-- First, count entries for every remote child PE |
---|
934 | DO i = 1, mychild%inter_npes |
---|
935 | ape => mychild%pes(i) |
---|
936 | ape%nrele = 0 |
---|
937 | ENDDO |
---|
938 | |
---|
939 | ! |
---|
940 | !-- Loop over number of coarse grid cells |
---|
941 | DO j = 1, nrp |
---|
942 | rempe = index_list(5,j) + 1 ! PE number on remote PE |
---|
943 | ape => mychild%pes(rempe) |
---|
944 | ape%nrele = ape%nrele + 1 ! Increment number of elements for this child PE |
---|
945 | ENDDO |
---|
946 | |
---|
947 | DO i = 1, mychild%inter_npes |
---|
948 | ape => mychild%pes(i) |
---|
949 | ALLOCATE( ape%locind(ape%nrele) ) |
---|
950 | ENDDO |
---|
951 | |
---|
952 | remind = 0 |
---|
953 | |
---|
954 | ! |
---|
955 | !-- Second, create lists |
---|
956 | !-- Loop over number of coarse grid cells |
---|
957 | DO j = 1, nrp |
---|
958 | rempe = index_list(5,j) + 1 |
---|
959 | ape => mychild%pes(rempe) |
---|
960 | remind(rempe) = remind(rempe)+1 |
---|
961 | ind = remind(rempe) |
---|
962 | ape%locind(ind)%i = index_list(1,j) |
---|
963 | ape%locind(ind)%j = index_list(2,j) |
---|
964 | ENDDO |
---|
965 | |
---|
966 | ! |
---|
967 | !-- Prepare number of elements for children PEs |
---|
968 | CALL pmc_alloc_mem( rldef, mychild%inter_npes*2 ) |
---|
969 | |
---|
970 | ! |
---|
971 | !-- Number of child PEs * size of INTEGER (i just arbitrary INTEGER) |
---|
972 | winsize = mychild%inter_npes*c_sizeof(i)*2 |
---|
973 | |
---|
974 | CALL MPI_WIN_CREATE( rldef, winsize, iwp, MPI_INFO_NULL, & |
---|
975 | mychild%intra_comm, indwin, ierr ) |
---|
976 | |
---|
977 | ! |
---|
978 | !-- Open window to set data |
---|
979 | CALL MPI_WIN_FENCE( 0, indwin, ierr ) |
---|
980 | |
---|
981 | rldef(1) = 0 ! index on remote PE 0 |
---|
982 | rldef(2) = remind(1) ! number of elements on remote PE 0 |
---|
983 | |
---|
984 | ! |
---|
985 | !-- Reserve buffer for index array |
---|
986 | DO i = 2, mychild%inter_npes |
---|
987 | i2 = (i-1) * 2 + 1 |
---|
988 | rldef(i2) = rldef(i2-2) + rldef(i2-1) * 2 ! index on remote PE |
---|
989 | rldef(i2+1) = remind(i) ! number of elements on remote PE |
---|
990 | ENDDO |
---|
991 | |
---|
992 | ! |
---|
993 | !-- Close window to allow child to access data |
---|
994 | CALL MPI_WIN_FENCE( 0, indwin, ierr ) |
---|
995 | |
---|
996 | ! |
---|
997 | !-- Child has retrieved data |
---|
998 | CALL MPI_WIN_FENCE( 0, indwin, ierr ) |
---|
999 | |
---|
1000 | i2 = 2 * mychild%inter_npes - 1 |
---|
1001 | winsize = ( rldef(i2) + rldef(i2+1) ) * 2 |
---|
1002 | |
---|
1003 | ! |
---|
1004 | !-- Make sure, MPI_ALLOC_MEM works |
---|
1005 | winsize = MAX( winsize, 1 ) |
---|
1006 | |
---|
1007 | CALL pmc_alloc_mem( remindw, INT( winsize ) ) |
---|
1008 | |
---|
1009 | CALL MPI_BARRIER( m_model_comm, ierr ) |
---|
1010 | CALL MPI_WIN_CREATE( remindw, winsize*c_sizeof(i), iwp, MPI_INFO_NULL, & |
---|
1011 | mychild%intra_comm, indwin2, ierr ) |
---|
1012 | ! |
---|
1013 | !-- Open window to set data |
---|
1014 | CALL MPI_WIN_FENCE( 0, indwin2, ierr ) |
---|
1015 | |
---|
1016 | ! |
---|
1017 | !-- Create the 2D index list |
---|
1018 | DO j = 1, nrp |
---|
1019 | rempe = index_list(5,j) + 1 ! PE number on remote PE |
---|
1020 | ape => mychild%pes(rempe) |
---|
1021 | i2 = rempe * 2 - 1 |
---|
1022 | ind = rldef(i2) + 1 |
---|
1023 | remindw(ind) = index_list(3,j) |
---|
1024 | remindw(ind+1) = index_list(4,j) |
---|
1025 | rldef(i2) = rldef(i2)+2 |
---|
1026 | ENDDO |
---|
1027 | |
---|
1028 | ! |
---|
1029 | !-- All data are set |
---|
1030 | CALL MPI_WIN_FENCE( 0, indwin2, ierr ) |
---|
1031 | |
---|
1032 | ! |
---|
1033 | !-- Don't know why, but this barrier is necessary before windows can be freed |
---|
1034 | !-- TODO: find out why this is required |
---|
1035 | CALL MPI_BARRIER( mychild%intra_comm, ierr ) |
---|
1036 | |
---|
1037 | CALL MPI_WIN_FREE( indwin, ierr ) |
---|
1038 | CALL MPI_WIN_FREE( indwin2, ierr ) |
---|
1039 | |
---|
1040 | ! |
---|
1041 | !-- TODO: check if the following idea needs to be done |
---|
1042 | !-- Sollte funktionieren, Problem mit MPI implementation |
---|
1043 | !-- https://www.lrz.de/services/software/parallel/mpi/onesided |
---|
1044 | !-- CALL MPI_Free_mem (remindw, ierr) |
---|
1045 | |
---|
1046 | END SUBROUTINE set_pe_index_list |
---|
1047 | |
---|
1048 | #endif |
---|
1049 | END MODULE pmc_parent |
---|