My Project
 All Classes Files Functions Variables Macros
mod_hallo.F90
Go to the documentation of this file.
1 module mod_hallo
2 USE parallel
3 implicit 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
54  ELSE
56  ENDIF
57 
58  end subroutine init_mod_hallo
59 
61  IMPLICIT NONE
62 
63  ALLOCATE(buffer(maxbuffersize))
64 
65  END SUBROUTINE create_standard_mpi_buffer
66 
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 
813 end module mod_hallo
814