source: LMDZ5/trunk/libf/dyn3dpar/mod_hallo.F90 @ 2194

Last change on this file since 2194 was 1907, checked in by lguez, 11 years ago

Added a copyright property to every file of the distribution, except
for the fcm files (which have their own copyright). Use svn propget on
a file to see the copyright. For instance:

$ svn propget copyright libf/phylmd/physiq.F90
Name of program: LMDZ
Creation date: 1984
Version: LMDZ5
License: CeCILL version 2
Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
See the license file in the root directory

Also added the files defining the CeCILL version 2 license, in French
and English, at the top of the LMDZ tree.

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 21.3 KB
Line 
1module mod_Hallo
2USE parallel_lmdz
3implicit none
4  logical,save :: use_mpi_alloc
5  integer, parameter :: MaxRequest=200
6  integer, parameter :: MaxProc=80
7  integer, parameter :: MaxBufferSize=1024*1024*16
8  integer, parameter :: ListSize=1000
9 
10  integer,save       :: MaxBufferSize_Used
11!$OMP THREADPRIVATE( MaxBufferSize_Used)
12
13   real,save,pointer,dimension(:) :: Buffer
14!$OMP THREADPRIVATE(Buffer)
15
16   integer,save,dimension(Listsize) :: Buffer_Pos
17   integer,save :: Index_Pos
18!$OMP THREADPRIVATE(Buffer_Pos,Index_pos)
19   
20  type Hallo
21    real, dimension(:,:),pointer :: Field
22    integer :: offset
23    integer :: size
24    integer :: NbLevel
25    integer :: Stride
26  end type Hallo
27 
28  type request_SR
29    integer :: NbRequest=0
30    integer :: Pos
31    integer :: Index
32    type(Hallo),dimension(MaxRequest) :: Hallo
33    integer :: MSG_Request
34  end type request_SR
35
36  type request
37    type(request_SR),dimension(0:MaxProc-1) :: RequestSend
38    type(request_SR),dimension(0:MaxProc-1) :: RequestRecv
39    integer :: tag=1
40  end type request
41 
42   
43  contains
44
45  subroutine Init_mod_hallo
46    implicit none
47
48    Index_Pos=1
49    Buffer_Pos(Index_Pos)=1
50    MaxBufferSize_Used=0
51
52    IF (use_mpi_alloc .AND. using_mpi) THEN
53      CALL create_global_mpi_buffer
54    ELSE
55      CALL create_standard_mpi_buffer
56    ENDIF
57     
58  end subroutine init_mod_hallo
59
60  SUBROUTINE create_standard_mpi_buffer
61  IMPLICIT NONE
62   
63    ALLOCATE(Buffer(MaxBufferSize))
64   
65  END SUBROUTINE create_standard_mpi_buffer
66 
67  SUBROUTINE create_global_mpi_buffer
68  IMPLICIT NONE
69#ifdef CPP_MPI
70  INCLUDE 'mpif.h'
71#endif 
72    POINTER (Pbuffer,MPI_Buffer(MaxBufferSize))
73    REAL :: MPI_Buffer
74#ifdef CPP_MPI
75    INTEGER(KIND=MPI_ADDRESS_KIND) :: BS
76#else
77    INTEGER(KIND=8) :: BS
78#endif
79    INTEGER :: i,ierr
80
81!  Allocation du buffer MPI
82      Bs=8*MaxBufferSize
83!$OMP CRITICAL (MPI)
84#ifdef CPP_MPI
85      CALL MPI_ALLOC_MEM(BS,MPI_INFO_NULL,Pbuffer,ierr)
86#endif
87!$OMP END CRITICAL (MPI)
88      DO i=1,MaxBufferSize
89        MPI_Buffer(i)=i
90      ENDDO
91     
92      CALL  Associate_buffer(MPI_Buffer)
93     
94  CONTAINS
95     
96     SUBROUTINE Associate_buffer(MPI_Buffer)
97     IMPLICIT NONE
98       REAL,DIMENSION(:),target :: MPI_Buffer 
99
100         Buffer=>MPI_Buffer
101 
102      END SUBROUTINE  Associate_buffer
103                                     
104  END SUBROUTINE create_global_mpi_buffer
105 
106     
107  subroutine allocate_buffer(Size,Index,Pos)
108  implicit none
109    integer :: Size
110    integer :: Index
111    integer :: Pos
112
113    if (Buffer_pos(Index_pos)+Size>MaxBufferSize_Used) MaxBufferSize_Used=Buffer_pos(Index_pos)+Size 
114    if (Buffer_pos(Index_pos)+Size>MaxBufferSize) then
115      print *,'STOP :: La taille de MaxBufferSize dans mod_hallo.F90 est trop petite !!!!'
116      stop
117    endif
118   
119    if (Index_pos>=ListSize) then
120      print *,'STOP :: La taille de ListSize dans mod_hallo.F90 est trop petite !!!!'
121      stop
122    endif
123     
124    Pos=Buffer_Pos(Index_Pos)
125    Buffer_Pos(Index_pos+1)=Buffer_Pos(Index_Pos)+Size
126    Index_Pos=Index_Pos+1
127    Index=Index_Pos
128   
129  end subroutine allocate_buffer
130     
131  subroutine deallocate_buffer(Index)
132  implicit none
133    integer :: Index
134   
135    Buffer_Pos(Index)=-1
136   
137    do while (Buffer_Pos(Index_Pos)==-1 .and. Index_Pos>1)
138      Index_Pos=Index_Pos-1
139    end do
140
141  end subroutine deallocate_buffer 
142 
143  subroutine SetTag(a_request,tag)
144  implicit none
145    type(request):: a_request
146    integer :: tag
147   
148    a_request%tag=tag
149  end subroutine SetTag
150 
151 
152  subroutine Init_Hallo(Field,Stride,NbLevel,offset,size,NewHallo)
153    integer :: Stride
154    integer :: NbLevel
155    integer :: size
156    integer :: offset
157    real, dimension(Stride,NbLevel),target :: Field
158    type(Hallo) :: NewHallo
159   
160    NewHallo%Field=>Field
161    NewHallo%Stride=Stride
162    NewHallo%NbLevel=NbLevel
163    NewHallo%size=size
164    NewHallo%offset=offset
165   
166   
167  end subroutine Init_Hallo
168 
169  subroutine Register_SendField(Field,ij,ll,offset,size,target,a_request)
170  implicit none
171
172#include "dimensions.h"
173#include "paramet.h"   
174   
175      INTEGER :: ij,ll,offset,size,target
176      REAL, dimension(ij,ll) :: Field
177      type(request),target :: a_request
178      type(request_SR),pointer :: Ptr_request
179
180      Ptr_Request=>a_request%RequestSend(target)
181      Ptr_Request%NbRequest=Ptr_Request%NbRequest+1
182      if (Ptr_Request%NbRequest>=MaxRequest) then
183        print *,'STOP :: La taille de MaxRequest dans mod_hallo.F90 est trop petite !!!!'
184        stop
185      endif     
186      call Init_Hallo(Field,ij,ll,offset,size,Ptr_request%Hallo(Ptr_Request%NbRequest))
187     
188   end subroutine Register_SendField     
189     
190  subroutine Register_RecvField(Field,ij,ll,offset,size,target,a_request)
191  implicit none
192
193#include "dimensions.h"
194#include "paramet.h"   
195   
196      INTEGER :: ij,ll,offset,size,target
197      REAL, dimension(ij,ll) :: Field
198      type(request),target :: a_request
199      type(request_SR),pointer :: Ptr_request
200
201      Ptr_Request=>a_request%RequestRecv(target)
202      Ptr_Request%NbRequest=Ptr_Request%NbRequest+1
203     
204      if (Ptr_Request%NbRequest>=MaxRequest) then
205        print *,'STOP :: La taille de MaxRequest dans mod_hallo.F90 est trop petite !!!!'
206        stop
207      endif   
208           
209      call Init_Hallo(Field,ij,ll,offset,size,Ptr_request%Hallo(Ptr_Request%NbRequest))
210
211     
212   end subroutine Register_RecvField     
213 
214  subroutine Register_SwapField(FieldS,FieldR,ij,ll,jj_Nb_New,a_request)
215 
216      implicit none
217#include "dimensions.h"
218#include "paramet.h"   
219   
220    INTEGER :: ij,ll
221    REAL, dimension(ij,ll) :: FieldS
222    REAL, dimension(ij,ll) :: FieldR
223    type(request) :: a_request
224    integer,dimension(0:MPI_Size-1) :: jj_Nb_New   
225    integer,dimension(0:MPI_Size-1) :: jj_Begin_New,jj_End_New
226   
227    integer ::i,jje,jjb
228   
229    jj_begin_New(0)=1
230    jj_End_New(0)=jj_Nb_New(0)
231    do i=1,MPI_Size-1
232      jj_begin_New(i)=jj_end_New(i-1)+1
233      jj_end_New(i)=jj_begin_new(i)+jj_Nb_New(i)-1
234    enddo
235   
236    do i=0,MPI_Size-1
237      if (i /= MPI_Rank) then
238        jjb=max(jj_begin_new(i),jj_begin)
239        jje=min(jj_end_new(i),jj_end)
240       
241        if (ij==ip1jm .and. jje==jjp1) jje=jjm
242       
243        if (jje >= jjb) then
244          call Register_SendField(FieldS,ij,ll,jjb,jje-jjb+1,i,a_request)
245        endif
246       
247        jjb=max(jj_begin_new(MPI_Rank),jj_begin_Para(i))
248        jje=min(jj_end_new(MPI_Rank),jj_end_Para(i))
249       
250        if (ij==ip1jm .and. jje==jjp1) jje=jjm
251       
252        if (jje >= jjb) then
253          call Register_RecvField(FieldR,ij,ll,jjb,jje-jjb+1,i,a_request)
254        endif
255       
256      endif
257    enddo
258   
259  end subroutine Register_SwapField   
260 
261 
262    subroutine Register_SwapFieldHallo(FieldS,FieldR,ij,ll,jj_Nb_New,Up,Down,a_request)
263 
264      implicit none
265#include "dimensions.h"
266#include "paramet.h"   
267   
268    INTEGER :: ij,ll,Up,Down
269    REAL, dimension(ij,ll) :: FieldS
270    REAL, dimension(ij,ll) :: FieldR
271    type(request) :: a_request
272    integer,dimension(0:MPI_Size-1) :: jj_Nb_New   
273    integer,dimension(0:MPI_Size-1) :: jj_Begin_New,jj_End_New
274   
275    integer ::i,jje,jjb
276   
277    jj_begin_New(0)=1
278    jj_End_New(0)=jj_Nb_New(0)
279    do i=1,MPI_Size-1
280      jj_begin_New(i)=jj_end_New(i-1)+1
281      jj_end_New(i)=jj_begin_new(i)+jj_Nb_New(i)-1
282    enddo
283   
284    do i=0,MPI_Size-1
285      jj_begin_New(i)=max(1,jj_begin_New(i)-Up)
286      jj_end_New(i)=min(jjp1,jj_end_new(i)+Down)
287    enddo
288   
289    do i=0,MPI_Size-1
290      if (i /= MPI_Rank) then
291        jjb=max(jj_begin_new(i),jj_begin)
292        jje=min(jj_end_new(i),jj_end)
293       
294        if (ij==ip1jm .and. jje==jjp1) jje=jjm
295       
296        if (jje >= jjb) then
297          call Register_SendField(FieldS,ij,ll,jjb,jje-jjb+1,i,a_request)
298        endif
299       
300        jjb=max(jj_begin_new(MPI_Rank),jj_begin_Para(i))
301        jje=min(jj_end_new(MPI_Rank),jj_end_Para(i))
302       
303        if (ij==ip1jm .and. jje==jjp1) jje=jjm
304       
305        if (jje >= jjb) then
306          call Register_RecvField(FieldR,ij,ll,jjb,jje-jjb+1,i,a_request)
307        endif
308       
309      endif
310    enddo
311   
312  end subroutine Register_SwapFieldHallo
313 
314  subroutine Register_Hallo(Field,ij,ll,RUp,Rdown,SUp,SDown,a_request)
315 
316      implicit none
317#include "dimensions.h"
318#include "paramet.h"   
319#ifdef CPP_MPI
320    include 'mpif.h'
321#endif   
322      INTEGER :: ij,ll
323      REAL, dimension(ij,ll) :: Field
324      INTEGER :: Sup,Sdown,rup,rdown
325      type(request) :: a_request
326      type(Hallo),pointer :: PtrHallo
327      LOGICAL :: SendUp,SendDown
328      LOGICAL :: RecvUp,RecvDown
329   
330 
331      SendUp=.TRUE.
332      SendDown=.TRUE.
333      RecvUp=.TRUE.
334      RecvDown=.TRUE.
335       
336      IF (pole_nord) THEN
337        SendUp=.FALSE.
338        RecvUp=.FALSE.
339      ENDIF
340 
341      IF (pole_sud) THEN
342        SendDown=.FALSE.
343        RecvDown=.FALSE.
344      ENDIF
345     
346      if (Sup.eq.0) then
347        SendUp=.FALSE.
348       endif
349     
350      if (Sdown.eq.0) then
351        SendDown=.FALSE.
352      endif
353
354      if (Rup.eq.0) then
355        RecvUp=.FALSE.
356      endif
357     
358      if (Rdown.eq.0) then
359        RecvDown=.FALSE.
360      endif
361     
362      IF (SendUp) THEN
363        call Register_SendField(Field,ij,ll,jj_begin,SUp,MPI_Rank-1,a_request)
364      ENDIF
365 
366      IF (SendDown) THEN
367        call Register_SendField(Field,ij,ll,jj_end-SDown+1,SDown,MPI_Rank+1,a_request)
368      ENDIF
369   
370 
371      IF (RecvUp) THEN
372        call Register_RecvField(Field,ij,ll,jj_begin-Rup,RUp,MPI_Rank-1,a_request)
373      ENDIF
374 
375      IF (RecvDown) THEN
376        call Register_RecvField(Field,ij,ll,jj_end+1,RDown,MPI_Rank+1,a_request)
377      ENDIF
378 
379    end subroutine Register_Hallo
380   
381    subroutine SendRequest(a_Request)
382      implicit none
383
384#include "dimensions.h"
385#include "paramet.h"
386#ifdef CPP_MPI
387      include 'mpif.h'
388#endif
389
390      type(request),target :: a_request
391      type(request_SR),pointer :: Req
392      type(Hallo),pointer :: PtrHallo
393      integer :: SizeBuffer
394      integer :: i,rank,l,ij,Pos,ierr
395      integer :: offset
396      real,dimension(:,:),pointer :: Field
397      integer :: Nb
398       
399      do rank=0,MPI_SIZE-1
400     
401        Req=>a_Request%RequestSend(rank)
402       
403        SizeBuffer=0
404        do i=1,Req%NbRequest
405          PtrHallo=>Req%Hallo(i)
406!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
407          DO l=1,PtrHallo%NbLevel
408            SizeBuffer=SizeBuffer+PtrHallo%size*iip1
409          ENDDO
410!$OMP ENDDO NOWAIT         
411        enddo
412     
413        if (SizeBuffer>0) then
414       
415          call allocate_buffer(SizeBuffer,Req%Index,Req%pos)
416
417          Pos=Req%Pos
418          do i=1,Req%NbRequest
419            PtrHallo=>Req%Hallo(i)
420            offset=(PtrHallo%offset-1)*iip1+1
421            Nb=iip1*PtrHallo%size-1
422            Field=>PtrHallo%Field
423
424!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)           
425            do l=1,PtrHallo%NbLevel
426!cdir NODEP
427              do ij=0,Nb
428                Buffer(Pos+ij)=Field(Offset+ij,l)
429              enddo
430             
431              Pos=Pos+Nb+1
432            enddo
433!$OMP END DO NOWAIT           
434          enddo
435   
436!$OMP CRITICAL (MPI)
437         
438#ifdef CPP_MPI
439         call MPI_ISSEND(Buffer(req%Pos),SizeBuffer,MPI_REAL_LMDZ,rank,a_request%tag+1000*omp_rank,     &
440                         COMM_LMDZ,Req%MSG_Request,ierr)
441#endif
442         IF (.NOT.using_mpi) THEN
443           PRINT *,'Erreur, echange MPI en mode sequentiel !!!'
444           STOP
445         ENDIF
446!         PRINT *,"-------------------------------------------------------------------"
447!         PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->"
448!         PRINT *,"Requete envoye au proc :",rank,"tag :",a_request%tag+1000*omp_rank
449!         PRINT *,"Taille du message :",SizeBuffer,"requete no :",Req%MSG_Request
450!         PRINT *,"-------------------------------------------------------------------"
451!$OMP END CRITICAL (MPI)
452        endif
453
454    enddo
455   
456           
457      do rank=0,MPI_SIZE-1
458         
459          Req=>a_Request%RequestRecv(rank)
460          SizeBuffer=0
461         
462          do i=1,Req%NbRequest
463            PtrHallo=>Req%Hallo(i)
464
465!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
466            DO l=1,PtrHallo%NbLevel
467              SizeBuffer=SizeBuffer+PtrHallo%size*iip1
468            ENDDO
469!$OMP ENDDO NOWAIT         
470          enddo
471       
472          if (SizeBuffer>0) then
473
474             call allocate_buffer(SizeBuffer,Req%Index,Req%Pos)
475!$OMP CRITICAL (MPI)
476
477#ifdef CPP_MPI
478             call MPI_IRECV(Buffer(Req%Pos),SizeBuffer,MPI_REAL_LMDZ,rank,a_request%tag+1000*omp_rank,     &
479                           COMM_LMDZ,Req%MSG_Request,ierr)
480#endif             
481             IF (.NOT.using_mpi) THEN
482               PRINT *,'Erreur, echange MPI en mode sequentiel !!!'
483               STOP
484             ENDIF
485
486!         PRINT *,"-------------------------------------------------------------------"
487!         PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->"
488!         PRINT *,"Requete en attente du proc :",rank,"tag :",a_request%tag+1000*omp_rank
489!         PRINT *,"Taille du message :",SizeBuffer,"requete no :",Req%MSG_Request
490!         PRINT *,"-------------------------------------------------------------------"
491
492!$OMP END CRITICAL (MPI)
493          endif
494     
495      enddo
496                       
497   end subroutine SendRequest
498   
499   subroutine WaitRequest(a_Request)
500   implicit none
501   
502#include "dimensions.h"
503#include "paramet.h"
504#ifdef CPP_MPI
505      include 'mpif.h'   
506#endif
507     
508      type(request),target :: a_request
509      type(request_SR),pointer :: Req
510      type(Hallo),pointer :: PtrHallo
511      integer, dimension(2*mpi_size) :: TabRequest
512#ifdef CPP_MPI
513      integer, dimension(MPI_STATUS_SIZE,2*mpi_size) :: TabStatus
514#else
515      integer, dimension(1,2*mpi_size) :: TabStatus
516#endif
517      integer :: NbRequest
518      integer :: i,rank,pos,ij,l,ierr
519      integer :: offset
520      integer :: Nb
521     
522     
523      NbRequest=0
524      do rank=0,MPI_SIZE-1
525        Req=>a_request%RequestSend(rank)
526        if (Req%NbRequest>0) then
527          NbRequest=NbRequest+1
528          TabRequest(NbRequest)=Req%MSG_Request
529        endif
530      enddo
531     
532      do rank=0,MPI_SIZE-1
533        Req=>a_request%RequestRecv(rank)
534        if (Req%NbRequest>0) then
535          NbRequest=NbRequest+1
536          TabRequest(NbRequest)=Req%MSG_Request
537        endif
538      enddo
539     
540      if (NbRequest>0) then
541!$OMP CRITICAL (MPI)
542!        PRINT *,"-------------------------------------------------------------------"
543!        PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->",NbRequest,"en attente"
544!        PRINT *,"No des requetes :",TabRequest(1:NbRequest)
545#ifdef CPP_MPI
546        call MPI_WAITALL(NbRequest,TabRequest,TabStatus,ierr)
547#endif
548!        PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->",NbRequest,"complete"
549!        PRINT *,"-------------------------------------------------------------------"
550!$OMP END CRITICAL (MPI)
551      endif
552      do rank=0,MPI_Size-1
553        Req=>a_request%RequestRecv(rank)
554        if (Req%NbRequest>0) then
555          Pos=Req%Pos
556          do i=1,Req%NbRequest
557            PtrHallo=>Req%Hallo(i)
558            offset=(PtrHallo%offset-1)*iip1+1
559            Nb=iip1*PtrHallo%size-1
560
561!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)           
562            do l=1,PtrHallo%NbLevel
563!cdir NODEP
564              do ij=0,Nb
565                PtrHallo%Field(offset+ij,l)=Buffer(Pos+ij)
566              enddo
567
568              Pos=Pos+Nb+1
569            enddo
570!$OMP ENDDO NOWAIT         
571          enddo
572        endif
573      enddo
574     
575      do rank=0,MPI_SIZE-1
576        Req=>a_request%RequestSend(rank)
577        if (Req%NbRequest>0) then
578          call deallocate_buffer(Req%Index)
579          Req%NbRequest=0
580        endif
581      enddo
582             
583      do rank=0,MPI_SIZE-1
584        Req=>a_request%RequestRecv(rank)
585        if (Req%NbRequest>0) then
586          call deallocate_buffer(Req%Index)
587          Req%NbRequest=0
588        endif
589      enddo
590     
591      a_request%tag=1
592    end subroutine WaitRequest
593     
594   subroutine WaitSendRequest(a_Request)
595   implicit none
596   
597#include "dimensions.h"
598#include "paramet.h"
599#ifdef CPP_MPI
600      include 'mpif.h'   
601#endif     
602      type(request),target :: a_request
603      type(request_SR),pointer :: Req
604      type(Hallo),pointer :: PtrHallo
605      integer, dimension(mpi_size) :: TabRequest
606#ifdef CPP_MPI
607      integer, dimension(MPI_STATUS_SIZE,mpi_size) :: TabStatus
608#else
609      integer, dimension(1,mpi_size) :: TabStatus
610#endif
611      integer :: NbRequest
612      integer :: i,rank,pos,ij,l,ierr
613      integer :: offset
614     
615     
616      NbRequest=0
617      do rank=0,MPI_SIZE-1
618        Req=>a_request%RequestSend(rank)
619        if (Req%NbRequest>0) then
620          NbRequest=NbRequest+1
621          TabRequest(NbRequest)=Req%MSG_Request
622        endif
623      enddo
624     
625
626      if (NbRequest>0) THEN
627!$OMP CRITICAL (MPI)     
628!        PRINT *,"-------------------------------------------------------------------"
629!        PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->",NbRequest,"en attente"
630!        PRINT *,"No des requetes :",TabRequest(1:NbRequest)
631#ifdef CPP_MPI
632        call MPI_WAITALL(NbRequest,TabRequest,TabStatus,ierr)
633#endif
634!        PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->",NbRequest,"complete"
635!        PRINT *,"-------------------------------------------------------------------"
636
637!$OMP END CRITICAL (MPI)
638      endif     
639     
640      do rank=0,MPI_SIZE-1
641        Req=>a_request%RequestSend(rank)
642        if (Req%NbRequest>0) then
643          call deallocate_buffer(Req%Index)
644          Req%NbRequest=0
645        endif
646      enddo
647             
648      a_request%tag=1
649    end subroutine WaitSendRequest
650   
651   subroutine WaitRecvRequest(a_Request)
652   implicit none
653   
654#include "dimensions.h"
655#include "paramet.h"
656#ifdef CPP_MPI
657      include 'mpif.h'   
658#endif
659     
660      type(request),target :: a_request
661      type(request_SR),pointer :: Req
662      type(Hallo),pointer :: PtrHallo
663      integer, dimension(mpi_size) :: TabRequest
664#ifdef CPP_MPI
665      integer, dimension(MPI_STATUS_SIZE,mpi_size) :: TabStatus
666#else
667      integer, dimension(1,mpi_size) :: TabStatus
668#endif
669      integer :: NbRequest
670      integer :: i,rank,pos,ij,l,ierr
671      integer :: offset,Nb
672     
673     
674      NbRequest=0
675     
676      do rank=0,MPI_SIZE-1
677        Req=>a_request%RequestRecv(rank)
678        if (Req%NbRequest>0) then
679          NbRequest=NbRequest+1
680          TabRequest(NbRequest)=Req%MSG_Request
681        endif
682      enddo
683     
684     
685      if (NbRequest>0) then
686!$OMP CRITICAL (MPI)     
687!        PRINT *,"-------------------------------------------------------------------"
688!        PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->",NbRequest,"en attente"
689!        PRINT *,"No des requetes :",TabRequest(1:NbRequest)
690#ifdef CPP_MPI
691        call MPI_WAITALL(NbRequest,TabRequest,TabStatus,ierr)
692#endif
693!        PRINT *,"Process de rang",mpi_rank,"Task : ",omp_rank,"--->",NbRequest,"complete"
694!        PRINT *,"-------------------------------------------------------------------"
695!$OMP END CRITICAL (MPI)     
696      endif
697     
698      do rank=0,MPI_Size-1
699        Req=>a_request%RequestRecv(rank)
700        if (Req%NbRequest>0) then
701          Pos=Req%Pos
702          do i=1,Req%NbRequest
703            PtrHallo=>Req%Hallo(i)
704            offset=(PtrHallo%offset-1)*iip1+1
705            Nb=iip1*PtrHallo%size-1
706!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)           
707            do l=1,PtrHallo%NbLevel
708!cdir NODEP
709              do ij=0,Nb
710                PtrHallo%Field(offset+ij,l)=Buffer(Pos+ij)
711              enddo
712                 Pos=Pos+Nb+1
713            enddo
714!$OMP END DO NOWAIT
715          enddo
716        endif
717      enddo
718     
719           
720      do rank=0,MPI_SIZE-1
721        Req=>a_request%RequestRecv(rank)
722        if (Req%NbRequest>0) then
723          call deallocate_buffer(Req%Index)
724          Req%NbRequest=0
725        endif
726      enddo
727     
728      a_request%tag=1
729    end subroutine WaitRecvRequest
730   
731   
732   
733    subroutine CopyField(FieldS,FieldR,ij,ll,jj_Nb_New)
734 
735      implicit none
736#include "dimensions.h"
737#include "paramet.h"   
738   
739    INTEGER :: ij,ll,l
740    REAL, dimension(ij,ll) :: FieldS
741    REAL, dimension(ij,ll) :: FieldR
742    integer,dimension(0:MPI_Size-1) :: jj_Nb_New   
743    integer,dimension(0:MPI_Size-1) :: jj_Begin_New,jj_End_New
744   
745    integer ::i,jje,jjb,ijb,ije
746   
747    jj_begin_New(0)=1
748    jj_End_New(0)=jj_Nb_New(0)
749    do i=1,MPI_Size-1
750      jj_begin_New(i)=jj_end_New(i-1)+1
751      jj_end_New(i)=jj_begin_new(i)+jj_Nb_New(i)-1
752    enddo
753   
754    jjb=max(jj_begin,jj_begin_new(MPI_Rank))
755    jje=min(jj_end,jj_end_new(MPI_Rank))
756    if (ij==ip1jm) jje=min(jje,jjm)
757
758    if (jje >= jjb) then
759      ijb=(jjb-1)*iip1+1
760      ije=jje*iip1
761
762!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
763      do l=1,ll
764        FieldR(ijb:ije,l)=FieldS(ijb:ije,l)
765      enddo
766!$OMP ENDDO NOWAIT
767    endif
768
769
770  end subroutine CopyField   
771
772  subroutine CopyFieldHallo(FieldS,FieldR,ij,ll,jj_Nb_New,Up,Down)
773 
774      implicit none
775#include "dimensions.h"
776#include "paramet.h"   
777   
778    INTEGER :: ij,ll,Up,Down
779    REAL, dimension(ij,ll) :: FieldS
780    REAL, dimension(ij,ll) :: FieldR
781    integer,dimension(0:MPI_Size-1) :: jj_Nb_New   
782    integer,dimension(0:MPI_Size-1) :: jj_Begin_New,jj_End_New
783
784    integer ::i,jje,jjb,ijb,ije,l
785
786     
787    jj_begin_New(0)=1
788    jj_End_New(0)=jj_Nb_New(0)
789    do i=1,MPI_Size-1
790      jj_begin_New(i)=jj_end_New(i-1)+1
791      jj_end_New(i)=jj_begin_new(i)+jj_Nb_New(i)-1
792    enddo
793
794       
795    jjb=max(jj_begin,jj_begin_new(MPI_Rank)-Up)
796    jje=min(jj_end,jj_end_new(MPI_Rank)+Down)
797    if (ij==ip1jm) jje=min(jje,jjm)
798   
799   
800    if (jje >= jjb) then
801      ijb=(jjb-1)*iip1+1
802      ije=jje*iip1
803
804!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
805      do l=1,ll
806        FieldR(ijb:ije,l)=FieldS(ijb:ije,l)
807      enddo
808!$OMP ENDDO NOWAIT
809
810    endif
811   end subroutine CopyFieldHallo       
812         
813end module mod_Hallo
814   
Note: See TracBrowser for help on using the repository browser.