source: LMDZ4/trunk/libf/dyn3dpar/mod_hallo.F90 @ 989

Last change on this file since 989 was 985, checked in by Laurent Fairhead, 17 years ago

Mise a jour de dyn3dpar par rapport a dyn3d, inclusion OpenMP et filtre FFT YM
LF

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