source: palm/trunk/SOURCE/pmc_server.f90 @ 1767

Last change on this file since 1767 was 1767, checked in by raasch, 8 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 26.8 KB
RevLine 
[1762]1MODULE pmc_server
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 terms
7! of the GNU General Public License as published by the Free Software Foundation,
8! either version 3 of the License, or (at your option) any later version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2015 Leibniz Universitaet Hannover
18!--------------------------------------------------------------------------------!
19!
20! Current revisions:
21! ------------------
22!
[1767]23!
[1762]24! Former revisions:
25! -----------------
26! $Id: pmc_server.f90 1767 2016-02-29 08:47:57Z raasch $
27!
[1767]28! 1766 2016-02-29 08:37:15Z raasch
29! modifications to allow for using PALM's pointer version
30! +new routine pmc_s_set_active_data_array
31!
[1765]32! 1764 2016-02-28 12:45:19Z raasch
33! cpp-statement added (nesting can only be used in parallel mode)
34!
[1763]35! 1762 2016-02-25 12:31:13Z hellstea
36! Initial revision by K. Ketelsen
[1762]37!
38! Description:
39! ------------
40!
41! Server part of Palm Model Coupler
42!------------------------------------------------------------------------------!
43
[1764]44#if defined( __parallel )
[1762]45   use, intrinsic :: iso_c_binding
46
[1764]47#if defined( __lc )
48    USE MPI
49#else
50    INCLUDE "mpif.h"
51#endif
52   USE  kinds
[1762]53   USE  PMC_general,               ONLY: ClientDef, PMC_MAX_MODELL,PMC_sort, DA_NameDef, DA_Desclen, DA_Namelen,       &
54                                         PMC_G_SetName, PMC_G_GetName, PeDef, ArrayDef
55   USE  PMC_handle_communicator,   ONLY: m_model_comm,m_model_rank,m_model_npes, m_to_client_comm,                     &
56                                         PMC_Server_for_Client, m_world_rank
57   USE   PMC_MPI_wrapper,          ONLY: PMC_Send_to_Client, PMC_Recv_from_Client, PMC_Bcast, PMC_Inter_Bcast,         &
58                                         PMC_Alloc_mem, PMC_Time
59
60   IMPLICIT none
61   PRIVATE
62   SAVE
63
64   TYPE ClientIndexDef
65      INTEGER                                        :: NrPoints
66      INTEGER,DIMENSION(:,:),allocatable             :: index_list_2d
67   END TYPE ClientIndexDef
68
69   TYPE(ClientDef),DIMENSION(PMC_MAX_MODELL)          :: Clients
70   TYPE(ClientIndexDef),DIMENSION(PMC_MAX_MODELL)     :: indClients
71
72   PUBLIC PMC_Server_for_Client
73
[1764]74!-- TO_DO: what is the meaning of this? Could variables declared in this module
75!--        also have single precision?
76!   INTEGER, PARAMETER :: dp = wp
[1762]77
78   ! INTERFACE section
79
80   INTERFACE PMC_ServerInit
81      MODULE procedure  PMC_ServerInit
82   END INTERFACE PMC_ServerInit
83
84    INTERFACE PMC_S_Set_2D_index_list
85        MODULE procedure PMC_S_Set_2D_index_list
86    END INTERFACE PMC_S_Set_2D_index_list
87
88    INTERFACE PMC_S_GetNextArray
89        MODULE procedure PMC_S_GetNextArray
90    END INTERFACE PMC_S_GetNextArray
91
92    INTERFACE PMC_S_Set_DataArray
93        MODULE procedure PMC_S_Set_DataArray_2d
94        MODULE procedure PMC_S_Set_DataArray_3d
95    END INTERFACE PMC_S_Set_DataArray
96
97    INTERFACE PMC_S_setInd_and_AllocMem
98        MODULE procedure PMC_S_setInd_and_AllocMem
99    END INTERFACE PMC_S_setInd_and_AllocMem
100
101    INTERFACE PMC_S_FillBuffer
102        MODULE procedure PMC_S_FillBuffer
103    END INTERFACE PMC_S_FillBuffer
104
105    INTERFACE PMC_S_GetData_from_Buffer
106        MODULE procedure PMC_S_GetData_from_Buffer
107    END INTERFACE PMC_S_GetData_from_Buffer
108
[1766]109    INTERFACE PMC_S_Set_Active_data_array
110        MODULE procedure PMC_S_Set_Active_data_array
111    END INTERFACE PMC_S_Set_Active_data_array
112
[1762]113    ! PUBLIC section
114
115    PUBLIC PMC_ServerInit, PMC_S_Set_2D_index_list, PMC_S_GetNextArray, PMC_S_Set_DataArray
[1766]116    PUBLIC PMC_S_setInd_and_AllocMem, PMC_S_FillBuffer, PMC_S_GetData_from_Buffer, PMC_S_Set_Active_data_array
[1762]117
118CONTAINS
119
120   SUBROUTINE PMC_ServerInit
121      IMPLICIT none
122      INTEGER                 :: i
123      INTEGER                 :: j
124      INTEGER                 :: ClientId
125      INTEGER                 :: istat
126
127      do i=1,size(PMC_Server_for_Client)-1
128         if(m_model_comm == 0) write(0,*) 'PMC_Server: Initialize client Id',PMC_Server_for_Client(i)
129
130         ClientId = PMC_Server_for_Client(i)
131
132         Clients(ClientId)%model_comm = m_model_comm
133         Clients(ClientId)%inter_comm = m_to_client_comm(ClientId)
134
135         ! Get rank and size
136         CALL MPI_Comm_rank (Clients(ClientId)%model_comm, Clients(ClientId)%model_rank, istat);
137         CALL MPI_Comm_size (Clients(ClientId)%model_comm, Clients(ClientId)%model_npes, istat);
138         CALL MPI_Comm_remote_size (Clients(ClientId)%inter_comm, Clients(ClientId)%inter_npes, istat);
139
140         ! Intra communicater is used for MPI_Get
141         CALL MPI_Intercomm_merge (Clients(ClientId)%inter_comm, .false., Clients(ClientId)%intra_comm, istat);
142         CALL MPI_Comm_rank (Clients(ClientId)%intra_comm, Clients(ClientId)%intra_rank, istat);
143
144         write(9,*) 'ClientId ',i,ClientId,m_world_rank, Clients(ClientId)%inter_npes
145
146         ALLOCATE (Clients(ClientId)%PEs(Clients(ClientId)%inter_npes))
147
148         do j=1,Clients(ClientId)%inter_npes                                ! Loop over all client PEs
149           NULLIFY(Clients(ClientId)%PEs(j)%Arrays)
150         end do
151
152         CALL Get_DA_names_from_client (ClientId)
153      end do
154
155      return
156   END SUBROUTINE PMC_ServerInit
157
158    SUBROUTINE PMC_S_Set_2D_index_list (ClientId, index_list)
159        IMPLICIT none
160        INTEGER,INTENT(IN)                         :: ClientId
161        INTEGER,DIMENSION(:,:),INTENT(INOUT)       :: index_list      !Index list will be modified in sort, therefore INOUT
162
163        !-- Local variables
164        INTEGER                 :: ip,is,ie,ian,ic,n
165        INTEGER                 :: istat
166
167        if(m_model_rank == 0)   then
168            CALL PMC_sort (index_list, 6)                       ! Sort to ascending Server PE
169            is = 1
170
171            do ip=0,m_model_npes-1
172
173                !       Split into Server PEs
174                ie = is-1                                     !there may be no entry for this PE
175                if(is <= size(index_list,2) .and. ie >= 0)  then
176                    do while ( index_list(6,ie+1) == ip)
177                        ie = ie+1
178                        if( ie == size(index_list,2)) EXIT
179                    end do
180
181                    ian = ie-is+1
182                else
183                    is  = -1
184                    ie  = -2
185                    ian = 0
186                end if
187
188                !       Send data to other server PEs
189
190                if(ip == 0)   then
191                    indClients(ClientId)%NrPoints = ian
192                    if(ian > 0)   then
193                        ALLOCATE (indClients(ClientId)%index_list_2d(6,ian))
194                        indClients(ClientId)%index_list_2d(:,1:ian) = index_list(:,is:ie)
195                    end if
196                else
197                    CALL MPI_Send (ian, 1, MPI_INTEGER, ip, 1000, m_model_comm, istat)
198                    if(ian > 0) then
199                        CALL MPI_Send (index_list(1,is), 6*ian, MPI_INTEGER, ip, 1001,                  &
200                            m_model_comm, istat)
201                    end if
202                end if
203                is = ie+1
204            end do
205        else
206            CALL MPI_Recv (indClients(ClientId)%NrPoints, 1, MPI_INTEGER, 0, 1000, m_model_comm,     &
207                MPI_STATUS_IGNORE, istat)
208            ian = indClients(ClientId)%NrPoints
209             if(ian > 0) then
210                ALLOCATE(indClients(ClientId)%index_list_2d(6,ian))
211                CALL MPI_RECV (indClients(ClientId)%index_list_2d, 6*ian, MPI_INTEGER, 0, 1001,        &
212                    m_model_comm, MPI_STATUS_IGNORE, istat)
213            end if
214        end if
215
216        CALL Set_PE_index_list (ClientId,Clients(ClientId),indClients(ClientId)%index_list_2d,indClients(ClientId)%NrPoints)
217
218        return
219    END SUBROUTINE PMC_S_Set_2D_index_list
220
221    logical function PMC_S_GetNextArray (ClientId, myName,Client_PeIndex)
222        INTEGER,INTENT(IN)                         :: ClientId
223        CHARACTER(len=*),INTENT(OUT)               :: myName
224
225        !-- local variables
226        INTEGER                      :: MyCoupleIndex
227        logical                      :: MyLast
228        CHARACTER(len=DA_Namelen)    :: loName
229        INTEGER,INTENT(IN),optional  :: Client_PeIndex
230
231        loName = ' '
232
233        CALL PMC_G_GetName (clients(ClientId), MyCoupleIndex, loName, MyLast, Client_PeIndex)
234
235        myName    = loName
236
237        PMC_S_GetNextArray = .NOT. MyLast                   ! Return true if valid array
238
239        return
240    END function PMC_S_GetNextArray
241
[1766]242    SUBROUTINE PMC_S_Set_DataArray_2d (ClientId, array, array_2 )
[1762]243        IMPLICIT none
244        INTEGER,INTENT(IN)                         :: ClientId
[1764]245!--   TO_DO: has array always to be of dp-kind, or can wp used here
246!--          this effects all respective declarations in this file
[1766]247        REAL(kind=dp),INTENT(IN),DIMENSION(:,:)           :: array
248        REAL(kind=dp),INTENT(IN),DIMENSION(:,:),OPTIONAL  :: array_2
[1762]249        !-- local variables
250        INTEGER                           :: NrDims
251        INTEGER,DIMENSION (4)             :: dims
252        INTEGER                           :: dim_order
253        TYPE(c_ptr)                       :: array_adr
[1766]254        TYPE(c_ptr)                       :: second_adr
[1762]255
256        dims = 1
257
258        NrDims    = 2
259        dims(1)   = size(array,1)
260        dims(2)   = size(array,2)
261        dim_order = 2
262
263        array_adr = c_loc(array)
264
[1766]265        IF ( PRESENT( array_2 ) )  THEN
266           second_adr = c_loc(array_2)
267           CALL PMC_S_SetArray (ClientId, NrDims, dims, dim_order, array_adr, second_adr = second_adr)
268        ELSE
269           CALL PMC_S_SetArray (ClientId, NrDims, dims, dim_order, array_adr)
270        ENDIF
[1762]271
272        return
273    END SUBROUTINE PMC_S_Set_DataArray_2d
274
[1766]275    SUBROUTINE PMC_S_Set_DataArray_3d (ClientId, array, nz_cl, nz, array_2 )
[1762]276        IMPLICIT none
277        INTEGER,INTENT(IN)                         :: ClientId
278        REAL(kind=dp),INTENT(IN),DIMENSION(:,:,:)  :: array
[1766]279        REAL(kind=dp),INTENT(IN),DIMENSION(:,:,:),OPTIONAL  :: array_2
[1762]280        INTEGER,INTENT(IN)                         :: nz_cl
281        INTEGER,INTENT(IN)                         :: nz
282        !-- local variables
283        INTEGER                           :: NrDims
284        INTEGER,DIMENSION (4)             :: dims
285        INTEGER                           :: dim_order
286        TYPE(c_ptr)                       :: array_adr
[1766]287        TYPE(c_ptr)                       :: second_adr
[1762]288
289        dims = 1
290
291        dims      = 0
292        NrDims    = 3
293        dims(1)   = size(array,1)
294        dims(2)   = size(array,2)
295        dims(3)   = size(array,3)
296        dim_order = 33
297        dims(4)   = nz_cl+dims(1)-nz                        ! works for first dimension 1:nz and 0:nz+1
298
299        array_adr = c_loc(array)
300
[1766]301!
302!--     In PALM's pointer version, two indices have to be stored internally.
303!--     The active address of the data array is set in swap_timelevel
304        IF ( PRESENT( array_2 ) )  THEN
305          second_adr = c_loc(array_2)
306          CALL PMC_S_SetArray (ClientId, NrDims, dims, dim_order, array_adr, second_adr = second_adr)
307        ELSE
308           CALL PMC_S_SetArray (ClientId, NrDims, dims, dim_order, array_adr)
309        ENDIF
[1762]310
311        return
312   END SUBROUTINE PMC_S_Set_DataArray_3d
313
314   SUBROUTINE PMC_S_setInd_and_AllocMem (ClientId)
315      IMPLICIT none
316      INTEGER,INTENT(IN)                      :: ClientId
317
318      INTEGER                                 :: i, istat, ierr
319      INTEGER                                 :: arlen, myIndex, tag
320      INTEGER                                 :: rCount                    ! count MPI requests
[1764]321      INTEGER(idp)                            :: bufsize                   ! Size of MPI data Window
[1762]322      TYPE(PeDef),POINTER                     :: aPE
323      TYPE(ArrayDef),POINTER                  :: ar
324      CHARACTER(len=DA_Namelen)               :: myName
325      INTEGER,DIMENSION(1024)                 :: req
326      Type(c_ptr)                             :: base_ptr
327      REAL(kind=wp),DIMENSION(:),POINTER      :: base_array
328      INTEGER(KIND=MPI_ADDRESS_KIND)          :: WinSize
329
330      myIndex = 1
331      rCount  = 0
332      bufsize = 8
333
334      ! First stride, Compute size and set index
335
336      do i=1,Clients(ClientId)%inter_npes
337         aPE => Clients(ClientId)%PEs(i)
338         tag = 200
339         do while (PMC_S_GetNextArray ( ClientId, myName,i))
340            ar  => aPE%Arrays
341            if(ar%dim_order == 2) then
342               arlen     = aPE%NrEle;                             ! 2D
343            else if(ar%dim_order == 33) then
344               arlen     = aPE%NrEle * ar%A_dim(4);               ! PALM 3D
345            else
346               arlen     = -1
347            end if
348            ar%BufIndex = myIndex
349
350            tag    = tag+1
351            rCount = rCount+1
352            CALL MPI_Isend (myIndex, 1, MPI_INTEGER, i-1, tag, Clients(ClientId)%inter_comm, req(rCount),ierr)
353
354            if(rCount == 1024) then                                  ! Maximum of 1024 outstanding requests
355               CALL MPI_Waitall (rCount, req, MPI_STATUSES_IGNORE, ierr)
356               rCount = 0;
357            end if
358
359            myIndex = myIndex+arlen
360            bufsize = bufsize+arlen
361            ar%BufSize = arlen
362
363         end do
364         if(rCount > 0) then                                          ! Wait for all send completed
365            CALL MPI_Waitall (rCount, req, MPI_STATUSES_IGNORE, ierr)
366         end if
367      end do
368
369      ! Create RMA (One Sided Communication) window for data buffer
370
371      CALL PMC_Alloc_mem (base_array, bufsize, base_ptr)
372      Clients(ClientId)%TotalBufferSize = bufsize*wp       !Total buffer size in Byte
373
374      WinSize = bufsize*wp
375!      write(9,*) 'PMC_S_SetInd_and_Mem ',m_model_rank,Clients(ClientId)%inter_npes,WinSize,bufsize
376      CALL MPI_Win_create (base_array, WinSize, wp, MPI_INFO_NULL, Clients(ClientId)%intra_comm, Clients(ClientId)%BufWin, ierr);
377      CALL MPI_Win_fence (0, Clients(ClientId)%BufWin, ierr);                    !  Open Window to set data
378      CALL MPI_Barrier(Clients(ClientId)%intra_comm, ierr)
379
380      ! Second stride, Set Buffer pointer
381
382      do i=1,Clients(ClientId)%inter_npes
383         aPE => Clients(ClientId)%PEs(i)
384         do while (PMC_S_GetNextArray ( ClientId, myName,i))
385            ar  => aPE%Arrays
[1764]386!--         TO_DO:  Adressrechnung ueberlegen?
[1762]387            ar%SendBuf = c_loc(base_array(ar%BufIndex))                         !kk Adressrechnung ueberlegen
388            if(ar%BufIndex+ar%BufSize > bufsize) then
[1764]389!--            TO_DO: can this error really happen, and what can be the reason?
[1762]390               write(0,'(a,i4,4i7,1x,a)') 'Buffer too small ',i,ar%BufIndex,ar%BufSize,ar%BufIndex+ar%BufSize,bufsize,trim(myName)
391               CALL MPI_Abort (MPI_COMM_WORLD, istat, ierr)
392            end if
393         end do
394      end do
395
396      return
397   END SUBROUTINE PMC_S_setInd_and_AllocMem
398
399   SUBROUTINE PMC_S_FillBuffer (ClientId, WaitTime)
400      IMPLICIT none
401      INTEGER,INTENT(IN)                      :: ClientId
402      REAL(kind=dp),INTENT(OUT),optional      :: WaitTime
403
404      !-- local variables
405      INTEGER                                 :: ip,ij,istat,ierr
406      INTEGER                                 :: myIndex
407      REAL(kind=dp)                           :: t1,t2
408      TYPE(PeDef),POINTER                     :: aPE
409      TYPE(ArrayDef),POINTER                  :: ar
410      CHARACTER(len=DA_Namelen)               :: myName
411      INTEGER,DIMENSION(1)                    :: buf_shape
412      REAL(kind=wp),POINTER,DIMENSION(:)      :: buf
413      REAL(kind=wp),POINTER,DIMENSION(:,:)    :: data_2d
414      REAL(kind=wp),POINTER,DIMENSION(:,:,:)  :: data_3d
415
416      t1 = PMC_Time()
417      CALL MPI_Barrier(Clients(ClientId)%intra_comm, ierr)              ! Wait for buffer empty
418      t2 = PMC_Time()
419      if(present(WaitTime)) WaitTime = t2-t1
420
421      do ip=1,Clients(ClientId)%inter_npes
422         aPE => Clients(ClientId)%PEs(ip)
423         do while (PMC_S_GetNextArray ( ClientId, myName,ip))
424            ar  => aPE%Arrays
425            myIndex=1
426            if(ar%dim_order == 2) then
427               buf_shape(1) = aPE%NrEle
428               CALL c_f_pointer(ar%SendBuf, buf, buf_shape)
429               CALL c_f_pointer(ar%data, data_2d, ar%A_dim(1:2))
430               do ij=1,aPE%NrEle
431                  buf(myIndex) = data_2d(aPE%locInd(ij)%j,aPE%locInd(ij)%i)
432                  myIndex = myIndex+1
433               end do
434            else if(ar%dim_order == 33) then
435               buf_shape(1) = aPE%NrEle*ar%A_dim(4)
436               CALL c_f_pointer(ar%SendBuf, buf, buf_shape)
437               CALL c_f_pointer(ar%data, data_3d, ar%A_dim(1:3))
438               do ij=1,aPE%NrEle
439                  buf(myIndex:myIndex+ar%A_dim(4)-1) = data_3d(1:ar%A_dim(4),aPE%locInd(ij)%j,aPE%locInd(ij)%i)
440                  myIndex = myIndex+ar%A_dim(4)
441               end do
442            else
[1764]443!--            TO_DO: can this error really happen, and what can be the reason?
[1762]444               write(0,*) "Illegal Order of Dimension ",ar%dim_order
445               CALL MPI_Abort (MPI_COMM_WORLD, istat, ierr);
446
447            end if
448          end do
449      end do
450
451      CALL MPI_Barrier(Clients(ClientId)%intra_comm, ierr)              ! buffer is full
452
453      return
454   END SUBROUTINE PMC_S_FillBuffer
455
456   SUBROUTINE PMC_S_GetData_from_Buffer (ClientId, WaitTime)
457      IMPLICIT none
458      INTEGER,INTENT(IN)                      :: ClientId
459      REAL(kind=dp),INTENT(OUT),optional      :: WaitTime
460
461      !-- local variables
462      INTEGER                                 :: ip,ij,istat,ierr
463      INTEGER                                 :: myIndex
464      REAL(kind=dp)                           :: t1,t2
465      TYPE(PeDef),POINTER                     :: aPE
466      TYPE(ArrayDef),POINTER                  :: ar
467      CHARACTER(len=DA_Namelen)               :: myName
468      INTEGER,DIMENSION(1)                    :: buf_shape
469      REAL(kind=wp),POINTER,DIMENSION(:)      :: buf
470      REAL(kind=wp),POINTER,DIMENSION(:,:)    :: data_2d
471      REAL(kind=wp),POINTER,DIMENSION(:,:,:)  :: data_3d
472
473      t1 = PMC_Time()
474      CALL MPI_Barrier(Clients(ClientId)%intra_comm, ierr)              ! Wait for MPI_Put from client
475      t2 = PMC_Time()
476      if(present(WaitTime)) WaitTime = t2-t1
477
478      do ip=1,Clients(ClientId)%inter_npes
479         aPE => Clients(ClientId)%PEs(ip)
480         do while (PMC_S_GetNextArray ( ClientId, myName,ip))
481            ar  => aPE%Arrays
482            myIndex=1
483            if(ar%dim_order == 2) then
484               buf_shape(1) = aPE%NrEle
485               CALL c_f_pointer(ar%SendBuf, buf, buf_shape)
486               CALL c_f_pointer(ar%data, data_2d, ar%A_dim(1:2))
487               do ij=1,aPE%NrEle
488                  data_2d(aPE%locInd(ij)%j,aPE%locInd(ij)%i) = buf(myIndex)
489                  myIndex = myIndex+1
490               end do
491            else if(ar%dim_order == 33) then
492               buf_shape(1) = aPE%NrEle*ar%A_dim(4)
493               CALL c_f_pointer(ar%SendBuf, buf, buf_shape)
494               CALL c_f_pointer(ar%data, data_3d, ar%A_dim(1:3))
495               do ij=1,aPE%NrEle
496                  data_3d(1:ar%A_dim(4),aPE%locInd(ij)%j,aPE%locInd(ij)%i) = buf(myIndex:myIndex+ar%A_dim(4)-1)
497                  myIndex = myIndex+ar%A_dim(4)
498               end do
499            else
[1764]500!--            TO_DO: can this error really happen, and what can be the reason?
[1762]501               write(0,*) "Illegal Order of Dimension ",ar%dim_order
502               CALL MPI_Abort (MPI_COMM_WORLD, istat, ierr);
503
504            end if
505          end do
506      end do
507
508      CALL MPI_Barrier(Clients(ClientId)%intra_comm, ierr)              ! data copy finished, buffer is free for use agein
509
510        return
511   END SUBROUTINE PMC_S_GetData_from_Buffer
512
513! Private SUBROUTINEs
514
515   SUBROUTINE Get_DA_names_from_client (ClientId)
516        IMPLICIT none
517        INTEGER,INTENT(IN)                    :: ClientId
518        !-- local variables
519        type(DA_NameDef)                      :: myName
520
521        !   Get Data Array Description and Name from Client
522
523        do
524            CALL PMC_Bcast ( myName%couple_index, 0, comm=m_to_client_comm(ClientId))
525            if(myName%couple_index == -1) EXIT
526            CALL PMC_Bcast ( myName%ServerDesc,   0, comm=m_to_client_comm(ClientId))
527            CALL PMC_Bcast ( myName%NameOnServer, 0, comm=m_to_client_comm(ClientId))
528            CALL PMC_Bcast ( myName%ClientDesc,   0, comm=m_to_client_comm(ClientId))
529            CALL PMC_Bcast ( myName%NameOnClient, 0, comm=m_to_client_comm(ClientId))
530
531            CALL PMC_G_SetName (clients(ClientID), myName%couple_index, myName%NameOnServer )
532        end do
533
534        return
535   END SUBROUTINE Get_DA_names_from_client
536
[1766]537   SUBROUTINE PMC_S_SetArray (ClientId, NrDims, dims, dim_order, array_adr, second_adr)
[1762]538      IMPLICIT none
539
540      INTEGER,INTENT(IN)                      :: ClientId
541      INTEGER,INTENT(IN)                      :: NrDims
542      INTEGER,INTENT(IN),DIMENSION(:)         :: dims
543      INTEGER,INTENT(IN)                      :: dim_order
544      TYPE(c_ptr),INTENT(IN)                  :: array_adr
[1766]545      TYPE(c_ptr),INTENT(IN),OPTIONAL         :: second_adr
[1762]546
547      INTEGER                                 :: i
548      TYPE(PeDef),POINTER                     :: aPE
549      TYPE(ArrayDef),POINTER                  :: ar
550      CHARACTER(len=DA_Namelen)               :: myName
551
552      !  Set Array for Client interPE 0
553
554       do i=1,Clients(ClientId)%inter_npes
555          aPE => Clients(ClientId)%PEs(i)
556          ar  => aPE%Arrays
557          ar%NrDims    = NrDims
558          ar%A_dim     = dims
559          ar%dim_order = dim_order
560          ar%data      = array_adr
[1766]561          if(present(second_adr)) then
562             ar%po_data(1) = array_adr
563             ar%po_data(2) = second_adr
564          else
565             ar%po_data(1) = C_NULL_PTR
566             ar%po_data(2) = C_NULL_PTR
567          end if
[1762]568       end do
569
570      return
571   END SUBROUTINE PMC_S_SetArray
572
573
[1766]574   SUBROUTINE PMC_S_Set_Active_data_array (ClientId,iactive)
575      IMPLICIT none
576
577      INTEGER,INTENT(IN)                      :: ClientId
578      INTEGER,INTENT(IN)                      :: iactive
579
580!--   local variables
581      INTEGER                                 :: i, ip
582      TYPE(PeDef),POINTER                     :: aPE
583      TYPE(ArrayDef),POINTER                  :: ar
584      CHARACTER(len=DA_Namelen)               :: myName
585
586      do ip=1,Clients(ClientId)%inter_npes
587         aPE => Clients(ClientId)%PEs(ip)
588         do while (PMC_S_GetNextArray ( ClientId, myName,ip))
589            ar  => aPE%Arrays
590            if(iactive == 1 .OR. iactive == 2)   then
591               ar%data = ar%po_data(iactive)
592            end if
593         end do
594      end do
595
596      return
597   END SUBROUTINE PMC_S_Set_Active_data_array
598
599
[1762]600    SUBROUTINE Set_PE_index_list (ClientId, myClient,index_list,NrP)
601       IMPLICIT none
602
603       INTEGER,INTENT(IN)                      :: ClientId
604       TYPE(ClientDef),INTENT(INOUT)           :: myClient
605       INTEGER,INTENT(IN),DIMENSION(:,:)       :: index_list
606       INTEGER,INTENT(IN)                      :: NrP
607
608!--    local variables
609       INTEGER                                 :: i,j,ind,ierr,i2
610       TYPE(PeDef),POINTER                     :: aPE
611       INTEGER                                 :: RemPE
612       INTEGER,DIMENSION(myClient%inter_npes)  :: RemInd
613       INTEGER,DIMENSION(:),POINTER            :: RemIndw
614       INTEGER,DIMENSION(:),POINTER            :: RLdef
615       INTEGER(KIND=MPI_ADDRESS_KIND)          :: WinSize
616       INTEGER                                 :: indWin,indWin2
617
618       ! First, count entries for every remote client PE
619
620       do i=1,myClient%inter_npes
621          aPE => myClient%PEs(i)
622          aPE%NrEle = 0
623       end do
624
625       do j=1,NrP                                ! loop over number of cells coarse grid
626          RemPE = index_list(5,j)+1              ! Pe number remote PE
627          aPE => myClient%PEs(RemPE)
628          aPE% NrEle = aPE% NrEle+1              ! Increment Number of elements for this client Pe
629       end do
630
631       do i=1,myClient%inter_npes
632          aPE => myClient%PEs(i)
633          ALLOCATE(aPE%locInd(aPE%NrEle))
634       end do
635
636       RemInd = 0
637
638       ! Second, Create lists
639
640       do j=1,NrP                                ! loop over number of cells coarse grid
641          RemPE = index_list(5,j)+1              ! Pe number remote PE
642          aPE => myClient%PEs(RemPE)
643          RemInd(RemPE)     = RemInd(RemPE)+1
644          ind               = RemInd(RemPE)
645          aPE%locInd(ind)%i = index_list(1,j)
646          aPE%locInd(ind)%j = index_list(2,j)
647       end do
648
649       !  Prepare Number of Elements for Client PEs
650       CALL PMC_Alloc_mem (RLdef, myClient%inter_npes*2)
651       WinSize = myClient%inter_npes*c_sizeof(i)*2   ! Number of Client PEs * size of INTEGER (i just arbitrary INTEGER)
652
653       CALL MPI_Win_create (RLdef, WinSize, iwp, MPI_INFO_NULL, myClient%intra_comm, indWin, ierr);
654       CALL MPI_Win_fence (0, indWin, ierr);         !  Open Window to set data
655
656       RLdef(1) = 0                                  ! Index on Remote PE 0
657       RLdef(2) = RemInd(1)                          ! Number of Elements on Rem PE 0
658
659       do i=2,myClient%inter_npes                    ! Reserve Buffer for index array
660          i2          = (i-1)*2+1
661          RLdef(i2)   = RLdef(i2-2) + RLdef(i2-1)*2  ! Index on Remote PE
662          RLdef(i2+1) = RemInd(i)                    ! Number of Elements on Remote PE
663       end do
664
665       CALL MPI_Win_fence (0, indWin, ierr);         ! Close Window to allow client to access data
666       CALL MPI_Win_fence (0, indWin, ierr);         ! Client has retrieved data
667
668       i2 = 2*myClient%inter_npes-1
669       WinSize = (RLdef(i2)+RLdef(i2+1))*2
670       WinSize = max(WinSize,1)                      ! Make sure, MPI_Alloc_mem works
671
672       CALL PMC_Alloc_mem (RemIndw, int(WinSize))
673
674       CALL MPI_Barrier (m_model_comm, ierr)
675       CALL MPI_Win_create (RemIndw, WinSize*c_sizeof(i), iwp, MPI_INFO_NULL, myClient%intra_comm, indWin2, ierr);
676
677       CALL MPI_Win_fence (0, indWin2, ierr);         !  Open Window to set data
678       do j=1,NrP                                ! this loop creates the 2D index list
679          RemPE = index_list(5,j)+1              ! Pe number remote PE
680          aPE => myClient%PEs(RemPE)
681          i2    = RemPE*2-1
682          ind   = RLdef(i2)+1
683          RemIndw(ind)   = index_list(3,j)
684          RemIndw(ind+1) = index_list(4,j)
685          RLdef(i2) = RLdef(i2)+2
686       end do
687       CALL MPI_Win_fence (0, indWin2, ierr);      !all data set
688
689       CALL MPI_Barrier(myClient%intra_comm, ierr) ! Dont know why, but this barrier is necessary before we can free the windows
690
691       CALL MPI_Win_free(indWin, ierr)
692       CALL MPI_Win_free(indWin2, ierr)
693
694!      Sollte funktionieren, Problem mit MPI implementation
695!      https://www.lrz.de/services/software/parallel/mpi/onesided
696!       CALL MPI_Free_mem (RemIndw, ierr)
697
698       return
699    END SUBROUTINE Set_PE_index_list
700
[1764]701#endif
[1762]702END MODULE pmc_server
Note: See TracBrowser for help on using the repository browser.