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