source: palm/trunk/SOURCE/pmc_child_mod.f90 @ 2965

Last change on this file since 2965 was 2841, checked in by knoop, 7 years ago

Bugfix: wrong placement of include 'mpif.h' corrected

  • Property svn:keywords set to Id
File size: 29.1 KB
RevLine 
[1933]1MODULE pmc_child
[1762]2
[2000]3!------------------------------------------------------------------------------!
[2696]4! This file is part of the PALM model system.
[1762]5!
[2000]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.
[1762]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!
[2718]18! Copyright 1997-2018 Leibniz Universitaet Hannover
[2000]19!------------------------------------------------------------------------------!
[1762]20!
21! Current revisions:
22! ------------------
[1834]23!
[2001]24!
[1762]25! Former revisions:
26! -----------------
27! $Id: pmc_child_mod.f90 2841 2018-02-27 15:02:57Z scharf $
[2841]28! Bugfix: wrong placement of include 'mpif.h' corrected
29!
30! 2809 2018-02-15 09:55:58Z schwenkel
[2809]31! Bugfix for gfortran: Replace the function C_SIZEOF with STORAGE_SIZE
32!
33! 2801 2018-02-14 16:01:55Z thiele
[2801]34! Introduce particle transfer in nested models.
35!
36! 2718 2018-01-02 08:49:38Z maronga
[2716]37! Corrected "Former revisions" section
38!
39! 2696 2017-12-14 17:12:51Z kanani
40! Change in file header (GPL part)
41!
42! 2599 2017-11-01 13:18:45Z hellstea
[2599]43! Some cleanup and commenting improvements only.
44!
45! 2101 2017-01-05 16:42:31Z suehring
[1762]46!
[2001]47! 2000 2016-08-20 18:09:15Z knoop
48! Forced header and separation lines into 80 columns
49!
[1933]50! 1897 2016-05-03 08:10:23Z raasch
51! Module renamed. Code clean up. The words server/client changed to parent/child.
52!
[1897]53! 1896 2016-05-03 08:06:41Z raasch
54! re-formatted to match PALM style
55!
[1851]56! 1850 2016-04-08 13:29:27Z maronga
57! Module renamed
58!
59!
[1834]60! 1833 2016-04-07 14:23:03Z raasch
61! gfortran requires pointer attributes for some array declarations,
62! long line wrapped
63!
[1809]64! 1808 2016-04-05 19:44:00Z raasch
65! MPI module used by default on all machines
66!
[1798]67! 1797 2016-03-21 16:50:28Z raasch
68! introduction of different datatransfer modes
69!
[1792]70! 1791 2016-03-11 10:41:25Z raasch
71! Debug write-statement commented out
72!
[1787]73! 1786 2016-03-08 05:49:27Z raasch
[1933]74! change in child-parent data transfer: parent now gets data from child
75! instead of that child puts it to the parent
[1787]76!
[1784]77! 1783 2016-03-06 18:36:17Z raasch
78! Bugfix: wrong data-type in MPI_WIN_CREATE replaced
79!
[1780]80! 1779 2016-03-03 08:01:28Z raasch
81! kind=dp replaced by wp, dim_order removed
82! array management changed from linked list to sequential loop
83!
[1765]84! 1764 2016-02-28 12:45:19Z raasch
85! cpp-statement added (nesting can only be used in parallel mode),
86! all kinds given in PALM style
87!
[1763]88! 1762 2016-02-25 12:31:13Z hellstea
89! Initial revision by K. Ketelsen
[1762]90!
91! Description:
92! ------------
93!
[1933]94! Child part of Palm Model Coupler
[2801]95!------------------------------------------------------------------------------!
[1762]96
[1764]97#if defined( __parallel )
[1762]98
[1896]99    USE, INTRINSIC ::  iso_c_binding
[1762]100
[2841]101#if !defined( __mpifh )
[1764]102    USE MPI
103#endif
[1896]104
105    USE kinds
[2841]106
[2801]107    USE pmc_general,                                                           &
108        ONLY:  arraydef, childdef, da_desclen, da_namedef, da_namelen, pedef,  &
[1896]109               pmc_da_name_err,  pmc_g_setname, pmc_max_array, pmc_status_ok
110
[2801]111    USE pmc_handle_communicator,                                               &
[1933]112        ONLY:  m_model_comm, m_model_npes, m_model_rank, m_to_parent_comm
[1896]113
[2801]114    USE pmc_mpi_wrapper,                                                       &
[1933]115        ONLY:  pmc_alloc_mem, pmc_bcast, pmc_inter_bcast, pmc_time
[1896]116
117    IMPLICIT NONE
118
[2841]119#if defined( __mpifh )
120    INCLUDE "mpif.h"
121#endif
122
[1762]123    PRIVATE
124    SAVE
125
[2801]126    TYPE(childdef), PUBLIC ::  me   !<
[1762]127
[2801]128    INTEGER(iwp) ::  myindex = 0         !< counter and unique number for data arrays
129    INTEGER(iwp) ::  next_array_in_list = 0   !<
[1762]130
131
[1933]132    INTERFACE pmc_childinit
133        MODULE PROCEDURE pmc_childinit
134    END INTERFACE pmc_childinit
[1762]135
[1896]136    INTERFACE pmc_c_clear_next_array_list
137        MODULE PROCEDURE pmc_c_clear_next_array_list
138    END INTERFACE pmc_c_clear_next_array_list
[1762]139
[1896]140    INTERFACE pmc_c_getbuffer
141        MODULE PROCEDURE pmc_c_getbuffer
142    END INTERFACE pmc_c_getbuffer
[1762]143
[1896]144    INTERFACE pmc_c_getnextarray
145        MODULE PROCEDURE pmc_c_getnextarray
146    END INTERFACE pmc_c_getnextarray
[1779]147
[1896]148    INTERFACE pmc_c_get_2d_index_list
149        MODULE PROCEDURE pmc_c_get_2d_index_list
150    END INTERFACE pmc_c_get_2d_index_list
[1762]151
[1896]152    INTERFACE pmc_c_putbuffer
153        MODULE PROCEDURE pmc_c_putbuffer
154    END INTERFACE pmc_c_putbuffer
[1762]155
[1896]156    INTERFACE pmc_c_setind_and_allocmem
157        MODULE PROCEDURE pmc_c_setind_and_allocmem
158    END INTERFACE pmc_c_setind_and_allocmem
[1762]159
[1896]160    INTERFACE pmc_c_set_dataarray
161        MODULE PROCEDURE pmc_c_set_dataarray_2d
162        MODULE PROCEDURE pmc_c_set_dataarray_3d
[2801]163        MODULE PROCEDURE pmc_c_set_dataarray_ip2d
[1896]164    END INTERFACE pmc_c_set_dataarray
[1762]165
[1896]166    INTERFACE pmc_set_dataarray_name
167        MODULE PROCEDURE pmc_set_dataarray_name
168        MODULE PROCEDURE pmc_set_dataarray_name_lastentry
169    END INTERFACE pmc_set_dataarray_name
[1762]170
171
[2801]172    PUBLIC pmc_childinit, pmc_c_clear_next_array_list, pmc_c_getbuffer,        &
173           pmc_c_getnextarray, pmc_c_putbuffer, pmc_c_setind_and_allocmem,     &
[1896]174           pmc_c_set_dataarray, pmc_set_dataarray_name, pmc_c_get_2d_index_list
[1762]175
[1896]176 CONTAINS
[1762]177
178
179
[1933]180 SUBROUTINE pmc_childinit
[1762]181
[1896]182     IMPLICIT NONE
[1762]183
[2801]184     INTEGER(iwp) ::  i        !<
185     INTEGER(iwp) ::  istat    !<
[1762]186
[1896]187!
188!--  Get / define the MPI environment
189     me%model_comm = m_model_comm
[1933]190     me%inter_comm = m_to_parent_comm
[1762]191
[1896]192     CALL MPI_COMM_RANK( me%model_comm, me%model_rank, istat )
193     CALL MPI_COMM_SIZE( me%model_comm, me%model_npes, istat )
194     CALL MPI_COMM_REMOTE_SIZE( me%inter_comm, me%inter_npes, istat )
195!
[2599]196!--  Intra-communicator is used for MPI_GET
[1896]197     CALL MPI_INTERCOMM_MERGE( me%inter_comm, .TRUE., me%intra_comm, istat )
198     CALL MPI_COMM_RANK( me%intra_comm, me%intra_rank, istat )
[1762]199
[1896]200     ALLOCATE( me%pes(me%inter_npes) )
[1779]201!
[2599]202!--  Allocate an array of type arraydef for all parent processes to store
203!--  information of then transfer array
[1896]204     DO  i = 1, me%inter_npes
205        ALLOCATE( me%pes(i)%array_list(pmc_max_array) )
206     ENDDO
[1762]207
[1933]208 END SUBROUTINE pmc_childinit
[1762]209
210
211
[2801]212 SUBROUTINE pmc_set_dataarray_name( parentarraydesc, parentarrayname,          &
[1933]213                                    childarraydesc, childarrayname, istat )
[1762]214
[1896]215    IMPLICIT NONE
[1762]216
[1933]217    CHARACTER(LEN=*), INTENT(IN) ::  parentarrayname  !<
218    CHARACTER(LEN=*), INTENT(IN) ::  parentarraydesc  !<
219    CHARACTER(LEN=*), INTENT(IN) ::  childarrayname   !<
220    CHARACTER(LEN=*), INTENT(IN) ::  childarraydesc   !<
[1762]221
[2801]222    INTEGER(iwp), INTENT(OUT) ::  istat  !<
[1896]223!
224!-- Local variables
225    TYPE(da_namedef) ::  myname  !<
[1762]226
[2801]227    INTEGER(iwp) ::  mype  !<
228    INTEGER(iwp) ::  my_addiarray = 0  !<
[1762]229
230
[1896]231    istat = pmc_status_ok
232!
233!-- Check length of array names
[2801]234    IF ( LEN( TRIM( parentarrayname) ) > da_namelen  .OR.                      &
[1933]235         LEN( TRIM( childarrayname) ) > da_namelen )  THEN
[1896]236       istat = pmc_da_name_err
237    ENDIF
[1762]238
[1896]239    IF ( m_model_rank == 0 )  THEN
240       myindex = myindex + 1
241       myname%couple_index = myIndex
[1933]242       myname%parentdesc   = TRIM( parentarraydesc )
243       myname%nameonparent = TRIM( parentarrayname )
244       myname%childdesc    = TRIM( childarraydesc )
245       myname%nameonchild  = TRIM( childarrayname )
[1896]246    ENDIF
[1762]247
[1896]248!
[2599]249!-- Broadcast to all child processes
[2801]250!
251!-- The complete description of an transfer names array is broadcasted
252
[1896]253    CALL pmc_bcast( myname%couple_index, 0, comm=m_model_comm )
[1933]254    CALL pmc_bcast( myname%parentdesc,   0, comm=m_model_comm )
255    CALL pmc_bcast( myname%nameonparent, 0, comm=m_model_comm )
256    CALL pmc_bcast( myname%childdesc,    0, comm=m_model_comm )
257    CALL pmc_bcast( myname%nameonchild,  0, comm=m_model_comm )
[1896]258!
[2599]259!-- Broadcast to all parent processes
[2801]260!-- The complete description of an transfer array names is broadcasted als to all parent processe
261!   Only the root PE of the broadcasts to parent using intra communicator
262
[1896]263    IF ( m_model_rank == 0 )  THEN
264        mype = MPI_ROOT
265    ELSE
266        mype = MPI_PROC_NULL
267    ENDIF
[1762]268
[1933]269    CALL pmc_bcast( myname%couple_index, mype, comm=m_to_parent_comm )
270    CALL pmc_bcast( myname%parentdesc,   mype, comm=m_to_parent_comm )
271    CALL pmc_bcast( myname%nameonparent, mype, comm=m_to_parent_comm )
272    CALL pmc_bcast( myname%childdesc,    mype, comm=m_to_parent_comm )
273    CALL pmc_bcast( myname%nameonchild,  mype, comm=m_to_parent_comm )
[1762]274
[1933]275    CALL pmc_g_setname( me, myname%couple_index, myname%nameonchild )
[1762]276
[1896]277 END SUBROUTINE pmc_set_dataarray_name
[1762]278
279
280
[1896]281 SUBROUTINE pmc_set_dataarray_name_lastentry( lastentry )
[1762]282
[1896]283    IMPLICIT NONE
[1762]284
[1896]285    LOGICAL, INTENT(IN), OPTIONAL ::  lastentry  !<
286!
287!-- Local variables
288    INTEGER ::  mype  !<
289    TYPE(dA_namedef) ::  myname  !<
[1762]290
[1896]291    myname%couple_index = -1
[1762]292
[1896]293    IF ( m_model_rank == 0 )  THEN
294       mype = MPI_ROOT
295    ELSE
296       mype = MPI_PROC_NULL
297    ENDIF
[1762]298
[1933]299    CALL pmc_bcast( myname%couple_index, mype, comm=m_to_parent_comm )
[1762]300
[1896]301 END SUBROUTINE pmc_set_dataarray_name_lastentry
[1762]302
303
304
[1896]305 SUBROUTINE pmc_c_get_2d_index_list
[1762]306
[1896]307    IMPLICIT NONE
[1762]308
[2801]309    INTEGER(iwp) :: dummy               !<
310    INTEGER(iwp) :: i, ierr, i2, j, nr  !<
311    INTEGER(iwp) :: indwin              !< MPI window object
312    INTEGER(iwp) :: indwin2             !< MPI window object
[1779]313
[1896]314    INTEGER(KIND=MPI_ADDRESS_KIND) :: win_size !< Size of MPI window 1 (in bytes)
315    INTEGER(KIND=MPI_ADDRESS_KIND) :: disp     !< Displacement unit (Integer = 4, floating poit = 8
316    INTEGER(KIND=MPI_ADDRESS_KIND) :: winsize  !< Size of MPI window 2 (in bytes)
[1779]317
[1896]318    INTEGER, DIMENSION(me%inter_npes*2) :: nrele  !< Number of Elements of a
319                                                  !< horizontal slice
320    INTEGER, DIMENSION(:), POINTER ::  myind  !<
[1779]321
[1896]322    TYPE(pedef), POINTER ::  ape  !> Pointer to pedef structure
[1762]323
324
[2809]325    win_size = STORAGE_SIZE( dummy )/8
[2801]326    CALL MPI_WIN_CREATE( dummy, win_size, iwp, MPI_INFO_NULL, me%intra_comm,   &
[1896]327                         indwin, ierr )
328!
[2801]329!-- Close window on child side and open on parent side
[1896]330    CALL MPI_WIN_FENCE( 0, indwin, ierr )
[2801]331
332!   Between the two MPI_WIN_FENCE calls, the parent can fill the RMA window
333
[1933]334!-- Close window on parent side and open on child side
[2801]335
[1896]336    CALL MPI_WIN_FENCE( 0, indwin, ierr )
[1762]337
[1896]338    DO  i = 1, me%inter_npes
339       disp = me%model_rank * 2
[2801]340       CALL MPI_GET( nrele((i-1)*2+1), 2, MPI_INTEGER, i-1, disp, 2,           &
[1896]341                     MPI_INTEGER, indwin, ierr )
342    ENDDO
343!
344!-- MPI_GET is non-blocking -> data in nrele is not available until MPI_FENCE is
345!-- called
346    CALL MPI_WIN_FENCE( 0, indwin, ierr )
347!
348!-- Allocate memory for index array
349    winsize = 0
350    DO  i = 1, me%inter_npes
351       ape => me%pes(i)
352       i2 = ( i-1 ) * 2 + 1
353       nr = nrele(i2+1)
354       IF ( nr > 0 )  THEN
355          ALLOCATE( ape%locind(nr) )
356       ELSE
357          NULLIFY( ape%locind )
358       ENDIF
359       winsize = MAX( nr, winsize )
360    ENDDO
[1762]361
[1896]362    ALLOCATE( myind(2*winsize) )
363    winsize = 1
364!
365!-- Local buffer used in MPI_GET can but must not be inside the MPI Window.
[2599]366!-- Here, we use a dummy for the MPI window because the parent processes do
367!-- not access the RMA window via MPI_GET or MPI_PUT
[2801]368    CALL MPI_WIN_CREATE( dummy, winsize, iwp, MPI_INFO_NULL, me%intra_comm,    &
[1896]369                         indwin2, ierr )
370!
371!-- MPI_GET is non-blocking -> data in nrele is not available until MPI_FENCE is
372!-- called
[2801]373
[1896]374    CALL MPI_WIN_FENCE( 0, indwin2, ierr )
[2801]375
376!   Between the two MPI_WIN_FENCE calls, the parent can fill the RMA window
377
[1896]378    CALL MPI_WIN_FENCE( 0, indwin2, ierr )
[1779]379
[1896]380    DO  i = 1, me%inter_npes
381       ape => me%pes(i)
382       nr = nrele(i*2)
383       IF ( nr > 0 )  THEN
384          disp = nrele(2*(i-1)+1)
385          CALL MPI_WIN_LOCK( MPI_LOCK_SHARED , i-1, 0, indwin2, ierr )
[2801]386          CALL MPI_GET( myind, 2*nr, MPI_INTEGER, i-1, disp, 2*nr,             &
[1896]387                        MPI_INTEGER, indwin2, ierr )
388          CALL MPI_WIN_UNLOCK( i-1, indwin2, ierr )
389          DO  j = 1, nr
390             ape%locind(j)%i = myind(2*j-1)
391             ape%locind(j)%j = myind(2*j)
392          ENDDO
393          ape%nrele = nr
394       ELSE
395          ape%nrele = -1
396       ENDIF
397    ENDDO
398!
399!-- Don't know why, but this barrier is necessary before we can free the windows
400    CALL MPI_BARRIER( me%intra_comm, ierr )
[1779]401
[1896]402    CALL MPI_WIN_FREE( indWin,  ierr )
403    CALL MPI_WIN_FREE( indwin2, ierr )
404    DEALLOCATE( myind )
[1762]405
[1896]406 END SUBROUTINE pmc_c_get_2d_index_list
[1762]407
[1779]408
409
[1896]410 SUBROUTINE pmc_c_clear_next_array_list
[1762]411
[1896]412    IMPLICIT NONE
[1762]413
[1896]414    next_array_in_list = 0
[1762]415
[1896]416 END SUBROUTINE pmc_c_clear_next_array_list
[1762]417
418
[1779]419
[1896]420 LOGICAL FUNCTION pmc_c_getnextarray( myname )
[1933]421
[1896]422!
423!--  List handling is still required to get minimal interaction with
424!--  pmc_interface
425     CHARACTER(LEN=*), INTENT(OUT) ::  myname  !<
426!
427!-- Local variables
428    TYPE(pedef), POINTER    :: ape
429    TYPE(arraydef), POINTER :: ar
[1779]430
431
[1896]432    next_array_in_list = next_array_in_list + 1
433!
[2599]434!-- Array names are the same on all child PEs, so take first process to
435!-- get the name   
[1896]436    ape => me%pes(1)
437!
438!-- Check if all arrays have been processed
439    IF ( next_array_in_list > ape%nr_arrays )  THEN
440       pmc_c_getnextarray = .FALSE.
441       RETURN
442    ENDIF
[1762]443
[1896]444    ar => ape%array_list( next_array_in_list )
[1762]445
[1896]446    myname = ar%name
447!
[2801]448!-- Return true if annother array
449!-- If all array have been processed, the RETURN statement a couple of lines above is active
450
[1896]451    pmc_c_getnextarray = .TRUE.
[1762]452
[2801]453 END FUNCTION pmc_c_getnextarray
[1764]454
[1762]455
[1786]456
[1896]457 SUBROUTINE pmc_c_set_dataarray_2d( array )
[1762]458
[1896]459    IMPLICIT NONE
[1762]460
[1896]461    REAL(wp), INTENT(IN) , DIMENSION(:,:), POINTER ::  array  !<
[1762]462
[2801]463    INTEGER(iwp)               ::  i       !<
464    INTEGER(iwp)               ::  nrdims  !<
465    INTEGER(iwp), DIMENSION(4) ::  dims    !<
[1786]466
[1896]467    TYPE(C_PTR)             ::  array_adr
468    TYPE(arraydef), POINTER ::  ar
469    TYPE(pedef), POINTER    ::  ape
[1762]470
471
[1896]472    dims    = 1
473    nrdims  = 2
474    dims(1) = SIZE( array, 1 )
475    dims(2) = SIZE( array, 2 )
[1762]476
[1896]477    array_adr = C_LOC( array )
[1786]478
[1896]479    DO  i = 1, me%inter_npes
480       ape => me%pes(i)
481       ar  => ape%array_list(next_array_in_list)
482       ar%nrdims = nrdims
[2801]483       ar%dimkey = nrdims
[1896]484       ar%a_dim  = dims
485       ar%data   = array_adr
486    ENDDO
[1786]487
[1896]488 END SUBROUTINE pmc_c_set_dataarray_2d
[1786]489
[2801]490 SUBROUTINE pmc_c_set_dataarray_ip2d( array )
[1786]491
[2801]492    IMPLICIT NONE
[1786]493
[2801]494    INTEGER(idp), INTENT(IN) , DIMENSION(:,:), POINTER ::  array  !<
495
496    INTEGER(iwp)               ::  i       !<
497    INTEGER(iwp)               ::  nrdims  !<
498    INTEGER(iwp), DIMENSION(4) ::  dims    !<
499
500    TYPE(C_PTR)             ::  array_adr
501    TYPE(arraydef), POINTER ::  ar
502    TYPE(pedef), POINTER    ::  ape
503
504    dims    = 1
505    nrdims  = 2
506    dims(1) = SIZE( array, 1 )
507    dims(2) = SIZE( array, 2 )
508
509    array_adr = C_LOC( array )
510
511    DO  i = 1, me%inter_npes
512       ape => me%pes(i)
513       ar  => ape%array_list(next_array_in_list)
514       ar%nrdims = nrdims
515       ar%dimkey = 22
516       ar%a_dim  = dims
517       ar%data   = array_adr
518    ENDDO
519
520 END SUBROUTINE pmc_c_set_dataarray_ip2d
521
[1896]522 SUBROUTINE pmc_c_set_dataarray_3d (array)
[1786]523
[1896]524    IMPLICIT NONE
[1786]525
[1896]526    REAL(wp), INTENT(IN), DIMENSION(:,:,:), POINTER ::  array  !<
[1786]527
[2801]528    INTEGER(iwp)                ::  i
529    INTEGER(iwp)                ::  nrdims
530    INTEGER(iwp), DIMENSION (4) ::  dims
531   
[1896]532    TYPE(C_PTR)             ::  array_adr
533    TYPE(pedef), POINTER    ::  ape
534    TYPE(arraydef), POINTER ::  ar
[1786]535
536
[1896]537    dims    = 1
538    nrdims  = 3
539    dims(1) = SIZE( array, 1 )
540    dims(2) = SIZE( array, 2 )
541    dims(3) = SIZE( array, 3 )
[1786]542
[1896]543    array_adr = C_LOC( array )
[1786]544
[1896]545    DO  i = 1, me%inter_npes
546       ape => me%pes(i)
547       ar  => ape%array_list(next_array_in_list)
548       ar%nrdims = nrdims
[2801]549       ar%dimkey = nrdims
[1896]550       ar%a_dim  = dims
551       ar%data   = array_adr
552    ENDDO
[1786]553
[1896]554 END SUBROUTINE pmc_c_set_dataarray_3d
[1762]555
556
557
[1896]558 SUBROUTINE pmc_c_setind_and_allocmem
[1762]559
[1896]560    IMPLICIT NONE
[1933]561
[1896]562!
[1933]563!-- Naming convention for appendices:  _pc  -> parent to child transfer
564!--                                    _cp  -> child to parent transfer
565!--                                    recv -> parent to child transfer
566!--                                    send -> child to parent transfer
[1896]567    CHARACTER(LEN=da_namelen) ::  myname  !<
[1786]568
[2801]569    INTEGER(iwp) ::  arlen    !<
570    INTEGER(iwp) ::  myindex  !<
571    INTEGER(iwp) ::  i        !<
572    INTEGER(iwp) ::  ierr     !<
573    INTEGER(iwp) ::  istat    !<
574    INTEGER(iwp) ::  j        !<
575    INTEGER(iwp) ::  rcount   !<
576    INTEGER(iwp) ::  tag      !<
[1762]577
[2801]578    INTEGER(iwp), PARAMETER ::  noindex = -1  !<
[1762]579
[1896]580    INTEGER(idp)                   ::  bufsize  !< size of MPI data window
581    INTEGER(KIND=MPI_ADDRESS_KIND) ::  winsize  !<
[1762]582
[2801]583    INTEGER(iwp),DIMENSION(1024) ::  req  !<
[1779]584
[1933]585    REAL(wp), DIMENSION(:), POINTER, SAVE ::  base_array_pc  !< base array
586    REAL(wp), DIMENSION(:), POINTER, SAVE ::  base_array_cp  !< base array
[1762]587
[1896]588    TYPE(pedef), POINTER    ::  ape       !<
589    TYPE(arraydef), POINTER ::  ar        !<
590    Type(C_PTR)             ::  base_ptr  !<
[1779]591
[1762]592
[1896]593    myindex = 0
594    bufsize = 8
[1797]595!
[1933]596!-- Parent to child direction.
[1896]597!-- First stride: compute size and set index
598    DO  i = 1, me%inter_npes
599       ape => me%pes(i)
600       tag = 200
601       DO  j = 1, ape%nr_arrays
602          ar => ape%array_list(j)
603!
[1933]604!--       Receive index from child
[1896]605          tag = tag + 1
[2801]606          CALL MPI_RECV( myindex, 1, MPI_INTEGER, i-1, tag, me%inter_comm,     &
[1896]607                         MPI_STATUS_IGNORE, ierr )
608          ar%recvindex = myindex
609!
[1933]610!--       Determine max, because child buffer is allocated only once
[2801]611!--       All 2D and 3d arrays use the same buffer
612
613          IF ( ar%nrdims == 3 )  THEN
[1896]614             bufsize = MAX( bufsize, ar%a_dim(1)*ar%a_dim(2)*ar%a_dim(3) )
615          ELSE
616             bufsize = MAX( bufsize, ar%a_dim(1)*ar%a_dim(2) )
617          ENDIF
618       ENDDO
619    ENDDO
[1779]620!
[1896]621!-- Create RMA (one sided communication) data buffer.
622!-- The buffer for MPI_GET can be PE local, i.e. it can but must not be part of
623!-- the MPI RMA window
[1933]624    CALL pmc_alloc_mem( base_array_pc, bufsize, base_ptr )
[1896]625    me%totalbuffersize = bufsize*wp  ! total buffer size in byte
626!
627!-- Second stride: set buffer pointer
628    DO  i = 1, me%inter_npes
629       ape => me%pes(i)
630       DO  j = 1, ape%nr_arrays
631          ar => ape%array_list(j)
632          ar%recvbuf = base_ptr
633       ENDDO
634    ENDDO
635!
[1933]636!-- Child to parent direction
[1896]637    myindex = 1
638    rcount  = 0
639    bufsize = 8
[1762]640
[1896]641    DO  i = 1, me%inter_npes
642       ape => me%pes(i)
643       tag = 300
644       DO  j = 1, ape%nr_arrays
645          ar => ape%array_list(j)
[2801]646          IF ( ar%nrdims == 2 )  THEN
[1896]647             arlen = ape%nrele
648          ELSEIF( ar%nrdims == 3 )  THEN
649             arlen = ape%nrele*ar%a_dim(1)
650          ENDIF
651          tag    = tag + 1
652          rcount = rcount + 1
653          IF ( ape%nrele > 0 )  THEN
[2801]654             CALL MPI_ISEND( myindex, 1, MPI_INTEGER, i-1, tag, me%inter_comm, &
[1896]655                             req(rcount), ierr )
656             ar%sendindex = myindex
657          ELSE
[2801]658             CALL MPI_ISEND( noindex, 1, MPI_INTEGER, i-1, tag, me%inter_comm, &
[1896]659                             req(rcount), ierr )
660             ar%sendindex = noindex
661          ENDIF
662!
[2801]663!--       Maximum of 1024 pending requests
664!         1024 is an arbitrary value just to make sure the number of pending
665!         requests is getting too large. It is possible that this value has to
666!         be adjusted in case of running the model on large number of cores.
667
[1896]668          IF ( rcount == 1024 )  THEN
669             CALL MPI_WAITALL( rcount, req, MPI_STATUSES_IGNORE, ierr )
670             rcount = 0
671          ENDIF
[1762]672
[1896]673          IF ( ape%nrele > 0 )  THEN
674             ar%sendsize = arlen
675             myindex     = myindex + arlen
676             bufsize     = bufsize + arlen
677          ENDIF
[1762]678
[1896]679       ENDDO
[1762]680
[1896]681       IF ( rcount > 0 )  THEN
682          CALL MPI_WAITALL( rcount, req, MPI_STATUSES_IGNORE, ierr )
683       ENDIF
[1762]684
[1896]685    ENDDO
686!
[1933]687!-- Create RMA (one sided communication) window for data buffer child to parent
[1896]688!-- transfer.
689!-- The buffer of MPI_GET (counter part of transfer) can be PE-local, i.e. it
690!-- can but must not be part of the MPI RMA window. Only one RMA window is
691!-- required to prepare the data
[1933]692!--        for parent -> child transfer on the parent side
[1896]693!-- and
[1933]694!--        for child -> parent transfer on the child side
695    CALL pmc_alloc_mem( base_array_cp, bufsize )
[1896]696    me%totalbuffersize = bufsize * wp  ! total buffer size in byte
697
698    winSize = me%totalbuffersize
699
[2801]700    CALL MPI_WIN_CREATE( base_array_cp, winsize, wp, MPI_INFO_NULL,            &
[1933]701                         me%intra_comm, me%win_parent_child, ierr )
702    CALL MPI_WIN_FENCE( 0, me%win_parent_child, ierr )
[1896]703    CALL MPI_BARRIER( me%intra_comm, ierr )
704!
705!-- Second stride: set buffer pointer
706    DO  i = 1, me%inter_npes
707       ape => me%pes(i)
708       DO  j = 1, ape%nr_arrays
[2599]709          ar => ape%array_list(j)         
[1896]710          IF ( ape%nrele > 0 )  THEN
[1933]711             ar%sendbuf = C_LOC( base_array_cp(ar%sendindex) )
712!
[1896]713!--          TODO: if this is an error to be really expected, replace the
714!--                following message by a meaningful standard PALM message using
715!--                the message-routine
716             IF ( ar%sendindex+ar%sendsize > bufsize )  THEN
[2801]717                WRITE( 0,'(a,i4,4i7,1x,a)') 'Child buffer too small ', i,      &
718                          ar%sendindex, ar%sendsize, ar%sendindex+ar%sendsize, &
[1896]719                          bufsize, TRIM( ar%name )
720                CALL MPI_ABORT( MPI_COMM_WORLD, istat, ierr )
721             ENDIF
722          ENDIF
723       ENDDO
724    ENDDO
725
726 END SUBROUTINE pmc_c_setind_and_allocmem
727
728
729
[2801]730 SUBROUTINE pmc_c_getbuffer( waittime, particle_transfer )
[1896]731
732    IMPLICIT NONE
733
734    REAL(wp), INTENT(OUT), OPTIONAL ::  waittime  !<
[2801]735    LOGICAL, INTENT(IN), OPTIONAL   ::  particle_transfer  !<
[1896]736
737    CHARACTER(LEN=da_namelen) ::  myname  !<
[2801]738   
739    LOGICAL                        ::  lo_ptrans!<
740   
741    INTEGER(iwp)                        ::  ierr    !<
742    INTEGER(iwp)                        ::  ij      !<
743    INTEGER(iwp)                        ::  ip      !<
744    INTEGER(iwp)                        ::  j       !<
745    INTEGER(iwp)                        ::  myindex !<
746    INTEGER(iwp)                        ::  nr      !< number of elements to get
747                                                    !< from parent
[1896]748    INTEGER(KIND=MPI_ADDRESS_KIND) ::  target_disp
749    INTEGER,DIMENSION(1)           ::  buf_shape
750
751    REAL(wp)                            ::  t1
752    REAL(wp)                            ::  t2
753
754    REAL(wp), POINTER, DIMENSION(:)     ::  buf
755    REAL(wp), POINTER, DIMENSION(:,:)   ::  data_2d
756    REAL(wp), POINTER, DIMENSION(:,:,:) ::  data_3d
757    TYPE(pedef), POINTER                ::  ape
758    TYPE(arraydef), POINTER             ::  ar
[2801]759    INTEGER(idp), POINTER, DIMENSION(:)     ::  ibuf      !<
760    INTEGER(idp), POINTER, DIMENSION(:,:)   ::  idata_2d  !<
[1896]761
762!
[1933]763!-- Synchronization of the model is done in pmci_synchronize.
764!-- Therefore the RMA window can be filled without
[1896]765!-- sychronization at this point and a barrier is not necessary.
[2801]766
767!-- In case waittime is present, the following barrier is necessary to
768!-- insure the same number of barrier calls on parent and child
769!-- This means, that here on child side two barriers are call successively
770!-- The parent is filling its buffer between the two barrier calls
771
[1896]772!-- Please note that waittime has to be set in pmc_s_fillbuffer AND
773!-- pmc_c_getbuffer
774    IF ( PRESENT( waittime ) )  THEN
775       t1 = pmc_time()
776       CALL MPI_BARRIER( me%intra_comm, ierr )
777       t2 = pmc_time()
778       waittime = t2 - t1
779    ENDIF
[2801]780
781    lo_ptrans = .FALSE.
782    IF ( PRESENT( particle_transfer))    lo_ptrans = particle_transfer
783
[1896]784!
[1933]785!-- Wait for buffer is filled.
[2801]786!
787!-- The parent side (in pmc_s_fillbuffer) is filling the buffer in the MPI RMA window
788!-- When the filling is complet, a MPI_BARRIER is called.
789!-- The child is not allowd to access the parent-buffer before it is completely filled
790!-- therefore the following barrier is required.
791
[1896]792    CALL MPI_BARRIER( me%intra_comm, ierr )
793
794    DO  ip = 1, me%inter_npes
795       ape => me%pes(ip)
796       DO  j = 1, ape%nr_arrays
797          ar => ape%array_list(j)
[2801]798
799          IF ( ar%dimkey == 2 .AND. .NOT.lo_ptrans)  THEN
[1896]800             nr = ape%nrele
[2801]801          ELSEIF ( ar%dimkey == 3 .AND. .NOT.lo_ptrans)  THEN
[1896]802             nr = ape%nrele * ar%a_dim(1)
[2801]803          ELSE IF ( ar%dimkey == 22 .AND. lo_ptrans)  THEN
804             nr = ape%nrele
805          ELSE
806             CYCLE                    ! Particle array ar not transferd here
[1896]807          ENDIF
808          buf_shape(1) = nr
[2801]809          IF ( lo_ptrans )   THEN
810             CALL C_F_POINTER( ar%recvbuf, ibuf, buf_shape )
811          ELSE
812             CALL C_F_POINTER( ar%recvbuf, buf, buf_shape )
813          ENDIF
[1896]814!
815!--       MPI passive target RMA
[2801]816!--       One data array is fetcht from MPI RMA window on parent
817
[1896]818          IF ( nr > 0 )  THEN
819             target_disp = ar%recvindex - 1
[2801]820             CALL MPI_WIN_LOCK( MPI_LOCK_SHARED , ip-1, 0,                     &
[1933]821                                me%win_parent_child, ierr )
[2801]822             IF ( lo_ptrans )   THEN
823                CALL MPI_GET( ibuf, nr*8, MPI_BYTE, ip-1, target_disp, nr*8, MPI_BYTE,  &               !There is no MPI_INTEGER8 datatype
824                                   me%win_parent_child, ierr )
825             ELSE
826                CALL MPI_GET( buf, nr, MPI_REAL, ip-1, target_disp, nr,        &
827                              MPI_REAL, me%win_parent_child, ierr )
828             ENDIF
[1933]829             CALL MPI_WIN_UNLOCK( ip-1, me%win_parent_child, ierr )
[1896]830          ENDIF
831          myindex = 1
[2801]832          IF ( ar%dimkey == 2 .AND. .NOT.lo_ptrans)  THEN
833
[1896]834             CALL C_F_POINTER( ar%data, data_2d, ar%a_dim(1:2) )
835             DO  ij = 1, ape%nrele
836                data_2d(ape%locind(ij)%j,ape%locind(ij)%i) = buf(myindex)
837                myindex = myindex + 1
838             ENDDO
[2801]839
840          ELSEIF ( ar%dimkey == 3 .AND. .NOT.lo_ptrans)  THEN
841
[1896]842             CALL C_F_POINTER( ar%data, data_3d, ar%a_dim(1:3) )
843             DO  ij = 1, ape%nrele
[2801]844                data_3d(:,ape%locind(ij)%j,ape%locind(ij)%i) =                 &
[1896]845                                              buf(myindex:myindex+ar%a_dim(1)-1)
846                myindex = myindex+ar%a_dim(1)
847             ENDDO
[2801]848
849          ELSEIF ( ar%dimkey == 22 .AND. lo_ptrans)  THEN
850             CALL C_F_POINTER( ar%data, idata_2d, ar%a_dim(1:2) )
851
852             DO  ij = 1, ape%nrele
853                idata_2d(ape%locind(ij)%j,ape%locind(ij)%i) = ibuf(myindex)
854                myindex = myindex + 1
855             ENDDO
856
[1896]857          ENDIF
858       ENDDO
859    ENDDO
860
861 END SUBROUTINE pmc_c_getbuffer
862
863
864
[2801]865 SUBROUTINE pmc_c_putbuffer( waittime , particle_transfer )
[1896]866
867    IMPLICIT NONE
868
869    REAL(wp), INTENT(OUT), OPTIONAL ::  waittime  !<
[2801]870    LOGICAL, INTENT(IN), OPTIONAL   ::  particle_transfer  !<
[1896]871
872    CHARACTER(LEN=da_namelen) ::  myname  !<
[2801]873   
874    LOGICAL ::  lo_ptrans!<
875   
876    INTEGER(iwp) ::  ierr         !<
877    INTEGER(iwp) ::  ij           !<
878    INTEGER(iwp) ::  ip           !<
879    INTEGER(iwp) ::  j            !<
880    INTEGER(iwp) ::  myindex      !<
881    INTEGER(iwp) ::  nr           !< number of elements to get from parent
882   
[1896]883    INTEGER(KIND=MPI_ADDRESS_KIND) ::  target_disp  !<
[2801]884   
[1896]885
[2801]886    INTEGER(iwp), DIMENSION(1) ::  buf_shape    !<
[1896]887
888    REAL(wp) ::  t1  !<
889    REAL(wp) ::  t2  !<
890
[2801]891    REAL(wp), POINTER, DIMENSION(:)         ::  buf      !<
892    REAL(wp), POINTER, DIMENSION(:,:)       ::  data_2d  !<
893    REAL(wp), POINTER, DIMENSION(:,:,:)     ::  data_3d  !<
894   
895    INTEGER(idp), POINTER, DIMENSION(:)     ::  ibuf      !<
896    INTEGER(idp), POINTER, DIMENSION(:,:)   ::  idata_2d  !<
[1896]897
[2801]898    TYPE(pedef), POINTER                    ::  ape  !<
899    TYPE(arraydef), POINTER                 ::  ar   !<
[1896]900
901!
902!-- Wait for empty buffer
[2801]903!-- Switch RMA epoche
904
[1896]905    t1 = pmc_time()
906    CALL MPI_BARRIER( me%intra_comm, ierr )
907    t2 = pmc_time()
908    IF ( PRESENT( waittime ) )  waittime = t2 - t1
909
[2801]910    lo_ptrans = .FALSE.
911    IF ( PRESENT( particle_transfer))    lo_ptrans = particle_transfer
912
[1896]913    DO  ip = 1, me%inter_npes
914       ape => me%pes(ip)
915       DO  j = 1, ape%nr_arrays
916          ar => aPE%array_list(j)
917          myindex = 1
[2801]918
919          IF ( ar%dimkey == 2 .AND. .NOT.lo_ptrans )  THEN
920
[1896]921             buf_shape(1) = ape%nrele
922             CALL C_F_POINTER( ar%sendbuf, buf,     buf_shape     )
923             CALL C_F_POINTER( ar%data,    data_2d, ar%a_dim(1:2) )
924             DO  ij = 1, ape%nrele
925                buf(myindex) = data_2d(ape%locind(ij)%j,ape%locind(ij)%i)
926                myindex = myindex + 1
927             ENDDO
[2801]928
929          ELSEIF ( ar%dimkey == 3 .AND. .NOT.lo_ptrans )  THEN
930
[1896]931             buf_shape(1) = ape%nrele*ar%a_dim(1)
932             CALL C_F_POINTER( ar%sendbuf, buf,     buf_shape     )
933             CALL C_F_POINTER( ar%data,    data_3d, ar%a_dim(1:3) )
934             DO  ij = 1, ape%nrele
[1933]935                buf(myindex:myindex+ar%a_dim(1)-1) =                            &
[1896]936                                    data_3d(:,ape%locind(ij)%j,ape%locind(ij)%i)
937                myindex = myindex + ar%a_dim(1)
938             ENDDO
[2801]939
940          ELSE IF ( ar%dimkey == 22 .AND. lo_ptrans)  THEN
941
942             buf_shape(1) = ape%nrele
943             CALL C_F_POINTER( ar%sendbuf, ibuf,     buf_shape     )
944             CALL C_F_POINTER( ar%data,    idata_2d, ar%a_dim(1:2) )
945
946             DO  ij = 1, ape%nrele
947                ibuf(myindex) = idata_2d(ape%locind(ij)%j,ape%locind(ij)%i)
948                myindex = myindex + 1
949             ENDDO
950
[1896]951          ENDIF
952       ENDDO
953    ENDDO
954!
955!-- Buffer is filled
[2801]956!-- Switch RMA epoche
957
[1896]958    CALL MPI_Barrier(me%intra_comm, ierr)
959
960 END SUBROUTINE pmc_c_putbuffer
961
[1764]962#endif
[1933]963 END MODULE pmc_child
Note: See TracBrowser for help on using the repository browser.