New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
modmpp.F in branches/UKMO/dev_r5518_AMM15_package/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES – NEMO

source: branches/UKMO/dev_r5518_AMM15_package/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES/modmpp.F @ 10251

Last change on this file since 10251 was 10251, checked in by kingr, 5 years ago

Rolled back to r10247 - i.e., undid merge of pkg br and 3.6_stable br

File size: 38.5 KB
Line 
1!
2! $Id: modmpp.F 2731 2011-04-08 12:05:35Z rblod $
3!
4C     AGRIF (Adaptive Grid Refinement In Fortran)
5C
6C     Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr)
7C                        Christophe Vouland (Christophe.Vouland@imag.fr)
8C
9C     This program is free software; you can redistribute it and/or modify
10C     it under the terms of the GNU General Public License as published by
11C     the Free Software Foundation; either version 2 of the License, or
12C     (at your option) any later version.
13C
14C     This program is distributed in the hope that it will be useful,
15C     but WITHOUT ANY WARRANTY; without even the implied warranty of
16C     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
17C     GNU General Public License for more details.
18C
19C     You should have received a copy of the GNU General Public License
20C     along with this program; if not, write to the Free Software
21C     Foundation, Inc., 59 Temple Place- Suite 330, Boston, MA 02111-1307, USA.
22C
23C
24C
25CCC   Module Agrif_mpp
26C
27      Module Agrif_mpp
28      Use Agrif_Types
29      Use Agrif_Arrays
30
31      Contains
32#ifdef key_mpp_mpi
33      Subroutine Get_External_Data_first(pttruetab,
34     &   cetruetab,pttruetabwhole,cetruetabwhole,nbdim,memberin,
35     &   memberout,memberoutall1,sendtoproc,recvfromproc,imin,imax,
36     &   imin_recv,imax_recv)
37
38      IMPLICIT NONE
39      INCLUDE 'mpif.h'
40      INTEGER :: nbdim
41      INTEGER,DIMENSION(nbdim,0:Agrif_NbProcs-1)    :: pttruetab,
42     &                                                 cetruetab
43      INTEGER,DIMENSION(nbdim,0:Agrif_NbProcs-1)    :: pttruetabwhole,
44     &                                                 cetruetabwhole
45      INTEGER :: k,i,k2
46      LOGICAL, DIMENSION(0:Agrif_Nbprocs-1) :: sendtoproc, recvfromproc
47      INTEGER,DIMENSION(nbdim,0:Agrif_NbProcs-1):: imin,imax,
48     &   imin_recv,imax_recv
49      LOGICAL :: memberin, memberout
50      INTEGER :: imintmp, imaxtmp,j,i1
51      INTEGER :: imin1,imax1
52      LOGICAL :: tochange,tochangebis
53      INTEGER,DIMENSION(nbdim,0:Agrif_NbProcs-1)    :: pttruetab2,
54     &                                                 cetruetab2
55      LOGICAL :: memberout1(1),memberoutall(0:Agrif_Nbprocs-1)
56      LOGICAL, OPTIONAL :: memberoutall1(0:Agrif_Nbprocs-1)
57      INTEGER :: code
58     
59C pttruetab2 and cetruetab2 are modified arraysin order to always
60C send the most inner points
61
62       
63        IF (present(memberoutall1)) THEN
64        memberoutall = memberoutall1
65        ELSE
66         memberout1(1) = memberout
67
68         CALL MPI_ALLGATHER(memberout1,1,MPI_LOGICAL,memberoutall,
69     &                  1,MPI_LOGICAL,MPI_COMM_AGRIF,code)
70        ENDIF
71         pttruetab2(:,Agrif_Procrank) = pttruetab(:,Agrif_Procrank)
72         cetruetab2(:,Agrif_Procrank) = cetruetab(:,Agrif_Procrank)
73         do k2=0,Agrif_Nbprocs-1
74            do i=1,nbdim
75
76           tochangebis=.TRUE.
77           DO i1=1,nbdim
78            IF (i .NE. i1) THEN
79              IF ((pttruetab(i1,Agrif_Procrank).NE.pttruetab(i1,k2)).OR.
80     &          (cetruetab(i1,Agrif_Procrank).NE.cetruetab(i1,k2))) THEN
81                   tochangebis = .FALSE.
82               EXIT
83              ENDIF
84             ENDIF
85           ENDDO
86
87           IF (tochangebis) THEN
88
89
90          imin1 = max(pttruetab(i,Agrif_Procrank),
91     &                    pttruetab(i,k2))
92          imax1 = min(cetruetab(i,Agrif_Procrank),
93     &                    cetruetab(i,k2))
94
95C Always send the most interior points
96
97           tochange = .false.
98           IF (cetruetab(i,Agrif_Procrank)> cetruetab(i,k2)) THEN
99
100           DO j=imin1,imax1
101             IF ((cetruetab(i,k2)-j) >
102     &             (j-pttruetab(i,Agrif_Procrank))) THEN
103             imintmp = j+1
104             tochange = .TRUE.
105             ELSE
106              EXIT
107             ENDIF
108          ENDDO
109          ENDIF
110
111           if (tochange) then
112C
113              pttruetab2(i,Agrif_Procrank) = imintmp
114C
115          endif
116
117           tochange = .FALSE.
118           imaxtmp=0
119           IF (pttruetab(i,Agrif_Procrank) < pttruetab(i,k2)) THEN
120          DO j=imax1,imin1,-1
121            IF ((j-pttruetab(i,k2)) >
122     &             (cetruetab(i,Agrif_Procrank)-j)) THEN
123             imaxtmp = j-1
124             tochange = .TRUE.
125            ELSE
126             EXIT
127            ENDIF
128          ENDDO
129
130          ENDIF
131
132                    if (tochange) then
133C
134              cetruetab2(i,Agrif_Procrank) = imaxtmp
135C
136          endif
137C
138
139          ENDIF
140           enddo
141         enddo
142
143
144       do k = 0,Agrif_NbProcs-1
145C
146        sendtoproc(k) = .true.
147C
148!CDIR SHORTLOOP
149        do i = 1,nbdim
150C
151          imin(i,k) = max(pttruetab2(i,Agrif_Procrank),
152     &                    pttruetabwhole(i,k))
153          imax(i,k) = min(cetruetab2(i,Agrif_Procrank),
154     &                    cetruetabwhole(i,k))
155C
156          if (imin(i,k) > imax(i,k)) then
157C
158              sendtoproc(k) = .false.
159C
160          endif
161C
162        enddo
163        IF (.NOT.memberoutall(k)) THEN
164           sendtoproc(k) = .FALSE.
165        ENDIF
166C
167      enddo
168     
169      Call Exchangesamelevel_first(sendtoproc,nbdim,imin,imax,
170     &     recvfromproc,imin_recv,imax_recv)
171     
172      End Subroutine Get_External_Data_first
173C
174      Subroutine Get_External_Data(tempC,tempCextend,pttruetab,
175     &   cetruetab,pttruetabwhole,cetruetabwhole,nbdim,memberin,
176     &   memberout,memberoutall1)
177
178      IMPLICIT NONE
179      INCLUDE 'mpif.h'
180      INTEGER :: nbdim
181      TYPE(Agrif_PVariable) :: tempC, tempCextend
182      INTEGER,DIMENSION(nbdim,0:Agrif_NbProcs-1)    :: pttruetab,
183     &                                                 cetruetab
184      INTEGER,DIMENSION(nbdim,0:Agrif_NbProcs-1)    :: pttruetabwhole,
185     &                                                 cetruetabwhole
186      INTEGER :: k,i,k2
187      LOGICAL :: sendtoproc(0:Agrif_Nbprocs-1)
188      INTEGER,DIMENSION(nbdim,0:Agrif_NbProcs-1)    :: imin,imax
189      LOGICAL :: memberin, memberout
190      INTEGER :: imintmp, imaxtmp,j,i1
191      INTEGER :: imin1,imax1
192      LOGICAL :: tochange,tochangebis
193      INTEGER,DIMENSION(nbdim,0:Agrif_NbProcs-1)    :: pttruetab2,
194     &                                                 cetruetab2
195      LOGICAL :: memberout1(1),memberoutall(0:Agrif_Nbprocs-1)
196      LOGICAL, OPTIONAL :: memberoutall1(0:Agrif_Nbprocs-1)
197      INTEGER :: code
198
199C pttruetab2 and cetruetab2 are modified arraysin order to always
200C send the most inner points
201
202       
203        IF (present(memberoutall1)) THEN
204        memberoutall = memberoutall1
205        ELSE
206         memberout1(1) = memberout
207
208         CALL MPI_ALLGATHER(memberout1,1,MPI_LOGICAL,memberoutall,
209     &                  1,MPI_LOGICAL,MPI_COMM_AGRIF,code)
210        ENDIF
211         pttruetab2(:,Agrif_Procrank) = pttruetab(:,Agrif_Procrank)
212         cetruetab2(:,Agrif_Procrank) = cetruetab(:,Agrif_Procrank)
213         do k2=0,Agrif_Nbprocs-1
214            do i=1,nbdim
215
216           tochangebis=.TRUE.
217           DO i1=1,nbdim
218            IF (i .NE. i1) THEN
219              IF ((pttruetab(i1,Agrif_Procrank).NE.pttruetab(i1,k2)).OR.
220     &          (cetruetab(i1,Agrif_Procrank).NE.cetruetab(i1,k2))) THEN
221                   tochangebis = .FALSE.
222               EXIT
223              ENDIF
224             ENDIF
225           ENDDO
226
227           IF (tochangebis) THEN
228
229
230          imin1 = max(pttruetab(i,Agrif_Procrank),
231     &                    pttruetab(i,k2))
232          imax1 = min(cetruetab(i,Agrif_Procrank),
233     &                    cetruetab(i,k2))
234
235C Always send the most interior points
236
237           tochange = .false.
238           IF (cetruetab(i,Agrif_Procrank)> cetruetab(i,k2)) THEN
239
240           DO j=imin1,imax1
241             IF ((cetruetab(i,k2)-j) >
242     &             (j-pttruetab(i,Agrif_Procrank))) THEN
243             imintmp = j+1
244             tochange = .TRUE.
245             ELSE
246              EXIT
247             ENDIF
248          ENDDO
249          ENDIF
250
251           if (tochange) then
252C
253              pttruetab2(i,Agrif_Procrank) = imintmp
254C
255          endif
256
257           tochange = .FALSE.
258           imaxtmp=0
259           IF (pttruetab(i,Agrif_Procrank) < pttruetab(i,k2)) THEN
260          DO j=imax1,imin1,-1
261            IF ((j-pttruetab(i,k2)) >
262     &             (cetruetab(i,Agrif_Procrank)-j)) THEN
263             imaxtmp = j-1
264             tochange = .TRUE.
265            ELSE
266             EXIT
267            ENDIF
268          ENDDO
269
270          ENDIF
271
272                    if (tochange) then
273C
274              cetruetab2(i,Agrif_Procrank) = imaxtmp
275C
276          endif
277C
278
279          ENDIF
280           enddo
281         enddo
282
283
284       do k = 0,Agrif_NbProcs-1
285C
286        sendtoproc(k) = .true.
287C
288!CDIR SHORTLOOP
289        do i = 1,nbdim
290C
291          imin(i,k) = max(pttruetab2(i,Agrif_Procrank),
292     &                    pttruetabwhole(i,k))
293          imax(i,k) = min(cetruetab2(i,Agrif_Procrank),
294     &                    cetruetabwhole(i,k))
295C
296          if (imin(i,k) > imax(i,k)) then
297C
298              sendtoproc(k) = .false.
299C
300          endif
301C
302        enddo
303        IF (.NOT.memberoutall(k)) THEN
304           sendtoproc(k) = .FALSE.
305        ENDIF
306C
307      enddo
308
309
310c       IF (.NOT.memberin) sendtoproc = .FALSE.
311
312      IF (memberout) THEN
313      Call Agrif_nbdim_allocation(tempCextend%var,
314     &                 pttruetabwhole(:,Agrif_ProcRank),
315     &                 cetruetabwhole(:,Agrif_ProcRank),nbdim)
316      call Agrif_nbdim_Full_VarEQreal(tempCextend%var,0.,nbdim)
317      ENDIF
318
319      IF (sendtoproc(Agrif_ProcRank)) THEN
320           Call Agrif_nbdim_VarEQvar(tempCextend%var,
321     &                               imin(:,Agrif_Procrank),
322     &                               imax(:,Agrif_Procrank),
323     &                               tempC%var,
324     &                               imin(:,Agrif_Procrank),
325     &                               imax(:,Agrif_Procrank),
326     &                               nbdim)
327      ENDIF
328
329      Call Exchangesamelevel(sendtoproc,nbdim,imin,imax,tempC,
330     &     tempCextend)
331
332      End Subroutine Get_External_Data
333
334       Subroutine ExchangeSameLevel(sendtoproc,nbdim,imin,imax,
335     &          tempC,tempCextend)
336      Implicit None
337      INTEGER :: nbdim
338      INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1) :: imin,imax
339      INTEGER,DIMENSION(nbdim,2,0:Agrif_Nbprocs-1) :: iminmax_temp
340      INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1) :: imin_recv,imax_recv
341      TYPE(Agrif_PVARIABLE) :: tempC,tempCextend
342      LOGICAL,DIMENSION(0:Agrif_Nbprocs-1)       :: sendtoproc
343      LOGICAL,DIMENSION(0:Agrif_Nbprocs-1)       :: recvfromproc
344      LOGICAL                                    :: res
345      TYPE(AGRIF_PVARIABLE), SAVE                      :: temprecv
346
347      INCLUDE 'mpif.h'
348          INTEGER :: i,k
349          INTEGER :: etiquette = 100
350          INTEGER :: code, datasize
351          INTEGER,DIMENSION(MPI_STATUS_SIZE)   :: statut
352
353
354      do k = 0,Agrif_ProcRank-1
355C
356C
357            Call MPI_SEND(sendtoproc(k),1,MPI_LOGICAL,k,etiquette,
358     &                    MPI_COMM_AGRIF,code)
359C
360            if (sendtoproc(k)) then
361C
362                iminmax_temp(:,1,k) = imin(:,k)
363                iminmax_temp(:,2,k) = imax(:,k)
364
365                Call MPI_SEND(iminmax_temp(:,:,k),
366     &                        2*nbdim,MPI_INTEGER,k,etiquette,
367     &                        MPI_COMM_AGRIF,code)
368C
369                datasize = 1
370C
371!CDIR SHORTLOOP
372                do i = 1,nbdim
373C
374                  datasize = datasize * (imax(i,k)-imin(i,k)+1)
375C
376                enddo
377C
378                SELECT CASE(nbdim)
379                CASE(1)
380                   Call MPI_SEND(tempC%var%array1(
381     &                        imin(1,k):imax(1,k)),
382     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
383     &                        MPI_COMM_AGRIF,code)
384                CASE(2)
385                   Call MPI_SEND(tempC%var%array2(
386     &                        imin(1,k):imax(1,k),
387     &                        imin(2,k):imax(2,k)),
388     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
389     &                        MPI_COMM_AGRIF,code)
390                CASE(3)
391               
392                  Call Agrif_Send_3Darray(tempC%var%array3,
393     &             lbound(tempC%var%array3),imin(:,k),imax(:,k),k)
394                CASE(4)
395                   Call MPI_SEND(tempC%var%array4(
396     &                        imin(1,k):imax(1,k),
397     &                        imin(2,k):imax(2,k),
398     &                        imin(3,k):imax(3,k),
399     &                        imin(4,k):imax(4,k)),
400     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
401     &                        MPI_COMM_AGRIF,code)
402                CASE(5)
403                   Call MPI_SEND(tempC%var%array5(
404     &                        imin(1,k):imax(1,k),
405     &                        imin(2,k):imax(2,k),
406     &                        imin(3,k):imax(3,k),
407     &                        imin(4,k):imax(4,k),
408     &                        imin(5,k):imax(5,k)),
409     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
410     &                        MPI_COMM_AGRIF,code)
411                CASE(6)
412                   Call MPI_SEND(tempC%var%array6(
413     &                        imin(1,k):imax(1,k),
414     &                        imin(2,k):imax(2,k),
415     &                        imin(3,k):imax(3,k),
416     &                        imin(4,k):imax(4,k),
417     &                        imin(5,k):imax(5,k),
418     &                        imin(6,k):imax(6,k)),
419     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
420     &                        MPI_COMM_AGRIF,code)
421                END SELECT
422C
423            endif
424
425C
426      enddo
427C
428C
429C     Reception from others processors of the necessary part of the parent grid
430      do k = Agrif_ProcRank+1,Agrif_Nbprocs-1
431C
432C
433            Call MPI_RECV(res,1,MPI_LOGICAL,k,etiquette,
434     &                    MPI_COMM_AGRIF,statut,code)
435C
436            recvfromproc(k) = res
437
438C
439            if (recvfromproc(k)) then
440C
441                Call MPI_RECV(iminmax_temp(:,:,k),
442     &                        2*nbdim,MPI_INTEGER,k,etiquette,
443     &                        MPI_COMM_AGRIF,statut,code)
444
445                imin_recv(:,k) = iminmax_temp(:,1,k)
446                imax_recv(:,k) = iminmax_temp(:,2,k)
447
448                datasize = 1
449C
450!CDIR SHORTLOOP
451                do i = 1,nbdim
452C
453                datasize = datasize * (imax_recv(i,k)-imin_recv(i,k)+1)
454C
455                enddo
456
457             IF (.Not.Associated(temprecv%var)) allocate(temprecv%var)
458             call Agrif_nbdim_allocation(temprecv%var,imin_recv(:,k),
459     &   imax_recv(:,k),nbdim)
460            SELECT CASE(nbdim)
461            CASE(1)
462              Call MPI_RECV(temprecv%var%array1,
463     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
464     &               MPI_COMM_AGRIF,statut,code)
465            CASE(2)
466              Call MPI_RECV(temprecv%var%array2,
467     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
468     &               MPI_COMM_AGRIF,statut,code)
469            CASE(3)
470              Call MPI_RECV(temprecv%var%array3,
471     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
472     &               MPI_COMM_AGRIF,statut,code)
473
474            CASE(4)
475              Call MPI_RECV(temprecv%var%array4,
476     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
477     &               MPI_COMM_AGRIF,statut,code)
478            CASE(5)
479              Call MPI_RECV(temprecv%var%array5,
480     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
481     &               MPI_COMM_AGRIF,statut,code)
482            CASE(6)
483              Call MPI_RECV(temprecv%var%array6,
484     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
485     &               MPI_COMM_AGRIF,statut,code)
486       END SELECT
487                       
488            Call where_valtabtotab_mpi(tempCextend%var,
489     &             temprecv%var,imin_recv(:,k),imax_recv(:,k),0.,nbdim)
490     
491                Call Agrif_nbdim_deallocation(temprecv%var,nbdim)
492C                deallocate(temprecv%var)
493
494            endif
495
496C
497      enddo
498
499C     Reception from others processors of the necessary part of the parent grid
500      do k = Agrif_ProcRank+1,Agrif_Nbprocs-1
501C
502C
503           
504            Call MPI_SEND(sendtoproc(k),1,MPI_LOGICAL,k,etiquette,
505     &                    MPI_COMM_AGRIF,code)
506C
507            if (sendtoproc(k)) then
508C
509                iminmax_temp(:,1,k) = imin(:,k)
510                iminmax_temp(:,2,k) = imax(:,k)
511
512                Call MPI_SEND(iminmax_temp(:,:,k),
513     &                        2*nbdim,MPI_INTEGER,k,etiquette,
514     &                        MPI_COMM_AGRIF,code)
515C
516                SELECT CASE(nbdim)
517                CASE(1)
518                datasize=SIZE(tempC%var%array1(
519     &                        imin(1,k):imax(1,k)))
520                   Call MPI_SEND(tempC%var%array1(
521     &                        imin(1,k):imax(1,k)),
522     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
523     &                        MPI_COMM_AGRIF,code)
524                CASE(2)
525                datasize=SIZE(tempC%var%array2(
526     &                        imin(1,k):imax(1,k),
527     &                        imin(2,k):imax(2,k)))
528                   Call MPI_SEND(tempC%var%array2(
529     &                        imin(1,k):imax(1,k),
530     &                        imin(2,k):imax(2,k)),
531     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
532     &                        MPI_COMM_AGRIF,code)
533                CASE(3)
534                datasize=SIZE(tempC%var%array3(
535     &                        imin(1,k):imax(1,k),
536     &                        imin(2,k):imax(2,k),
537     &                        imin(3,k):imax(3,k)))
538                   Call MPI_SEND(tempC%var%array3(
539     &                        imin(1,k):imax(1,k),
540     &                        imin(2,k):imax(2,k),
541     &                        imin(3,k):imax(3,k)),
542     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
543     &                        MPI_COMM_AGRIF,code)
544                CASE(4)
545                datasize=SIZE(tempC%var%array4(
546     &                        imin(1,k):imax(1,k),
547     &                        imin(2,k):imax(2,k),
548     &                        imin(3,k):imax(3,k),
549     &                        imin(4,k):imax(4,k)))
550                   Call MPI_SEND(tempC%var%array4(
551     &                        imin(1,k):imax(1,k),
552     &                        imin(2,k):imax(2,k),
553     &                        imin(3,k):imax(3,k),
554     &                        imin(4,k):imax(4,k)),
555     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
556     &                        MPI_COMM_AGRIF,code)
557                CASE(5)
558                datasize=SIZE(tempC%var%array5(
559     &                        imin(1,k):imax(1,k),
560     &                        imin(2,k):imax(2,k),
561     &                        imin(3,k):imax(3,k),
562     &                        imin(4,k):imax(4,k),
563     &                        imin(5,k):imax(5,k)))
564                   Call MPI_SEND(tempC%var%array5(
565     &                        imin(1,k):imax(1,k),
566     &                        imin(2,k):imax(2,k),
567     &                        imin(3,k):imax(3,k),
568     &                        imin(4,k):imax(4,k),
569     &                        imin(5,k):imax(5,k)),
570     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
571     &                        MPI_COMM_AGRIF,code)
572                CASE(6)
573                datasize=SIZE(tempC%var%array6(
574     &                        imin(1,k):imax(1,k),
575     &                        imin(2,k):imax(2,k),
576     &                        imin(3,k):imax(3,k),
577     &                        imin(4,k):imax(4,k),
578     &                        imin(5,k):imax(5,k),
579     &                        imin(6,k):imax(6,k)))
580                   Call MPI_SEND(tempC%var%array6(
581     &                        imin(1,k):imax(1,k),
582     &                        imin(2,k):imax(2,k),
583     &                        imin(3,k):imax(3,k),
584     &                        imin(4,k):imax(4,k),
585     &                        imin(5,k):imax(5,k),
586     &                        imin(6,k):imax(6,k)),
587     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
588     &                        MPI_COMM_AGRIF,code)
589                END SELECT
590C
591            endif
592
593C
594      enddo
595C
596C
597C     Reception from others processors of the necessary part of the parent grid
598      do k = Agrif_ProcRank-1,0,-1
599C
600C
601            Call MPI_RECV(res,1,MPI_LOGICAL,k,etiquette,
602     &                    MPI_COMM_AGRIF,statut,code)
603C
604            recvfromproc(k) = res
605
606C
607            if (recvfromproc(k)) then
608C
609                Call MPI_RECV(iminmax_temp(:,:,k),
610     &                        2*nbdim,MPI_INTEGER,k,etiquette,
611     &                        MPI_COMM_AGRIF,statut,code)
612
613C                imin_recv(:,k) = iminmax_temp(:,1,k)
614C                imax_recv(:,k) = iminmax_temp(:,2,k)
615
616C                datasize = 1
617C
618C                do i = 1,nbdim
619C
620C                datasize = datasize * (imax_recv(i,k)-imin_recv(i,k)+1)
621C
622C                enddo
623             IF (.Not.Associated(temprecv%var)) allocate(temprecv%var)
624             call Agrif_nbdim_allocation(temprecv%var,
625     &   iminmax_temp(:,1,k),iminmax_temp(:,2,k),nbdim)
626            SELECT CASE(nbdim)
627            CASE(1)
628              datasize=SIZE(temprecv%var%array1)
629              Call MPI_RECV(temprecv%var%array1,
630     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
631     &               MPI_COMM_AGRIF,statut,code)
632            CASE(2)
633              datasize=SIZE(temprecv%var%array2)
634              Call MPI_RECV(temprecv%var%array2,
635     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
636     &               MPI_COMM_AGRIF,statut,code)
637            CASE(3)
638              datasize=SIZE(temprecv%var%array3)
639              Call MPI_RECV(temprecv%var%array3,
640     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
641     &               MPI_COMM_AGRIF,statut,code)
642
643            CASE(4)
644              datasize=SIZE(temprecv%var%array4)
645              Call MPI_RECV(temprecv%var%array4,
646     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
647     &               MPI_COMM_AGRIF,statut,code)
648            CASE(5)
649              datasize=SIZE(temprecv%var%array5)
650              Call MPI_RECV(temprecv%var%array5,
651     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
652     &               MPI_COMM_AGRIF,statut,code)
653            CASE(6)
654              datasize=SIZE(temprecv%var%array6)
655              Call MPI_RECV(temprecv%var%array6,
656     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
657     &               MPI_COMM_AGRIF,statut,code)
658          END SELECT
659           
660            Call where_valtabtotab_mpi(tempCextend%var,
661     &             temprecv%var,iminmax_temp(:,1,k),iminmax_temp(:,2,k)
662     &            ,0.,nbdim)
663     
664                Call Agrif_nbdim_deallocation(temprecv%var,nbdim)
665C                deallocate(temprecv%var)
666            endif
667
668C
669      enddo
670
671          End Subroutine ExchangeSamelevel
672         
673       Subroutine ExchangeSameLevel_first(sendtoproc,nbdim,imin,imax,
674     &          recvfromproc,imin_recv,imax_recv)
675      Implicit None
676      INTEGER :: nbdim
677      INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1) :: imin,imax
678      INTEGER,DIMENSION(nbdim,2,0:Agrif_Nbprocs-1) :: iminmax_temp
679      INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1) :: imin_recv,imax_recv
680      LOGICAL,DIMENSION(0:Agrif_Nbprocs-1)       :: sendtoproc
681      LOGICAL,DIMENSION(0:Agrif_Nbprocs-1)       :: recvfromproc
682      LOGICAL                                    :: res
683
684      INCLUDE 'mpif.h'
685          INTEGER :: i,k
686          INTEGER :: etiquette = 100
687          INTEGER :: code, datasize
688          INTEGER,DIMENSION(MPI_STATUS_SIZE)   :: statut
689
690     
691      do k = 0,Agrif_ProcRank-1
692C
693C
694            Call MPI_SEND(sendtoproc(k),1,MPI_LOGICAL,k,etiquette,
695     &                    MPI_COMM_AGRIF,code)
696C
697            if (sendtoproc(k)) then
698C
699                iminmax_temp(:,1,k) = imin(:,k)
700                iminmax_temp(:,2,k) = imax(:,k)
701
702                Call MPI_SEND(iminmax_temp(:,:,k),
703     &                        2*nbdim,MPI_INTEGER,k,etiquette,
704     &                        MPI_COMM_AGRIF,code)
705C
706            endif
707
708C
709      enddo
710C
711C
712C     Reception from others processors of the necessary part of the parent grid
713      do k = Agrif_ProcRank+1,Agrif_Nbprocs-1
714C
715C
716            Call MPI_RECV(res,1,MPI_LOGICAL,k,etiquette,
717     &                    MPI_COMM_AGRIF,statut,code)
718C
719            recvfromproc(k) = res
720
721C
722            if (recvfromproc(k)) then
723C
724                Call MPI_RECV(iminmax_temp(:,:,k),
725     &                        2*nbdim,MPI_INTEGER,k,etiquette,
726     &                        MPI_COMM_AGRIF,statut,code)
727
728                imin_recv(:,k) = iminmax_temp(:,1,k)
729                imax_recv(:,k) = iminmax_temp(:,2,k)
730            endif
731
732C
733      enddo
734
735C     Reception from others processors of the necessary part of the parent grid
736      do k = Agrif_ProcRank+1,Agrif_Nbprocs-1
737C
738C
739           
740            Call MPI_SEND(sendtoproc(k),1,MPI_LOGICAL,k,etiquette,
741     &                    MPI_COMM_AGRIF,code)
742C
743            if (sendtoproc(k)) then
744C
745                iminmax_temp(:,1,k) = imin(:,k)
746                iminmax_temp(:,2,k) = imax(:,k)
747
748                Call MPI_SEND(iminmax_temp(:,:,k),
749     &                        2*nbdim,MPI_INTEGER,k,etiquette,
750     &                        MPI_COMM_AGRIF,code)
751C
752            endif
753
754C
755      enddo
756C
757C
758C     Reception from others processors of the necessary part of the parent grid
759      do k = Agrif_ProcRank-1,0,-1
760C
761C
762            Call MPI_RECV(res,1,MPI_LOGICAL,k,etiquette,
763     &                    MPI_COMM_AGRIF,statut,code)
764C
765            recvfromproc(k) = res
766
767C
768            if (recvfromproc(k)) then
769C
770                Call MPI_RECV(iminmax_temp(:,:,k),
771     &                        2*nbdim,MPI_INTEGER,k,etiquette,
772     &                        MPI_COMM_AGRIF,statut,code)
773
774                imin_recv(:,k) = iminmax_temp(:,1,k)
775                imax_recv(:,k) = iminmax_temp(:,2,k)
776            endif
777
778C
779      enddo
780
781          End Subroutine ExchangeSamelevel_first         
782
783       Subroutine ExchangeSameLevel2(sendtoproc,recvfromproc,
784     &   nbdim,
785     &          pttruetabwhole,cetruetabwhole,imin,imax,
786     &          imin_recv,imax_recv,memberout,tempC,tempCextend)
787      Implicit None
788      INTEGER :: nbdim
789      INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1) :: imin,imax
790      INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1) :: pttruetabwhole,
791     &     cetruetabwhole
792      INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1) :: imin_recv,imax_recv
793      TYPE(Agrif_PVARIABLE) :: tempC,tempCextend
794      LOGICAL,DIMENSION(0:Agrif_Nbprocs-1)       :: sendtoproc
795      LOGICAL,DIMENSION(0:Agrif_Nbprocs-1)       :: recvfromproc
796      LOGICAL                                    :: res
797      LOGICAL :: memberout
798      TYPE(AGRIF_PVARIABLE), SAVE                      :: temprecv
799
800      INCLUDE 'mpif.h'
801          INTEGER :: i,k
802          INTEGER :: etiquette = 100
803          INTEGER :: code, datasize
804          INTEGER,DIMENSION(MPI_STATUS_SIZE)   :: statut
805
806      IF (memberout) THEN
807      Call Agrif_nbdim_allocation(tempCextend%var,
808     &                 pttruetabwhole(:,Agrif_ProcRank),
809     &                 cetruetabwhole(:,Agrif_ProcRank),nbdim)
810      call Agrif_nbdim_Full_VarEQreal(tempCextend%var,0.,nbdim)
811      ENDIF
812
813      IF (sendtoproc(Agrif_ProcRank)) THEN
814           Call Agrif_nbdim_VarEQvar(tempCextend%var,
815     &                               imin(:,Agrif_Procrank),
816     &                               imax(:,Agrif_Procrank),
817     &                               tempC%var,
818     &                               imin(:,Agrif_Procrank),
819     &                               imax(:,Agrif_Procrank),
820     &                               nbdim)
821      ENDIF
822
823      do k = 0,Agrif_ProcRank-1
824C
825C
826C
827            if (sendtoproc(k)) then
828C
829                datasize = 1
830C
831!CDIR SHORTLOOP
832                do i = 1,nbdim
833C
834                  datasize = datasize * (imax(i,k)-imin(i,k)+1)
835C
836                enddo
837C
838
839                SELECT CASE(nbdim)
840                CASE(1)
841                   Call MPI_SEND(tempC%var%array1(
842     &                        imin(1,k):imax(1,k)),
843     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
844     &                        MPI_COMM_AGRIF,code)
845                CASE(2)               
846                   Call MPI_SEND(tempC%var%array2(
847     &                        imin(1,k):imax(1,k),
848     &                        imin(2,k):imax(2,k)),
849     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
850     &                        MPI_COMM_AGRIF,code)
851                CASE(3)
852     
853                  Call Agrif_Send_3Darray(tempC%var%array3,
854     &             lbound(tempC%var%array3),imin(:,k),imax(:,k),k)
855                CASE(4)
856                   Call MPI_SEND(tempC%var%array4(
857     &                        imin(1,k):imax(1,k),
858     &                        imin(2,k):imax(2,k),
859     &                        imin(3,k):imax(3,k),
860     &                        imin(4,k):imax(4,k)),
861     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
862     &                        MPI_COMM_AGRIF,code)
863                CASE(5)
864                   Call MPI_SEND(tempC%var%array5(
865     &                        imin(1,k):imax(1,k),
866     &                        imin(2,k):imax(2,k),
867     &                        imin(3,k):imax(3,k),
868     &                        imin(4,k):imax(4,k),
869     &                        imin(5,k):imax(5,k)),
870     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
871     &                        MPI_COMM_AGRIF,code)
872                CASE(6)
873                   Call MPI_SEND(tempC%var%array6(
874     &                        imin(1,k):imax(1,k),
875     &                        imin(2,k):imax(2,k),
876     &                        imin(3,k):imax(3,k),
877     &                        imin(4,k):imax(4,k),
878     &                        imin(5,k):imax(5,k),
879     &                        imin(6,k):imax(6,k)),
880     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
881     &                        MPI_COMM_AGRIF,code)
882                END SELECT
883C
884            endif
885
886C
887      enddo
888C
889C
890C     Reception from others processors of the necessary part of the parent grid
891      do k = Agrif_ProcRank+1,Agrif_Nbprocs-1
892
893C
894            if (recvfromproc(k)) then
895
896                datasize = 1
897C
898!CDIR SHORTLOOP
899                do i = 1,nbdim
900C
901                datasize = datasize * (imax_recv(i,k)-imin_recv(i,k)+1)
902C
903                enddo
904
905             IF (.Not.Associated(temprecv%var)) allocate(temprecv%var)
906             call Agrif_nbdim_allocation(temprecv%var,imin_recv(:,k),
907     &   imax_recv(:,k),nbdim)
908            SELECT CASE(nbdim)
909            CASE(1)
910              Call MPI_RECV(temprecv%var%array1,
911     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
912     &               MPI_COMM_AGRIF,statut,code)
913            CASE(2)           
914              Call MPI_RECV(temprecv%var%array2,
915     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
916     &               MPI_COMM_AGRIF,statut,code)
917            CASE(3)       
918              Call MPI_RECV(temprecv%var%array3,
919     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
920     &               MPI_COMM_AGRIF,statut,code)
921
922            CASE(4)
923              Call MPI_RECV(temprecv%var%array4,
924     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
925     &               MPI_COMM_AGRIF,statut,code)
926            CASE(5)
927              Call MPI_RECV(temprecv%var%array5,
928     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
929     &               MPI_COMM_AGRIF,statut,code)
930            CASE(6)
931              Call MPI_RECV(temprecv%var%array6,
932     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
933     &               MPI_COMM_AGRIF,statut,code)
934       END SELECT
935                       
936            Call where_valtabtotab_mpi(tempCextend%var,
937     &             temprecv%var,imin_recv(:,k),imax_recv(:,k),0.,nbdim)
938     
939                Call Agrif_nbdim_deallocation(temprecv%var,nbdim)
940C                deallocate(temprecv%var)
941
942            endif
943
944C
945      enddo
946
947C     Reception from others processors of the necessary part of the parent grid
948      do k = Agrif_ProcRank+1,Agrif_Nbprocs-1
949C
950C
951            if (sendtoproc(k)) then
952C
953                SELECT CASE(nbdim)
954                CASE(1)
955                datasize=SIZE(tempC%var%array1(
956     &                        imin(1,k):imax(1,k)))
957                   Call MPI_SEND(tempC%var%array1(
958     &                        imin(1,k):imax(1,k)),
959     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
960     &                        MPI_COMM_AGRIF,code)
961                CASE(2)               
962                datasize=SIZE(tempC%var%array2(
963     &                        imin(1,k):imax(1,k),
964     &                        imin(2,k):imax(2,k)))
965                   Call MPI_SEND(tempC%var%array2(
966     &                        imin(1,k):imax(1,k),
967     &                        imin(2,k):imax(2,k)),
968     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
969     &                        MPI_COMM_AGRIF,code)
970                CASE(3)
971                datasize=SIZE(tempC%var%array3(
972     &                        imin(1,k):imax(1,k),
973     &                        imin(2,k):imax(2,k),
974     &                        imin(3,k):imax(3,k)))
975                   Call MPI_SEND(tempC%var%array3(
976     &                        imin(1,k):imax(1,k),
977     &                        imin(2,k):imax(2,k),
978     &                        imin(3,k):imax(3,k)),
979     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
980     &                        MPI_COMM_AGRIF,code)
981                CASE(4)
982                datasize=SIZE(tempC%var%array4(
983     &                        imin(1,k):imax(1,k),
984     &                        imin(2,k):imax(2,k),
985     &                        imin(3,k):imax(3,k),
986     &                        imin(4,k):imax(4,k)))
987                   Call MPI_SEND(tempC%var%array4(
988     &                        imin(1,k):imax(1,k),
989     &                        imin(2,k):imax(2,k),
990     &                        imin(3,k):imax(3,k),
991     &                        imin(4,k):imax(4,k)),
992     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
993     &                        MPI_COMM_AGRIF,code)
994                CASE(5)
995                datasize=SIZE(tempC%var%array5(
996     &                        imin(1,k):imax(1,k),
997     &                        imin(2,k):imax(2,k),
998     &                        imin(3,k):imax(3,k),
999     &                        imin(4,k):imax(4,k),
1000     &                        imin(5,k):imax(5,k)))
1001                   Call MPI_SEND(tempC%var%array5(
1002     &                        imin(1,k):imax(1,k),
1003     &                        imin(2,k):imax(2,k),
1004     &                        imin(3,k):imax(3,k),
1005     &                        imin(4,k):imax(4,k),
1006     &                        imin(5,k):imax(5,k)),
1007     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
1008     &                        MPI_COMM_AGRIF,code)
1009                CASE(6)
1010                datasize=SIZE(tempC%var%array6(
1011     &                        imin(1,k):imax(1,k),
1012     &                        imin(2,k):imax(2,k),
1013     &                        imin(3,k):imax(3,k),
1014     &                        imin(4,k):imax(4,k),
1015     &                        imin(5,k):imax(5,k),
1016     &                        imin(6,k):imax(6,k)))
1017                   Call MPI_SEND(tempC%var%array6(
1018     &                        imin(1,k):imax(1,k),
1019     &                        imin(2,k):imax(2,k),
1020     &                        imin(3,k):imax(3,k),
1021     &                        imin(4,k):imax(4,k),
1022     &                        imin(5,k):imax(5,k),
1023     &                        imin(6,k):imax(6,k)),
1024     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
1025     &                        MPI_COMM_AGRIF,code)
1026                END SELECT
1027C
1028            endif
1029
1030C
1031      enddo
1032C
1033C
1034C     Reception from others processors of the necessary part of the parent grid
1035      do k = Agrif_ProcRank-1,0,-1
1036C
1037
1038C
1039            if (recvfromproc(k)) then
1040C
1041             IF (.Not.Associated(temprecv%var)) allocate(temprecv%var)
1042             call Agrif_nbdim_allocation(temprecv%var,
1043     &   imin_recv(:,k),imax_recv(:,k),nbdim)
1044            SELECT CASE(nbdim)
1045            CASE(1)
1046              datasize=SIZE(temprecv%var%array1)
1047              Call MPI_RECV(temprecv%var%array1,
1048     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
1049     &               MPI_COMM_AGRIF,statut,code)
1050            CASE(2)           
1051              datasize=SIZE(temprecv%var%array2)
1052              Call MPI_RECV(temprecv%var%array2,
1053     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
1054     &               MPI_COMM_AGRIF,statut,code)
1055            CASE(3)           
1056              datasize=SIZE(temprecv%var%array3)
1057              Call MPI_RECV(temprecv%var%array3,
1058     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
1059     &               MPI_COMM_AGRIF,statut,code)
1060
1061            CASE(4)
1062              datasize=SIZE(temprecv%var%array4)
1063              Call MPI_RECV(temprecv%var%array4,
1064     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
1065     &               MPI_COMM_AGRIF,statut,code)
1066            CASE(5)
1067              datasize=SIZE(temprecv%var%array5)
1068              Call MPI_RECV(temprecv%var%array5,
1069     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
1070     &               MPI_COMM_AGRIF,statut,code)
1071            CASE(6)
1072              datasize=SIZE(temprecv%var%array6)
1073              Call MPI_RECV(temprecv%var%array6,
1074     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette,
1075     &               MPI_COMM_AGRIF,statut,code)
1076          END SELECT
1077           
1078            Call where_valtabtotab_mpi(tempCextend%var,
1079     &             temprecv%var,imin_recv(:,k),imax_recv(:,k)
1080     &            ,0.,nbdim)
1081     
1082                Call Agrif_nbdim_deallocation(temprecv%var,nbdim)
1083C                deallocate(temprecv%var)
1084            endif
1085
1086C
1087      enddo
1088                 
1089          End Subroutine ExchangeSamelevel2
1090         
1091          Subroutine Agrif_Send_3Darray(tab3D,bounds,imin,imax,k)
1092          integer, dimension(3) :: bounds, imin, imax
1093          real,dimension(bounds(1):,bounds(2):,bounds(3):),target
1094     &                             :: tab3D
1095          integer :: k
1096          integer :: etiquette = 100
1097          integer :: datasize, code
1098        INCLUDE 'mpif.h'
1099
1100          datasize = SIZE(tab3D(
1101     &                       imin(1):imax(1),
1102     &                        imin(2):imax(2),
1103     &                        imin(3):imax(3)))
1104       
1105                   Call MPI_SEND(tab3D(
1106     &                        imin(1):imax(1),
1107     &                        imin(2):imax(2),
1108     &                        imin(3):imax(3)),
1109     &                        datasize,MPI_DOUBLE_PRECISION,k,etiquette,
1110     &                        MPI_COMM_AGRIF,code)
1111     
1112         End Subroutine Agrif_Send_3Darray
1113
1114#else
1115      Subroutine Agrif_mpp_empty()
1116      End Subroutine Agrif_mpp_empty
1117#endif
1118
1119      End Module Agrif_mpp
Note: See TracBrowser for help on using the repository browser.