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

Last change on this file since 814 was 807, checked in by lsce, 17 years ago

ACo + YM : amelioration de la gestion du buffer mpi

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