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.
Changeset 898 for trunk/AGRIF/AGRIF_FILES/modmpp.F – NEMO

Ignore:
Timestamp:
2008-04-22T17:35:20+02:00 (16 years ago)
Author:
rblod
Message:

Correct some bugs in agrif optimization and add MPP optimization, see #42

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/AGRIF/AGRIF_FILES/modmpp.F

    r779 r898  
    3131      Contains 
    3232#ifdef AGRIF_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_WORLD,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 
    33173C 
    34174      Subroutine Get_External_Data(tempC,tempCextend,pttruetab, 
     
    249389     &                        MPI_COMM_WORLD,code) 
    250390                CASE(3) 
     391                 
    251392                  Call Agrif_Send_3Darray(tempC%var%array3, 
    252393     &             lbound(tempC%var%array3),imin(:,k),imax(:,k),k) 
     
    529670 
    530671          End Subroutine ExchangeSamelevel 
    531  
     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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,statut,code) 
     913            CASE(2)            
     914              Call MPI_RECV(temprecv%var%array2, 
     915     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette, 
     916     &               MPI_COMM_WORLD,statut,code) 
     917            CASE(3)        
     918              Call MPI_RECV(temprecv%var%array3, 
     919     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette, 
     920     &               MPI_COMM_WORLD,statut,code) 
     921 
     922            CASE(4) 
     923              Call MPI_RECV(temprecv%var%array4, 
     924     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette, 
     925     &               MPI_COMM_WORLD,statut,code) 
     926            CASE(5) 
     927              Call MPI_RECV(temprecv%var%array5, 
     928     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette, 
     929     &               MPI_COMM_WORLD,statut,code) 
     930            CASE(6) 
     931              Call MPI_RECV(temprecv%var%array6, 
     932     &               datasize,MPI_DOUBLE_PRECISION,k,etiquette, 
     933     &               MPI_COMM_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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_WORLD,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           
    5321091          Subroutine Agrif_Send_3Darray(tab3D,bounds,imin,imax,k) 
    5331092          integer, dimension(3) :: bounds, imin, imax 
Note: See TracChangeset for help on using the changeset viewer.