File: /nfs/home/0/users/jenkins/mfix.git/model/dmp_modules/sendrecv3_mod.f

1             module sendrecv3
2             use parallel_mpi
3             use debug
4             use geometry
5             use compar
6             use indices
7             implicit none
8     
9             integer, pointer, dimension(:) :: &
10                     recvproc1, recvtag1, xrecv1, recvijk1, &
11                     sendproc1, sendtag1, xsend1, sendijk1, &
12                     recvproc2, recvtag2, xrecv2, recvijk2, &
13                     sendproc2, sendtag2, xsend2, sendijk2
14     
15     
16             integer,pointer, dimension(:) :: &
17                 send_persistent_request, recv_persistent_request,     &
18                 send_persistent_request1, send_persistent_request2,   &
19                 recv_persistent_request1, recv_persistent_request2
20     
21             integer :: nrecv1,nsend1, nrecv2,nsend2
22     
23     !---------------
24     !EFD extra layer
25     !---------------
26             integer, pointer, dimension(:) :: &
27                     recvproc3, recvtag3, xrecv3, recvijk3, &
28                     sendproc3, sendtag3, xsend3, sendijk3, &
29                     send_persistent_request3,recv_persistent_request3
30             integer :: nrecv3,nsend3
31             integer, parameter :: nlayers = 3
32     
33             logical,parameter :: localfunc=.false.
34     
35             logical,parameter :: use_persistent_message=.true.
36     
37     
38             double precision, dimension(:), pointer :: &
39                     dsendbuffer, drecvbuffer
40             integer, dimension(:), pointer :: &
41                     isendbuffer, irecvbuffer
42             character, dimension(:), pointer :: &
43                     csendbuffer, crecvbuffer
44     
45             integer :: nrecv,nsend
46             integer, pointer, dimension(:) :: &
47                     recvrequest, sendrequest, &
48                     xrecv,recvproc, recvijk, recvtag, &
49                     xsend,sendproc, sendijk, sendtag
50     
51             integer :: &
52                kstart_all_myPE, jstart_all_myPE, istart_all_myPE, &
53                kend_all_myPE, jend_all_myPE, iend_all_myPE
54     
55             integer :: communicator
56     
57     !       -----------------
58     !       generic interface
59     !       -----------------
60             interface sendrecv3_begin
61             module procedure &
62                     sendrecv3_begin_1d, &
63                     sendrecv3_begin_1i, &
64                     sendrecv3_begin_1c
65             end interface
66     
67             interface sendrecv3_end
68             module procedure &
69                     sendrecv3_end_1d, &
70                     sendrecv3_end_1i, &
71                     sendrecv3_end_1c
72             end interface
73     
74             interface send_recv3
75             module procedure &
76                     send_recv3_1d, send_recv3_2d, send_recv3_3d, &
77                     send_recv3_1i, &
78                     send_recv3_1c
79             end interface
80     
81     
82             contains
83     
84             subroutine ijk_of( ijkp, i,j,k )
85             integer, intent(in) :: ijkp
86             integer, intent(out) :: i,j,k
87     
88             integer :: k1,k2, j1,j2, i1,i2, &
89                     ijk, isize,jsize,ksize, gijk
90     
91             character(len=32), parameter :: name = "ijk_of"
92             logical :: isok_k, isok_j, isok_i, is_same, isok
93     
94             ijk = ijkp
95     
96     
97     !---------------
98     !EFD extra layer
99     !---------------
100             i1 = istart4_all(myPE)
101             i2 = iend4_all(myPE)
102             j1 = jstart4_all(myPE)
103             j2 = jend4_all(myPE)
104             k1 = kstart4_all(myPE)
105             k2 = kend4_all(myPE)
106     
107     
108     
109     
110     
111     
112             ksize = (k2-k1+1)
113             jsize = (j2-j1+1)
114             isize = (i2-i1+1)
115     
116     
117             if (mod(ijk,isize*jsize).ne.0) then
118                     k = int( ijk/(isize*jsize) ) + k1
119             else
120                     k = int( ijk/(isize*jsize) ) + k1 -1
121             endif
122             ijk = ijk - (k-k1)*(isize*jsize)
123     
124             if (mod(ijk,isize).ne.0) then
125                     j = int( ijk/isize ) + j1
126             else
127                     j = int( ijk/isize ) + j1 - 1
128             endif
129             ijk = ijk - (j-j1)*isize
130     
131             i = (ijk-1) + i1
132     !       ------------
133     !       double check
134     !       ------------
135             isok_i = (i1 <= i) .and. (i <= i2)
136             isok_j = (j1 <= j) .and. (j <= j2)
137             isok_k = (k1 <= k) .and. (k <= k2)
138             gijk = 1 + (i-i1) + (j-j1)*(i2-i1+1) + &
139                     (k-k1)*(j2-j1+1)*(i2-i1+1)
140             is_same = (gijk .eq. ijkp)
141             isok = isok_i .and. isok_j .and. isok_k .and. is_same
142             if (.not.isok) then
143                 call write_debug( name, 'i,j,k ', i,j,k )
144                 call write_debug( name, 'ijkp, gijk ', ijkp, gijk )
145             endif
146     
147     
148             return
149             end subroutine ijk_of
150     
151     
152             subroutine ijk_of_gl( ijkp, i,j,k )
153             integer, intent(in) :: ijkp
154             integer, intent(out) :: i,j,k
155     
156             integer :: k1,k2, j1,j2, i1,i2, &
157                     ijk, isize,jsize,ksize, gijk
158     
159             character(len=32), parameter :: name = "ijk_of_gl"
160             logical :: isok_k, isok_j, isok_i, is_same, isok
161     
162             ijk = ijkp
163     
164     
165     !---------------
166     !EFD extra layer
167     !---------------
168             k1 = minval( kstart4_all(:) )
169             k2 = maxval( kend4_all(:) )
170     
171             j1 = minval( jstart4_all(:) )
172             j2 = maxval( jend4_all(:) )
173     
174             i1 = minval( istart4_all(:) )
175             i2 = maxval( iend4_all(:) )
176     
177     
178     
179             ksize = (k2-k1+1)
180             jsize = (j2-j1+1)
181             isize = (i2-i1+1)
182     
183     
184             if (mod(ijk,isize*jsize).ne.0) then
185                     k = int( ijk/(isize*jsize) ) + k1
186             else
187                     k = int( ijk/(isize*jsize) ) + k1 -1
188             endif
189             ijk = ijk - (k-k1)*(isize*jsize)
190     
191             if (mod(ijk,isize).ne.0) then
192                     j = int( ijk/isize ) + j1
193             else
194                     j = int( ijk/isize ) + j1 - 1
195             endif
196             ijk = ijk - (j-j1)*isize
197     
198             i = (ijk-1) + i1
199     !       ------------
200     !       double check
201     !       ------------
202             isok_i = (i1 <= i) .and. (i <= i2)
203             isok_j = (j1 <= j) .and. (j <= j2)
204             isok_k = (k1 <= k) .and. (k <= k2)
205             gijk = 1 + (i-i1) + (j-j1)*(i2-i1+1) + &
206                     (k-k1)*(j2-j1+1)*(i2-i1+1)
207             is_same = (gijk .eq. ijkp)
208             isok = isok_i .and. isok_j .and. isok_k .and. is_same
209             if (.not.isok) then
210                 call write_debug( name, 'i,j,k ', i,j,k )
211                 call write_debug( name, 'ijkp, gijk ', ijkp, gijk )
212             endif
213     
214     
215             return
216             end subroutine ijk_of_gl
217     
218             subroutine sendrecv3_init(        &
219                     comm,                    &
220                     cyclic_i,cyclic_j,cyclic_k, idebug )
221     
222                     use functions
223     
224             implicit none
225     
226             integer, intent(in) :: comm
227             logical,intent(in) :: cyclic_i,cyclic_j,cyclic_k
228     
229             integer, intent(in), optional :: idebug
230     
231     !       -------------------------------------
232     !       set up tables and data structures for
233     !       exchanging ghost regions
234     !       -------------------------------------
235     
236     !       ---------------
237     !       local variables
238     !       ---------------
239             character(len=80), parameter :: name = 'sendrecv3_init'
240     
241             character(len=80), pointer, dimension(:) :: line
242             integer :: ip, lmax
243     
244             integer :: layer,request, source, tag, datatype
245     
246             integer :: lidebug
247             integer :: isize,jsize,ksize, ijksize
248             integer :: recvsize1, recvsize2, &
249                        sendsize1, sendsize2
250     
251             integer :: iter, i,j,k, ii, jj,kk, &
252                     ntotal, icount,ipos, &
253                     ilayer,        i1,i2,  j1,j2, k1,k2,  &
254                     ijk, ijk2, iproc, jproc, src,dest, &
255                     ierror
256     
257             logical :: isok, isvalid, ismine, is_halobc
258     
259             integer, dimension(:,:,:), pointer :: ijk2proc
260             integer, pointer, dimension(:) :: &
261                     istartx,iendx, jstartx,jendx, kstartx,kendx, &
262                     ncount, &
263                     recvproc, recvtag, xrecv, recvijk,  &
264                     sendproc, sendtag, xsend, sendijk
265     
266             logical, parameter :: jfastest = .true.
267     
268     
269             integer, parameter :: message_tag_offset = 1000
270     
271     
272     !       ----------------
273     !       inline functions
274     !       ----------------
275             integer :: message_tag
276     
277             message_tag(src,dest) = message_tag_offset + (1+src + dest*10*numPEs)
278     
279             nullify( &
280                     recvproc1, recvtag1, xrecv1, recvijk1, &
281                     sendproc1, sendtag1, xsend1, sendijk1, &
282                     recvproc2, recvtag2, xrecv2, recvijk2, &
283                     sendproc2, sendtag2, xsend2, sendijk2)
284     
285             nullify( &
286                     send_persistent_request, recv_persistent_request,   &
287                     send_persistent_request1, send_persistent_request2, &
288                     recv_persistent_request1, recv_persistent_request2 )
289     
290             nullify( dsendbuffer, drecvbuffer )
291             nullify( isendbuffer, irecvbuffer )
292             nullify( csendbuffer, crecvbuffer )
293     
294             nullify( &
295                     recvrequest, sendrequest, &
296                     xrecv,recvproc, recvijk, recvtag, &
297                     xsend,sendproc, sendijk, sendtag )
298     
299     
300     
301     
302     !       --------------------
303     !       initialize variables
304     !       --------------------
305             lidebug = 0
306             if (present(idebug)) then
307               lidebug = idebug
308             endif
309     
310             communicator = comm
311             call MPI_COMM_SIZE( comm, numPEs, ierror )
312             call MPI_Check( 'sendrecv3_init:MPI_COMM_SIZE ', ierror )
313     
314             call MPI_COMM_RANK( comm, myPE, ierror )
315             call MPI_Check( 'sendrecv3_init:MPI_COMM_RANK ', ierror )
316     
317     !       ---------------------------
318     !       check obtain bounds of domain
319     !       ---------------------------
320     
321     
322     !       -----------------------
323     !       check bounds for k-axis
324     !       -----------------------
325             call assert( kmin1 .eq. minval( kstart1_all(:) ), &
326                     '** sendrecv3_init: invalid kmin1, ' // &
327                     ' kmin1, minval(kstart1_all(:)) ', &
328                       kmin1, minval(kstart1_all(:)) )
329     
330             call assert( kmin2 .eq. minval( kstart2_all(:) ), &
331                     '** sendrecv3_init: invalid kmin2, ' // &
332                     ' kmin2, minval(kstart2_all(:)) ', &
333                       kmin2, minval(kstart2_all(:)) )
334     
335             call assert( kmin3 .eq. minval( kstart3_all(:) ), &
336                     '** sendrecv3_init: invalid kmin3, ' // &
337                     ' kmin3, minval(kstart3_all(:)) ', &
338                       kmin3, minval(kstart3_all(:)) )
339     
340     
341             call assert( kmax1 .eq. maxval( kend1_all(:) ), &
342                     '** sendrecv3_init: invalid kmax1, ' // &
343                     ' kmax1, maxval(kend1_all(:)) ', &
344                       kmax1, maxval(kend1_all(:)) )
345     
346             call assert( kmax2 .eq. maxval( kend2_all(:) ), &
347                     '** sendrecv3_init: invalid kmax2, ' // &
348                     ' kmax2, maxval(kend2_all(:)) ', &
349                       kmax2, maxval(kend2_all(:)) )
350     
351             call assert( kmax3 .eq. maxval( kend3_all(:) ), &
352                     '** sendrecv3_init: invalid kmax3, ' // &
353                     ' kmax3, maxval(kend3_all(:)) ', &
354                       kmax3, maxval(kend3_all(:)) )
355     
356     !---------------
357     !EFD extra layer
358     !---------------
359             call assert( kmin4 .eq. minval( kstart4_all(:) ), &
360                     '** sendrecv3_init: invalid kmin4, ' // &
361                     ' kmin4, minval(kstart4_all(:)) ', &
362                       kmin4, minval(kstart4_all(:)) )
363             call assert( kmax4 .eq. maxval( kend4_all(:) ), &
364                     '** sendrecv3_init: invalid kmax4, ' // &
365                     ' kmax4, maxval(kend4_all(:)) ', &
366                       kmax4, maxval(kend4_all(:)) )
367     
368     !       -----------------------
369     !       check bounds for j-axis
370     !       -----------------------
371             call assert( jmin1 .eq. minval( jstart1_all(:) ), &
372                     '** sendrecv3_init: invalid jmin1, ' // &
373                     ' jmin1, minval(jstart1_all(:)) ', &
374                       jmin1, minval(jstart1_all(:)) )
375     
376             call assert( jmin2 .eq. minval( jstart2_all(:) ), &
377                     '** sendrecv3_init: invalid jmin2, ' // &
378                     ' jmin2, minval(jstart2_all(:)) ', &
379                       jmin2, minval(jstart2_all(:)) )
380     
381             call assert( jmin3 .eq. minval( jstart3_all(:) ), &
382                     '** sendrecv3_init: invalid jmin3, ' // &
383                     ' jmin3, minval(jstart3_all(:)) ', &
384                       jmin3, minval(jstart3_all(:)) )
385     
386     
387             call assert( jmax1 .eq. maxval( jend1_all(:) ), &
388                     '** sendrecv3_init: invalid jmax1, ' // &
389                     ' jmax1, maxval(jend1_all(:)) ', &
390                       jmax1, maxval(jend1_all(:)) )
391     
392             call assert( jmax2 .eq. maxval( jend2_all(:) ), &
393                     '** sendrecv3_init: invalid jmax2, ' // &
394                     ' jmax2, maxval(jend2_all(:)) ', &
395                       jmax2, maxval(jend2_all(:)) )
396     
397             call assert( jmax3 .eq. maxval( jend3_all(:) ), &
398                     '** sendrecv3_init: invalid jmax3, ' // &
399                     ' jmax3, maxval(jend3_all(:)) ', &
400                       jmax3, maxval(jend3_all(:)) )
401     !---------------
402     !EFD extra layer
403     !---------------
404             call assert( jmin4 .eq. minval( jstart4_all(:) ), &
405                     '** sendrecv3_init: invalid jmin4, ' // &
406                     ' jmin4, minval(jstart4_all(:)) ', &
407                       jmin4, minval(jstart4_all(:)) )
408             call assert( jmax4 .eq. maxval( jend4_all(:) ), &
409                     '** sendrecv3_init: invalid jmax4, ' // &
410                     ' jmax4, maxval(jend4_all(:)) ', &
411                       jmax4, maxval(jend4_all(:)) )
412     
413     !       -----------------------
414     !       check bounds for i-axis
415     !       -----------------------
416             call assert( imin1 .eq. minval( istart1_all(:) ), &
417                     '** sendrecv3_init: invalid imin1, ' // &
418                     ' imin1, minval(istart1_all(:)) ', &
419                       imin1, minval(istart1_all(:)) )
420     
421             call assert( imin2 .eq. minval( istart2_all(:) ), &
422                     '** sendrecv3_init: invalid imin2, ' // &
423                     ' imin2, minval(istart2_all(:)) ', &
424                       imin2, minval(istart2_all(:)) )
425     
426             call assert( imin3 .eq. minval( istart3_all(:) ), &
427                     '** sendrecv3_init: invalid imin3, ' // &
428                     ' imin3, minval(istart3_all(:)) ', &
429                       imin3, minval(istart3_all(:)) )
430     
431     
432             call assert( imax1 .eq. maxval( iend1_all(:) ), &
433                     '** sendrecv3_init: invalid imax1, ' // &
434                     ' imax1, maxval(iend1_all(:)) ', &
435                       imax1, maxval(iend1_all(:)) )
436     
437             call assert( imax2 .eq. maxval( iend2_all(:) ), &
438                     '** sendrecv3_init: invalid imax2, ' // &
439                     ' imax2, maxval(iend2_all(:)) ', &
440                       imax2, maxval(iend2_all(:)) )
441     
442             call assert( imax3 .eq. maxval( iend3_all(:) ), &
443                     '** sendrecv3_init: invalid imax3, ' // &
444                     ' imax3, maxval(iend3_all(:)) ', &
445                       imax3, maxval(iend3_all(:)) )
446     
447     !---------------
448     !EFD extra layer
449     !---------------
450             call assert( imin4 .eq. minval( istart4_all(:) ), &
451                     '** sendrecv3_init: invalid imin4, ' // &
452                     ' imin4, minval(istart4_all(:)) ', &
453                       imin4, minval(istart4_all(:)) )
454             call assert( imax4 .eq. maxval( iend4_all(:) ), &
455                     '** sendrecv3_init: invalid imax4, ' // &
456                     ' imax4, maxval(iend4_all(:)) ', &
457                       imax4, maxval(iend4_all(:)) )
458     
459     
460     
461             call assert( jmin1 .le. jmax1, &
462                     '** sendrecv3_init: jmin1,jmax1 ', jmin1,jmax1 )
463             call assert( jmin2 .le. jmax2, &
464                     '** sendrecv3_init: jmin2,jmax2 ', jmin2,jmax2 )
465             call assert( jmin3 .le. jmax3, &
466                     '** sendrecv3_init: jmin3,jmax3 ', jmin3,jmax3 )
467     
468             call assert( kmin1 .le. kmax1, &
469                     '** sendrecv3_init: kmin1,kmax1 ', kmin1,kmax1 )
470             call assert( kmin2 .le. kmax2, &
471                     '** sendrecv3_init: kmin2,kmax2 ', kmin2,kmax2 )
472             call assert( kmin3 .le. kmax3, &
473                     '** sendrecv3_init: kmin3,kmax3 ', kmin3,kmax3 )
474     
475             call assert( imin1 .le. imax1, &
476                     '** sendrecv3_init: imin1,imax1 ', imin1,imax1 )
477             call assert( imin2 .le. imax2, &
478                     '** sendrecv3_init: imin2,imax2 ', imin2,imax2 )
479             call assert( imin3 .le. imax3, &
480                     '** sendrecv3_init: imin3,imax3 ', imin3,imax3 )
481     
482     
483     !---------------
484     !EFD extra layer
485     !---------------
486             call assert( jmin4 .le. jmax4, &
487                     '** sendrecv3_init: jmin4,jmax4 ', jmin4,jmax4 )
488             call assert( imin4 .le. imax4, &
489                     '** sendrecv3_init: imin4,imax4 ', imin4,imax4 )
490             call assert( kmin4 .le. kmax4, &
491                     '** sendrecv3_init: kmin4,kmax4 ', kmin4,kmax4 )
492     
493     
494     
495     
496     
497     
498             k1 = min( kmin1, min(kmin2, kmin3) )
499             k2 = max( kmax1, max(kmax2, kmax3) )
500             j1 = min( jmin1, min(jmin2, jmin3) )
501             j2 = max( jmax1, max(jmax2, jmax3) )
502             i1 = min( imin1, min(imin2, imin3) )
503             i2 = max( imax1, max(imax2, imax3) )
504     
505     !----------------
506     !EFD extra layer
507     !----------------
508             k1 = min(k1,kmin4)
509             k2 = max(k2,kmax4)
510             j1 = min(j1,jmin4)
511             j2 = max(j2,jmax4)
512             i1 = min(i1,imin4)
513             i2 = max(i2,imax4)
514     
515     
516             allocate( ijk2proc( i1:i2, j1:j2, k1:k2 ) )
517     
518     
519             if(localfunc) then
520     !       --------------------------------------
521     !       double check ijk_of()
522     !       --------------------------------------
523     
524     
525     !---------------
526     !EFD extra layer
527     !---------------
528             kstart_all_myPE = min( kstart3_all(myPE), kstart4_all(myPE))
529             kend_all_myPE = max( kend3_all(myPE), kend4_all(myPE))
530     
531             jstart_all_myPE = min( jstart3_all(myPE), jstart4_all(myPE))
532             jend_all_myPE = max( jend3_all(myPE), jend4_all(myPE))
533     
534             istart_all_myPE = min( istart3_all(myPE), istart4_all(myPE))
535             iend_all_myPE = max( iend3_all(myPE), iend4_all(myPE))
536     
537     
538             do k=kstart_all_myPE,kend_all_myPE
539             do j=jstart_all_myPE,jend_all_myPE
540             do i=istart_all_myPE,iend_all_myPE
541                ijk = funijk(i,j,k)
542                call ijk_of(ijk, ii,jj,kk)
543                ijk2 = funijk( ii,jj,kk)
544     
545                isvalid = (ii.eq.i).and.(jj.eq.j).and.(kk.eq.k).and.(ijk.eq.ijk2)
546                if (.not.isvalid) then
547                     call write_debug( name, 'error with ijk_of ')
548     
549                     call write_debug( name, 'istart3_all(myPE),iend3_all(myPE) ', &
550                                     istart3_all(myPE),iend3_all(myPE) )
551                     call write_debug( name, 'jstart3_all(myPE),jend3_all(myPE) ', &
552                                     jstart3_all(myPE),jend3_all(myPE) )
553                     call write_debug( name, 'kstart3_all(myPE),kend3_all(myPE) ', &
554                                     kstart3_all(myPE),kend3_all(myPE) )
555     
556                     call write_debug( name, 'i,j,k, ijk ', i,j,k, ijk )
557                     call write_debug( name, 'ii,jj,kk,  ijk2 ',ii,jj,kk,ijk2 )
558     
559                endif
560             enddo
561             enddo
562             enddo
563             endif ! Local Function
564     
565     
566             if (lidebug.ge.1) then
567                call write_debug( name, 'imap ', imap )
568                call write_debug( name, 'jmap ', jmap )
569                call write_debug( name, 'kmap ', kmap )
570             endif
571     
572     
573     !       ----------------------------
574     !       set up table ijk2proc(:,:,:)
575     !
576     !       ijk2proc(i,j,k) maps (i,j,k) index to
577     !       unique processor that 'owns' that node.
578     !       ----------------------------
579     
580     
581     
582             ijk2proc( :,:,: ) = 0
583     
584     !       --------------------------------------------------
585     !       double check domain decomposition that
586     !       each interior node is assigned to UNIQUE processor
587     !       --------------------------------------------------
588             do iproc=0,numPEs-1
589     
590     
591     
592     
593                i1 = istart1_all(iproc)
594                i2 = iend1_all(iproc)
595                j1 = jstart1_all(iproc)
596                j2 = jend1_all(iproc)
597                k1 = kstart1_all(iproc)
598                k2 = kend1_all(iproc)
599     
600                if(istart3_all(iproc).eq.imin3) i1 = istart3_all(iproc)
601                if(iend3_all(iproc).eq.imax3) i2 = iend3_all(iproc)
602                if(jstart3_all(iproc).eq.jmin3) j1 = jstart3_all(iproc)
603                if(jend3_all(iproc).eq.jmax3) j2 = jend3_all(iproc)
604                if(kstart3_all(iproc).eq.kmin3) k1 = kstart3_all(iproc)
605                if(kend3_all(iproc).eq.kmax3) k2 = kend3_all(iproc)
606     
607     
608     
609     
610                do k=k1,k2
611                do j=j1,j2
612                do i=i1,i2
613                     ijk2proc(i,j,k) = ijk2proc(i,j,k) + 1
614                enddo
615                enddo
616                enddo
617     
618             enddo
619     
620             do k=lbound(ijk2proc,3),ubound(ijk2proc,3)
621             do j=lbound(ijk2proc,2),ubound(ijk2proc,2)
622             do i=lbound(ijk2proc,1),ubound(ijk2proc,1)
623               isvalid = (ijk2proc(i,j,k) .eq. 1)
624               if (.not.isvalid) then
625     
626                     call write_debug(name, ' invalid decomposition ')
627                     call write_debug(name, 'i,j,k ',i,j,k )
628                     call write_debug(name, 'ijk2proc(i,j,k) ', ijk2proc(i,j,k))
629     
630                     call mfix_exit( myPE )
631               endif
632             enddo
633             enddo
634             enddo
635     
636             ijk2proc(:,:,:) = -1
637             do iproc=0,numPEs-1
638     
639     
640     
641                i1 = istart1_all(iproc)
642                i2 = iend1_all(iproc)
643                j1 = jstart1_all(iproc)
644                j2 = jend1_all(iproc)
645                k1 = kstart1_all(iproc)
646                k2 = kend1_all(iproc)
647     
648                if(istart3_all(iproc).eq.imin3) i1 = istart3_all(iproc)
649                if(iend3_all(iproc).eq.imax3) i2 = iend3_all(iproc)
650                if(jstart3_all(iproc).eq.jmin3) j1 = jstart3_all(iproc)
651                if(jend3_all(iproc).eq.jmax3) j2 = jend3_all(iproc)
652                if(kstart3_all(iproc).eq.kmin3) k1 = kstart3_all(iproc)
653                if(kend3_all(iproc).eq.kmax3) k2 = kend3_all(iproc)
654     
655     
656     !---------------
657     !EFD extra layer
658     !but not sure if this is needed
659     !---------------
660     !           if(istart4_all(iproc).eq.imin4) i1 = istart4_all(iproc)
661     !           if(iend4_all(iproc).eq.imax4) i2 = iend4_all(iproc)
662     !           if(jstart4_all(iproc).eq.jmin4) j1 = jstart4_all(iproc)
663     !           if(jend4_all(iproc).eq.jmax4) j2 = jend4_all(iproc)
664     !           if(kstart4_all(iproc).eq.kmin4) k1 = kstart4_all(iproc)
665     !           if(kend4_all(iproc).eq.kmax4) k2 = kend4_all(iproc)
666     
667     
668                do k=k1,k2
669                do j=j1,j2
670                do i=i1,i2
671                     ijk2proc(i,j,k) = iproc
672                enddo
673                enddo
674                enddo
675     
676             enddo
677     
678     
679     
680        allocate( ncount(0:numPEs-1) )
681     
682        allocate( istartx(0:numPEs-1) )
683        allocate( jstartx(0:numPEs-1) )
684        allocate( kstartx(0:numPEs-1) )
685     
686        allocate( iendx(0:numPEs-1) )
687        allocate( jendx(0:numPEs-1) )
688        allocate( kendx(0:numPEs-1) )
689     
690     !---------------
691     !EFD extra layer
692     !---------------
693     do ilayer=nlayers,nlayers
694     
695        if (ilayer.eq.1) then
696             kstartx(:) = kstart2_all(:)
697             kendx(:) = kend2_all(:)
698             jstartx(:) = jstart2_all(:)
699             jendx(:) = jend2_all(:)
700             istartx(:) = istart2_all(:)
701             iendx(:) = iend2_all(:)
702        else if (ilayer.eq.2) then
703             kstartx(:) = kstart3_all(:)
704             kendx(:) = kend3_all(:)
705             jstartx(:) = jstart3_all(:)
706             jendx(:) = jend3_all(:)
707             istartx(:) = istart3_all(:)
708             iendx(:) = iend3_all(:)
709        else if (ilayer.eq.3) then
710     !
711     !---------------
712     !EFD extra layer
713     !---------------
714             kstartx(:) = kstart4_all(:)
715             kendx(:) = kend4_all(:)
716             jstartx(:) = jstart4_all(:)
717             jendx(:) = jend4_all(:)
718             istartx(:) = istart4_all(:)
719             iendx(:) = iend4_all(:)
720     
721        endif
722     
723     
724     
725        if (lidebug.ge.1) then
726             call write_debug(name, 'determine send schedule ', myPE )
727        endif
728     
729     !  -----------------------
730     !  determine send schedule
731     !
732     !  examine all neighboring processors
733     !  to see if they need my data
734     !  -----------------------
735     
736     
737     ! -----------------------------------
738     ! first pass to determine array sizes
739     ! -----------------------------------
740     
741     
742     
743        ncount(:) = 0
744     
745        do iproc=0,numPEs-1
746         if (iproc.ne.myPE) then
747     
748     
749             k1 = lbound(ijk2proc,3)
750             k2 = ubound(ijk2proc,3)
751             j1 = lbound(ijk2proc,2)
752             j2 = ubound(ijk2proc,2)
753             i1 = lbound(ijk2proc,1)
754             i2 = ubound(ijk2proc,1)
755     
756     
757             do k=kstartx(iproc),kendx(iproc)
758             do j=jstartx(iproc),jendx(iproc)
759             do i=istartx(iproc),iendx(iproc)
760     
761               ii = imap(i)
762               jj = jmap(j)
763               kk = kmap(k)
764     
765               isvalid  = (k1.le.kk).and.(kk.le.k2)
766               call assert( isvalid, '** sendrecv3_init: invalid kk ', kk )
767     
768               isvalid  = (j1.le.jj).and.(jj.le.j2)
769               call assert( isvalid, '** sendrecv3_init: invalid jj ', jj )
770     
771               isvalid  = (i1.le.ii).and.(ii.le.i2)
772               call assert( isvalid, '** sendrecv3_init: invalid ii ', ii )
773     
774               jproc = ijk2proc( ii,jj,kk )
775     
776               ismine = (jproc .eq. myPE)
777               if (ismine) then
778                     ncount(iproc) = ncount(iproc) + 1
779               endif
780     
781             enddo
782             enddo
783             enddo
784     
785         endif
786        enddo
787     
788     
789     !  --------------
790     !  prepare arrays
791     !  --------------
792        ntotal = 0
793        nsend = 0
794        do iproc=0,numPEs-1
795           ntotal = ntotal + ncount(iproc)
796           if (ncount(iproc).ge.1) then
797             nsend = nsend + 1
798           endif
799        enddo
800     
801        if (lidebug.ge.1) then
802             call write_debug( name, 'ncount = ', ncount )
803             call write_debug( name, 'nsend, ntotal ', nsend, ntotal )
804        endif
805     
806     
807        allocate( xsend(nsend+1) )
808     
809        allocate( sendijk( max(1,ntotal) ) )
810        allocate( sendproc(max(1,nsend)) )
811     
812        nsend = 0
813        do iproc=0,numPEs-1
814          if (ncount(iproc).ne.0) then
815             nsend = nsend + 1
816             sendproc(nsend) = iproc
817          endif
818        enddo
819     
820        xsend(1) = 1
821        do i=1,nsend
822          iproc = sendproc(i)
823          xsend(i+1) = xsend(i) + ncount(iproc)
824        enddo
825     
826        allocate( sendtag( max(1,nsend) ) )
827        do ii=1,nsend
828          iproc = sendproc(ii)
829          src = myPE
830          dest = iproc
831          sendtag(ii) = message_tag( src, dest )
832        enddo
833     
834     ! -----------------------------
835     ! second pass to fill in arrays
836     ! -----------------------------
837     
838     
839       ipos = 1
840       do iter=1,nsend
841          iproc = sendproc(iter)
842          icount = 0
843     
844             do k=kstartx(iproc),kendx(iproc)
845     
846             if (jfastest) then
847     
848                do i=istartx(iproc),iendx(iproc)
849                do j=jstartx(iproc),jendx(iproc)
850     
851     
852                  ii = imap(i)
853                  jj = jmap(j)
854                  kk = kmap(k)
855                  jproc = ijk2proc(ii,jj,kk)
856                  ismine = (jproc.eq.myPE)
857                  if (ismine) then
858                      icount = icount + 1
859                   ijk = funijk(ii,jj,kk)
860     
861                   ipos = xsend(iter)-1 + icount
862                   sendijk( ipos ) = ijk
863                  endif
864     
865               enddo
866               enddo
867     
868     
869             else
870     
871                do j=jstartx(iproc),jendx(iproc)
872                do i=istartx(iproc),iendx(iproc)
873     
874     
875                  ii = imap(i)
876                  jj = jmap(j)
877                  kk = kmap(k)
878                  jproc = ijk2proc(ii,jj,kk)
879                  ismine = (jproc.eq.myPE)
880                  if (ismine) then
881                      icount = icount + 1
882                   ijk = funijk(ii,jj,kk)
883     
884                   ipos = xsend(iter)-1 + icount
885                   sendijk( ipos ) = ijk
886                  endif
887     
888               enddo
889               enddo
890     
891            endif
892            enddo
893     
894          isvalid = (icount .eq. ncount(iproc))
895          call assert( isvalid, &
896                '** sendrecv3_init: icount != ncount(iproc) ', iproc)
897     
898        enddo
899     
900     
901     
902        if (lidebug.ge.1) then
903             call write_debug(name, 'determine recv schedule ', myPE )
904        endif
905     
906     ! ---------------------------
907     ! determine recv schedule
908     !
909     ! examine nodes in my ghost region and
910     ! see what data is needed from my neighbors
911     ! ---------------------------
912     
913     ! -----------------------------------
914     ! first pass to determine array sizes
915     ! -----------------------------------
916     
917       ncount(:) = 0
918     
919             k1 = lbound(ijk2proc,3)
920             k2 = ubound(ijk2proc,3)
921             j1 = lbound(ijk2proc,2)
922             j2 = ubound(ijk2proc,2)
923             i1 = lbound(ijk2proc,1)
924             i2 = ubound(ijk2proc,1)
925     
926     
927             do k=kstartx(myPE),kendx(myPE)
928             do j=jstartx(myPE),jendx(myPE)
929             do i=istartx(myPE),iendx(myPE)
930     
931               ii = imap(i)
932               jj = jmap(j)
933               kk = kmap(k)
934     
935               isvalid  = (k1.le.kk).and.(kk.le.k2)
936               call assert( isvalid, '** sendrecv3_init: invalid kk ', kk )
937     
938               isvalid  = (j1.le.jj).and.(jj.le.j2)
939               call assert( isvalid, '** sendrecv3_init: invalid jj ', jj )
940     
941               isvalid  = (i1.le.ii).and.(ii.le.i2)
942               call assert( isvalid, '** sendrecv3_init: invalid ii ', ii )
943     
944     
945     
946     
947     
948               iproc = ijk2proc(ii,jj,kk)
949               is_halobc = (iproc.eq.-1)
950               ismine = (iproc.eq.myPE)
951               if (.not.ismine) then
952     
953                   isvalid = (0 .le. iproc) .and. &
954                             (iproc.le.numPEs-1) .and. &
955                             (iproc.ne.myPE)
956     
957                   call assert( isvalid, &
958                         '** sendrecv3_init: invalid iproc ',iproc)
959     
960                   ncount(iproc) = ncount(iproc) + 1
961               endif
962             enddo
963             enddo
964             enddo
965     
966        ncount(myPE) = 0
967     
968        ntotal = 0
969        do iproc=0,numPEs-1
970          ntotal = ntotal + ncount(iproc)
971        enddo
972     
973        nrecv = count( ncount(:) .ne. 0)
974     
975        allocate( recvproc( max(1,nrecv) ) )
976     
977        nrecv = 0
978        do iproc=0,numPEs-1
979             if (ncount(iproc).ne.0) then
980                nrecv = nrecv + 1
981                recvproc(nrecv) = iproc
982             endif
983        enddo
984     
985        allocate( xrecv(nrecv+1) )
986        allocate( recvijk(max(1,ntotal)) )
987     
988        xrecv(1) = 1
989        do iter=1,nrecv
990          iproc = recvproc(iter)
991          xrecv(iter+1) = xrecv(iter) + ncount(iproc)
992        enddo
993     
994        allocate( recvtag( max(1,nrecv) ) )
995     
996        do iter=1,nrecv
997           iproc = recvproc(iter)
998           src = iproc
999           dest = myPE
1000           recvtag(iter) = message_tag( src, dest )
1001        enddo
1002     
1003     
1004     ! ----------------------------
1005     ! second pass to fill in array
1006     ! ----------------------------
1007        if (lidebug.ge.1) then
1008             call write_debug( name, 'recv second pass ', myPE )
1009        endif
1010     
1011             ipos = 1
1012     
1013       do iter=1,nrecv
1014           jproc = recvproc(iter)
1015     
1016             do k=kstartx(myPE),kendx(myPE)
1017     
1018             if (jfastest) then
1019     
1020                do i=istartx(myPE),iendx(myPE)
1021                do j=jstartx(myPE),jendx(myPE)
1022     
1023                  ii = imap(i)
1024                  jj = jmap(j)
1025                  kk = kmap(k)
1026     
1027                  iproc = ijk2proc(ii,jj,kk)
1028                  is_halobc = (iproc.eq.-1)
1029     
1030                  ismine = (iproc.eq.myPE)
1031                  if ((.not.ismine) .and. (iproc.eq.jproc)) then
1032     
1033     
1034                     ijk = funijk( i,j,k)
1035                     recvijk( ipos ) = ijk
1036                     ipos = ipos + 1
1037                  endif
1038                enddo
1039                enddo
1040     
1041     
1042     
1043             else
1044     
1045                do j=jstartx(myPE),jendx(myPE)
1046                do i=istartx(myPE),iendx(myPE)
1047     
1048                  ii = imap(i)
1049                  jj = jmap(j)
1050                  kk = kmap(k)
1051     
1052                  iproc = ijk2proc(ii,jj,kk)
1053                  is_halobc = (iproc.eq.-1)
1054     
1055                  ismine = (iproc.eq.myPE)
1056                  if ((.not.ismine) .and. (iproc.eq.jproc)) then
1057     
1058     
1059                     ijk = funijk( i,j,k)
1060                     recvijk( ipos ) = ijk
1061                     ipos = ipos + 1
1062                  endif
1063                enddo
1064                enddo
1065     
1066             endif
1067             enddo
1068     
1069        enddo
1070     
1071         if (ilayer.eq.1) then
1072     
1073             nsend1 = nsend
1074             xsend1 => xsend
1075             sendijk1 => sendijk
1076             sendproc1 => sendproc
1077             sendtag1 => sendtag
1078     
1079             nrecv1 = nrecv
1080             xrecv1 => xrecv
1081             recvijk1 => recvijk
1082             recvproc1 => recvproc
1083             recvtag1 => recvtag
1084     
1085         else if (ilayer.eq.2) then
1086     
1087             nsend2 = nsend
1088             xsend2 => xsend
1089             sendijk2 => sendijk
1090             sendproc2 => sendproc
1091             sendtag2 => sendtag
1092     
1093             nrecv2 = nrecv
1094             xrecv2 => xrecv
1095             recvijk2 => recvijk
1096             recvproc2 => recvproc
1097             recvtag2 => recvtag
1098     
1099         else if (ilayer.eq.3) then
1100     !---------------
1101     !EFD extra layer
1102     !---------------
1103             nsend3 = nsend
1104             xsend3 => xsend
1105             sendijk3 => sendijk
1106             sendproc3 => sendproc
1107             sendtag3 => sendtag
1108     
1109             nrecv3 = nrecv
1110             xrecv3 => xrecv
1111             recvijk3 => recvijk
1112             recvproc3 => recvproc
1113             recvtag3 => recvtag
1114     
1115     
1116         endif
1117     
1118     
1119         nullify( xsend )
1120         nullify( sendijk )
1121         nullify( sendproc )
1122         nullify( sendtag )
1123     
1124         nullify( xrecv )
1125         nullify( recvijk )
1126         nullify( recvproc )
1127         nullify( recvtag )
1128     
1129     enddo ! do ilayer
1130     
1131     
1132             deallocate( ncount )
1133             deallocate( ijk2proc )
1134     
1135             deallocate( istartx )
1136             deallocate( jstartx )
1137             deallocate( kstartx )
1138             deallocate( iendx )
1139             deallocate( jendx )
1140             deallocate( kendx )
1141     
1142             nullify( ncount )
1143             nullify( ijk2proc )
1144     
1145             nullify( istartx )
1146             nullify( jstartx )
1147             nullify( kstartx )
1148             nullify( iendx )
1149             nullify( jendx )
1150             nullify( kendx )
1151     
1152     
1153     
1154         if (lidebug.ge.1) then
1155     
1156             call write_debug( name, ' allocate message buffers ' )
1157             call write_debug( name, 'nrecv1 ', nrecv1 )
1158             call write_debug( name, 'recvproc1 ', recvproc1 )
1159             call write_debug( name, 'recvtag1 ', recvtag1 )
1160             call write_debug( name, 'xrecv1 ', xrecv1 )
1161     
1162     
1163             lmax = size(recvijk1)
1164             allocate( line(lmax) )
1165             line(:) = " "
1166     
1167             ip = 1
1168             do ii=lbound(recvijk1,1),ubound(recvijk1,1)
1169                ijk = recvijk1(ii)
1170                if(localfunc) then
1171                call ijk_of(ijk,i,j,k)
1172                else
1173                i = i_of(ijk)
1174                j = j_of(ijk)
1175                k = k_of(ijk)
1176                endif
1177                write(line(ip),9001) ii,ijk, i,j,k
1178      9001      format('recvijk1( ', i6,') = ', &
1179                          i6, '( ', i6,',',i6,',',i6,') ')
1180                ip = ip + 1
1181             enddo
1182             call write_error( name, line, lmax )
1183             deallocate( line )
1184             nullify( line )
1185     
1186             lmax = size(recvijk2)
1187             allocate( line(lmax) )
1188             line(:) = " "
1189     
1190             ip = 1
1191             do ii=lbound(recvijk2,1),ubound(recvijk2,1)
1192                ijk = recvijk2(ii)
1193                if(localfunc) then
1194                call ijk_of(ijk,i,j,k)
1195                else
1196                i = i_of(ijk)
1197                j = j_of(ijk)
1198                k = k_of(ijk)
1199                endif
1200     
1201                write(line(ip),9101) ii,ijk, i,j,k
1202      9101      format('recvijk2( ', i6,') = ', &
1203                           i6, '( ', i6,',',i6,',',i6,') ')
1204                ip = ip + 1
1205             enddo
1206             call write_error( name, line, lmax )
1207             deallocate( line )
1208             nullify( line )
1209     
1210     
1211             call write_debug( name, ' allocate message buffers ' )
1212             call write_debug( name, 'nsend1 ', nsend1 )
1213             call write_debug( name, 'sendproc1 ', sendproc1 )
1214             call write_debug( name, 'sendtag1 ', sendtag1 )
1215             call write_debug( name, 'xsend1 ', xsend1 )
1216     
1217     
1218     
1219             lmax = size(sendijk1)
1220             allocate(line(lmax))
1221             line(:) = " "
1222     
1223             ip = 1
1224             do ii=lbound(sendijk1,1),ubound(sendijk1,1)
1225                ijk = sendijk1(ii)
1226                if(localfunc) then
1227                call ijk_of(ijk,i,j,k)
1228                else
1229                i = i_of(ijk)
1230                j = j_of(ijk)
1231                k = k_of(ijk)
1232                endif
1233     
1234                write(line(ip),9002) ii,ijk,   i,j,k
1235      9002      format('sendijk1( ', i6,') = ', &
1236                             i6, '( ', i6,',',i6,',',i6,') ')
1237                ip = ip + 1
1238             enddo
1239     
1240             call write_error( name, line, lmax )
1241             deallocate( line )
1242             nullify( line )
1243     
1244     
1245             lmax = size(sendijk2)
1246             allocate(line(lmax))
1247             line(:) = " "
1248     
1249             ip = 1
1250             do ii=lbound(sendijk2,1),ubound(sendijk2,1)
1251                ijk = sendijk2(ii)
1252                if(localfunc) then
1253                call ijk_of(ijk,i,j,k)
1254                else
1255                i = i_of(ijk)
1256                j = j_of(ijk)
1257                k = k_of(ijk)
1258                endif
1259     
1260                write(line(ip),9102) ii,ijk,   i,j,k
1261      9102      format('sendijk2( ', i6,') = ', &
1262                                i6, '( ', i6,',',i6,',',i6,') ')
1263                ip = ip + 1
1264             enddo
1265     
1266             call write_error( name, line, lmax )
1267             deallocate( line )
1268             nullify( line )
1269     
1270     
1271     
1272         endif
1273     
1274     
1275     
1276     ! ------------------------
1277     ! allocate message buffers
1278     ! ------------------------
1279     
1280     
1281         isize = max(1,max(nsend1,nsend2))
1282         allocate( sendrequest( isize ) )
1283         allocate( send_persistent_request1( isize ) )
1284         allocate( send_persistent_request2( isize ) )
1285     
1286         isize = max(1,max(nrecv1,nrecv2))
1287         allocate( recvrequest( isize ) )
1288         allocate( recv_persistent_request1( isize ) )
1289         allocate( recv_persistent_request2( isize ) )
1290     
1291     ! -----------------------------------
1292     ! preallocate buffers for common case
1293     ! -----------------------------------
1294             recvsize1 = xrecv1( nrecv1+1)-1
1295             recvsize2 = xrecv2( nrecv2+1)-1
1296     
1297             isize = max(1,max(recvsize1,recvsize2))
1298             allocate( drecvbuffer( isize ) )
1299     
1300             sendsize1 = xsend1( nsend1+1)-1
1301             sendsize2 = xsend2( nsend2+1)-1
1302     
1303             isize = max(1,max(sendsize1,sendsize2))
1304             allocate( dsendbuffer( isize ) )
1305     
1306     
1307      if (use_persistent_message) then
1308     
1309        datatype = MPI_DOUBLE_PRECISION
1310     
1311        do layer=nlayers,nlayers
1312     
1313     
1314             if (layer.eq.1) then
1315               nrecv = nrecv1
1316               recvtag =>recvtag1
1317               recvproc => recvproc1
1318               recvijk => recvijk1
1319               xrecv => xrecv1
1320     
1321               nsend = nsend1
1322               sendtag => sendtag1
1323               sendproc => sendproc1
1324               sendijk => sendijk1
1325               xsend => xsend1
1326     
1327               send_persistent_request => send_persistent_request1
1328               recv_persistent_request => recv_persistent_request1
1329     
1330             else if (layer.eq.2) then
1331               nrecv = nrecv2
1332               recvtag =>recvtag2
1333               recvproc => recvproc2
1334               recvijk => recvijk2
1335               xrecv => xrecv2
1336     
1337               nsend = nsend2
1338               sendtag => sendtag2
1339               sendproc => sendproc2
1340               sendijk => sendijk2
1341               xsend => xsend2
1342     
1343               send_persistent_request => send_persistent_request2
1344               recv_persistent_request => recv_persistent_request2
1345     
1346             else if (layer.eq.3) then
1347     !----------------
1348     !EFD extra layer
1349     !----------------
1350     
1351               nrecv = nrecv3
1352               recvtag =>recvtag3
1353               recvproc => recvproc3
1354               recvijk => recvijk3
1355               xrecv => xrecv3
1356     
1357               nsend = nsend3
1358               sendtag => sendtag3
1359               sendproc => sendproc3
1360               sendijk => sendijk3
1361               xsend => xsend3
1362     
1363               send_persistent_request => send_persistent_request3
1364               recv_persistent_request => recv_persistent_request3
1365     
1366             endif
1367     
1368     
1369     
1370                  do ii=1,nrecv
1371                     j1 = xrecv(ii)
1372                     j2 = xrecv(ii+1)-1
1373                     icount = j2-j1+1
1374                     source = recvproc( ii )
1375                     tag = recvtag( ii )
1376     
1377     
1378                     if (lidebug.ge.2) then
1379     
1380                          call write_debug(name, 'mpi_recv_init: ii,j1,j2 ', &
1381                                                  ii,j1,j2 )
1382                          call write_debug(name, 'icount, source, tag ', &
1383                                                  icount,source,tag )
1384                     endif
1385     
1386     
1387                     call MPI_RECV_INIT( drecvbuffer(j1), icount, datatype, &
1388                                             source, tag, comm, request, ierror )
1389                     call MPI_Check( 'sendrecv3_begin_1d:MPI_IRECV ', ierror )
1390     
1391                     recv_persistent_request(ii) = request
1392                  enddo
1393     
1394     
1395                  do ii=1,nsend
1396                      j1 = xsend(ii)
1397                      j2 = xsend(ii+1)-1
1398                       dest = sendproc( ii )
1399                       tag = sendtag( ii )
1400                       icount = j2-j1+1
1401     
1402                       if (lidebug.ge.2) then
1403     
1404                           call write_debug(name, 'mpi_send_init: ii,j1,j2 ', &
1405                                              ii,j1,j2)
1406                           call write_debug(name, 'icount, dest, tag ', &
1407                                              icount,dest,tag )
1408                       endif
1409     
1410     
1411                       call MPI_SEND_INIT( dsendbuffer(j1), icount, datatype, &
1412                                             dest, tag, &
1413                                             comm, request, ierror )
1414                      call MPI_Check( 'sendrecv3_begin_1d:MPI_SEND_INIT ', ierror )
1415     
1416                      send_persistent_request( ii ) = request
1417                   enddo
1418              enddo
1419     
1420             endif
1421     
1422             if (lidebug.ge.1) then
1423                call write_debug(name, ' end of sendrecv3_init ', myPE )
1424             endif
1425     
1426             end subroutine sendrecv3_init
1427     
1428     
1429     
1430             subroutine sendrecv3_begin_1d( XX, ilayer, idebug )
1431     
1432             use functions
1433     
1434             implicit none
1435     
1436             integer, intent(in),optional :: ilayer
1437             double precision, intent(inout), dimension(:) :: XX
1438             integer, intent(in), optional :: idebug
1439     
1440     !       interface
1441     !
1442     !       subroutine MPI_ISEND( buffer, count, datatype, dest, tag, &
1443     !                        comm, request, ierror )
1444     !       double precision buffer(*)
1445     !       integer count,datatype,dest,tag,comm,request,ierror
1446     !       end subroutine MPI_ISEND
1447     !
1448     !        subroutine MPI_IRECV( buffer, count, datatype, source, tag, &
1449     !                       comm, request, ierror )
1450     !       double precision buffer(*)
1451     !       integer count,datatype,source,tag,comm,request,ierror
1452     !       end subroutine MPI_IRECV
1453     !
1454     !       end interface
1455     
1456     
1457     
1458     
1459     !       ---------------
1460     !       local variables
1461     !       ---------------
1462             character(len=80), parameter :: name = 'sendrecv3_begin_1d'
1463     
1464             integer :: lidebug
1465     
1466             integer ::  layer, datatype, comm, recvsize, sendsize, &
1467                     ijk,jj,j1,j2, request, ii,count,source,dest, tag, ierror
1468     
1469              lidebug = 0
1470     
1471             if (present(idebug)) then
1472                lidebug = idebug
1473             endif
1474     
1475             layer = 1
1476             if (present(ilayer)) then
1477                layer = ilayer
1478             endif
1479     
1480             if (layer.eq.1) then
1481               nrecv = nrecv1
1482               recvtag =>recvtag1
1483               recvproc => recvproc1
1484               recvijk => recvijk1
1485               xrecv => xrecv1
1486     
1487               nsend = nsend1
1488               sendtag => sendtag1
1489               sendproc => sendproc1
1490               sendijk => sendijk1
1491               xsend => xsend1
1492     
1493               send_persistent_request => send_persistent_request1
1494               recv_persistent_request => recv_persistent_request1
1495     
1496             else if (layer.eq.2) then
1497               nrecv = nrecv2
1498               recvtag =>recvtag2
1499               recvproc => recvproc2
1500               recvijk => recvijk2
1501               xrecv => xrecv2
1502     
1503               nsend = nsend2
1504               sendtag => sendtag2
1505               sendproc => sendproc2
1506               sendijk => sendijk2
1507               xsend => xsend2
1508     
1509               send_persistent_request => send_persistent_request2
1510               recv_persistent_request => recv_persistent_request2
1511     
1512             else if (layer.eq.3) then
1513     !---------------
1514     !EFD extra layer
1515     !---------------
1516     
1517               nrecv = nrecv3
1518               recvtag =>recvtag3
1519               recvproc => recvproc3
1520               recvijk => recvijk3
1521               xrecv => xrecv3
1522     
1523               nsend = nsend3
1524               sendtag => sendtag3
1525               sendproc => sendproc3
1526               sendijk => sendijk3
1527               xsend => xsend3
1528     
1529               send_persistent_request => send_persistent_request3
1530               recv_persistent_request => recv_persistent_request3
1531     
1532             endif
1533     
1534     
1535     !   --------------------------
1536     !   post asynchronous receives
1537     !   --------------------------
1538     
1539         if (lidebug.ge.1) then
1540             call write_debug(name, 'post asynchronous receives, nrecv = ', nrecv )
1541         endif
1542     
1543         if (nrecv.ge.1) then
1544             recvsize = xrecv( nrecv+1)-1
1545     
1546     
1547             if (lidebug.ge.1) then
1548                     call write_debug( name, 'recvsize, ubound(drecvbuffer,1) ', &
1549                             recvsize, ubound(drecvbuffer,1) )
1550     
1551                     call write_debug( name, 'ubound(xrecv,1) ', &
1552                                             ubound(xrecv,1) )
1553                     call write_debug( name, 'ubound(recvproc,1) ', &
1554                                             ubound(recvproc,1) )
1555                     call write_debug( name, 'ubound(recvtag,1) ', &
1556                                             ubound(recvtag,1) )
1557     
1558             endif
1559     
1560     !       -------------
1561     !       post receives
1562     !       -------------
1563             datatype = MPI_DOUBLE_PRECISION
1564             comm = communicator
1565     
1566          if (use_persistent_message) then
1567     
1568     !          ---------------------------------------
1569     !          persistent request already established
1570     !          ---------------------------------------
1571                if (lidebug.ge.2) then
1572                    call write_debug( name,'before startall for recv ',&
1573                             recv_persistent_request)
1574                endif
1575     
1576                call MPI_STARTALL( nrecv, recv_persistent_request, ierror )
1577     
1578                if (lidebug.ge.2) then
1579                    call write_debug( name,'after startall for recv, ierror',&
1580                             ierror)
1581                endif
1582     
1583                call MPI_Check( 'sendrecv3_begin: MPI_STARTALL ', ierror )
1584     
1585     
1586          else
1587     !    ---------
1588     !    use irecv
1589     !    ---------
1590     
1591             do ii=1,nrecv
1592                j1 = xrecv(ii)
1593                j2 = xrecv(ii+1)-1
1594                count = j2-j1+1
1595                source = recvproc( ii )
1596                tag = recvtag( ii )
1597     
1598     
1599                if (lidebug.ge.2) then
1600     
1601                     call write_debug(name, 'mpi_irecv: ii,j1,j2 ', ii,j1,j2 )
1602                     call write_debug(name, 'count, source, tag ', &
1603                                             count,source,tag )
1604                endif
1605     
1606                call MPI_IRECV( drecvbuffer(j1), count, datatype, source, tag, &
1607                             comm, request, ierror )
1608     
1609                call MPI_Check( 'sendrecv3_begin_1d:MPI_IRECV ', ierror )
1610     
1611                 recvrequest( ii ) = request
1612             enddo
1613           endif
1614     
1615        endif
1616     
1617     !  -----------------------
1618     !  post asynchronous sends
1619     !  -----------------------
1620     
1621        if (lidebug.ge.1) then
1622     
1623             call write_debug(name, 'post asynchronous sends ')
1624        endif
1625     
1626        if (nsend.ge.1) then
1627             sendsize = xsend( nsend+1)-1
1628     
1629     
1630             if (lidebug.ge.1) then
1631     
1632                     call write_debug( name, &
1633                             'sendsize, ubound(dsendbuffer,1) ', &
1634                             sendsize, ubound(dsendbuffer,1) )
1635     
1636                     call write_debug( name, 'ubound(xsend,1) ', &
1637                                             ubound(xsend,1) )
1638                     call write_debug( name, 'ubound(sendproc,1) ', &
1639                                             ubound(sendproc,1) )
1640                     call write_debug( name, 'ubound(sendtag,1) ', &
1641                                             ubound(sendtag,1) )
1642     
1643             endif
1644     
1645     
1646     
1647     
1648     !       -------------
1649     !       perform sends
1650     !       -------------
1651             datatype = MPI_DOUBLE_PRECISION
1652             comm = communicator
1653     
1654         if (use_persistent_message) then
1655     
1656     !       -----------------------------
1657     !       perform copy into dsendbuffer
1658     !       -----------------------------
1659             j1 = xsend(1)
1660             j2 = xsend(nsend+1)-1
1661     
1662             do jj=j1,j2
1663                   ijk = sendijk( jj )
1664                   dsendbuffer( jj )  = XX(ijk)
1665             enddo
1666     
1667     
1668     
1669              if (lidebug.ge.2) then
1670                 call write_debug(name,'before mpi_startall send ',&
1671                     send_persistent_request )
1672              endif
1673     
1674              call MPI_STARTALL( nsend, send_persistent_request, ierror )
1675     
1676              if (lidebug .ge.2) then
1677                call write_debug(name,'after mpi_startall send ',&
1678                     send_persistent_request )
1679              endif
1680     
1681              call MPI_Check( 'sendrecv3_begin_1d:MPI_STARTALL ', ierror )
1682     
1683     
1684         else
1685     
1686             do ii=1,nsend
1687     
1688     !       ----------------------------
1689     !       perform copy into dsendbuffer
1690     !       ----------------------------
1691     
1692                 j1 = xsend(ii)
1693                 j2 = xsend(ii+1)-1
1694                 count = j2-j1+1
1695     
1696                 do jj=j1,j2
1697                     ijk = sendijk( jj )
1698                     dsendbuffer(jj) = XX(ijk)
1699                 enddo
1700     
1701                 dest = sendproc( ii )
1702                 tag = sendtag( ii )
1703     
1704                 if (lidebug.ge.2) then
1705     
1706                     call write_debug(name, 'mpi_isend: ii,j1,j2 ', ii,j1,j2)
1707                     call write_debug(name, 'count, dest, tag ', count,dest,tag )
1708                 endif
1709     
1710     
1711                 call MPI_ISEND( dsendbuffer(j1), count, datatype, dest, tag, &
1712                             comm, request, ierror )
1713                call MPI_Check( 'sendrecv3_begin_1d:MPI_ISEND ', ierror )
1714     
1715                sendrequest( ii ) = request
1716             enddo
1717     
1718          endif
1719     
1720         endif
1721     
1722             return
1723             end subroutine sendrecv3_begin_1d
1724     
1725     
1726             subroutine sendrecv3_begin_1i( XX, ilayer, idebug )
1727     
1728             use functions
1729     
1730             implicit none
1731     
1732             integer, intent(in),optional :: ilayer
1733             integer, intent(inout), dimension(:) :: XX
1734             integer, intent(in), optional :: idebug
1735     
1736     !       interface
1737     !
1738     !       subroutine MPI_ISEND( buffer, count, datatype, dest, tag, &
1739     !                        comm, request, ierror )
1740     !       integer buffer(*)
1741     !       integer count,datatype,dest,tag,comm,request,ierror
1742     !       end subroutine MPI_ISEND
1743     !
1744     !        subroutine MPI_IRECV( buffer, count, datatype, source, tag, &
1745     !                       comm, request, ierror )
1746     !       integer buffer(*)
1747     !       integer count,datatype,source,tag,comm,request,ierror
1748     !       end subroutine MPI_IRECV
1749     !
1750     !       end interface
1751     
1752     
1753     
1754     
1755     !       ---------------
1756     !       local variables
1757     !       ---------------
1758             character(len=80), parameter :: name = 'sendrecv3_begin_1i'
1759     
1760             integer :: lidebug
1761     
1762             integer ::  layer, datatype, comm, recvsize, sendsize, &
1763                     ijk,jj,j1,j2, request, ii,count,source,dest, tag, ierror
1764     
1765             lidebug = 0
1766             if (present(idebug)) then
1767                lidebug = idebug
1768             endif
1769     
1770             layer = 1
1771             if (present(ilayer)) then
1772                layer = ilayer
1773             endif
1774     
1775             if (layer.eq.1) then
1776               nrecv = nrecv1
1777               recvtag =>recvtag1
1778               recvproc => recvproc1
1779               recvijk => recvijk1
1780               xrecv => xrecv1
1781     
1782               nsend = nsend1
1783               sendtag => sendtag1
1784               sendproc => sendproc1
1785               sendijk => sendijk1
1786               xsend => xsend1
1787             else if (layer.eq.2) then
1788               nrecv = nrecv2
1789               recvtag =>recvtag2
1790               recvproc => recvproc2
1791               recvijk => recvijk2
1792               xrecv => xrecv2
1793     
1794               nsend = nsend2
1795               sendtag => sendtag2
1796               sendproc => sendproc2
1797               sendijk => sendijk2
1798               xsend => xsend2
1799             else if (layer.eq.3) then
1800     !---------------
1801     !EFD extra layer
1802     !---------------
1803               nrecv = nrecv3
1804               recvtag =>recvtag3
1805               recvproc => recvproc3
1806               recvijk => recvijk3
1807               xrecv => xrecv3
1808     
1809               nsend = nsend3
1810               sendtag => sendtag3
1811               sendproc => sendproc3
1812               sendijk => sendijk3
1813               xsend => xsend3
1814             endif
1815     
1816     
1817     !   --------------------------
1818     !   post asynchronous receives
1819     !   --------------------------
1820     
1821         if (lidebug.ge.1) then
1822             call write_debug(name, &
1823                   'post asynchronous receives, nrecv = ', nrecv )
1824         endif
1825     
1826         if (nrecv.ge.1) then
1827             recvsize = xrecv( nrecv+1)-1
1828             allocate( irecvbuffer( recvsize ) )
1829     
1830             if (lidebug.ge.1) then
1831                     call write_debug( name, &
1832                             'recvsize, ubound(irecvbuffer,1) ', &
1833                             recvsize, ubound(irecvbuffer,1) )
1834     
1835                     call write_debug( name, 'ubound(xrecv,1) ', &
1836                                             ubound(xrecv,1) )
1837                     call write_debug( name, 'ubound(recvproc,1) ', &
1838                                             ubound(recvproc,1) )
1839                     call write_debug( name, 'ubound(recvtag,1) ', &
1840                                             ubound(recvtag,1) )
1841     
1842             endif
1843     
1844     !       -------------
1845     !       post receives
1846     !       -------------
1847             datatype = MPI_INTEGER
1848             comm = communicator
1849     
1850             do ii=1,nrecv
1851                j1 = xrecv(ii)
1852                j2 = xrecv(ii+1)-1
1853                count = j2-j1+1
1854                source = recvproc( ii )
1855                tag = recvtag( ii )
1856     
1857     
1858                if (lidebug.ge.2) then
1859     
1860                     call write_debug(name, 'mpi_irecv: ii,j1,j2 ', ii,j1,j2 )
1861                     call write_debug(name, 'count, source, tag ', &
1862                                             count,source,tag )
1863                endif
1864     
1865                call MPI_IRECV( irecvbuffer(j1), count, datatype, source, tag, &
1866                             comm, request, ierror )
1867     
1868                call MPI_Check( 'sendrecv3_begin_1i:MPI_IRECV ', ierror )
1869     
1870                 recvrequest( ii ) = request
1871             enddo
1872     
1873        endif
1874     
1875     !  -----------------------
1876     !  post asynchronous sends
1877     !  -----------------------
1878     
1879        if (lidebug.ge.1) then
1880     
1881             call write_debug(name, 'post asynchronous sends ')
1882        endif
1883     
1884        if (nsend.ge.1) then
1885             sendsize = xsend( nsend+1)-1
1886             allocate( isendbuffer( sendsize ) )
1887     
1888             if (lidebug.ge.1) then
1889     
1890                     call write_debug( name, 'sendsize, ubound(isendbuffer,1) ', &
1891                             sendsize, ubound(isendbuffer,1) )
1892     
1893                     call write_debug( name, 'ubound(xsend,1) ', &
1894                                             ubound(xsend,1) )
1895                     call write_debug( name, 'ubound(sendproc,1) ', &
1896                                             ubound(sendproc,1) )
1897                     call write_debug( name, 'ubound(sendtag,1) ', &
1898                                             ubound(sendtag,1) )
1899     
1900             endif
1901     
1902     
1903     
1904     
1905     !       -------------
1906     !       perform sends
1907     !       -------------
1908             datatype = MPI_INTEGER
1909             comm = communicator
1910     
1911             do ii=1,nsend
1912     
1913     !       ----------------------------
1914     !       perform copy into sendbuffer
1915     !       ----------------------------
1916     
1917                 j1 = xsend(ii)
1918                 j2 = xsend(ii+1)-1
1919                 count = j2-j1+1
1920     
1921                 do jj=j1,j2
1922                     ijk = sendijk( jj )
1923                     isendbuffer(jj) = XX(ijk)
1924                 enddo
1925     
1926                 dest = sendproc( ii )
1927                 tag = sendtag( ii )
1928     
1929                 if (lidebug.ge.2) then
1930     
1931                     call write_debug(name, 'mpi_isend: ii,j1,j2 ', ii,j1,j2)
1932                     call write_debug(name, 'count, dest, tag ', count,dest,tag )
1933                 endif
1934     
1935     
1936                 call MPI_ISEND( isendbuffer(j1), count, datatype, dest, tag, &
1937                             comm, request, ierror )
1938                call MPI_Check( 'sendrecv3_begin_1i:MPI_ISEND ', ierror )
1939     
1940                sendrequest( ii ) = request
1941             enddo
1942     
1943         endif
1944     
1945             return
1946             end subroutine sendrecv3_begin_1i
1947     
1948     
1949             subroutine sendrecv3_begin_1c( XX, ilayer, idebug )
1950     
1951               use functions
1952     
1953             implicit none
1954     
1955             integer, intent(in),optional :: ilayer
1956             character(len=*), intent(inout), dimension(:) :: XX
1957             integer, intent(in), optional :: idebug
1958     
1959     !       interface
1960     !
1961     !       subroutine MPI_ISEND( buffer, count, datatype, dest, tag, &
1962     !                        comm, request, ierror )
1963     !       character(len=*) buffer(*)
1964     !       integer count,datatype,dest,tag,comm,request,ierror
1965     !       end subroutine MPI_ISEND
1966     !
1967     !        subroutine MPI_IRECV( buffer, count, datatype, source, tag, &
1968     !                       comm, request, ierror )
1969     !       character(len=*) buffer(*)
1970     !       integer count,datatype,source,tag,comm,request,ierror
1971     !       end subroutine MPI_IRECV
1972     !
1973     !       end interface
1974     
1975     
1976     
1977     
1978     !       ---------------
1979     !       local variables
1980     !       ---------------
1981             character(len=80), parameter :: name = 'sendrecv3_begin_1c'
1982     
1983             integer :: lidebug
1984     
1985             integer ::  layer, datatype, comm, recvsize, sendsize, &
1986                     ijk,jj,j1,j2, request, ii,count,source,dest, tag, ierror
1987     
1988             integer :: ic, clen, jpos
1989     
1990             lidebug = 0
1991             if (present(idebug)) then
1992                lidebug = idebug
1993             endif
1994     
1995             layer = 1
1996             if (present(ilayer)) then
1997                layer = ilayer
1998             endif
1999     
2000             jpos = lbound(XX,1)
2001             clen = len( XX( jpos ) )
2002     
2003             if (layer.eq.1) then
2004               nrecv = nrecv1
2005               recvtag =>recvtag1
2006               recvproc => recvproc1
2007               recvijk => recvijk1
2008               xrecv => xrecv1
2009     
2010               nsend = nsend1
2011               sendtag => sendtag1
2012               sendproc => sendproc1
2013               sendijk => sendijk1
2014               xsend => xsend1
2015             else if (layer.eq.2) then
2016               nrecv = nrecv2
2017               recvtag =>recvtag2
2018               recvproc => recvproc2
2019               recvijk => recvijk2
2020               xrecv => xrecv2
2021     
2022               nsend = nsend2
2023               sendtag => sendtag2
2024               sendproc => sendproc2
2025               sendijk => sendijk2
2026               xsend => xsend2
2027             else if (layer.eq.3) then
2028     !----------------
2029     !EFD extra layer
2030     !----------------
2031     
2032               nrecv = nrecv3
2033               recvtag =>recvtag3
2034               recvproc => recvproc3
2035               recvijk => recvijk3
2036               xrecv => xrecv3
2037     
2038               nsend = nsend3
2039               sendtag => sendtag3
2040               sendproc => sendproc3
2041               sendijk => sendijk3
2042               xsend => xsend3
2043             endif
2044     
2045     
2046     !   --------------------------
2047     !   post asynchronous receives
2048     !   --------------------------
2049     
2050         if (lidebug.ge.1) then
2051             call write_debug(name, 'post asynchronous receives, nrecv = ', nrecv )
2052         endif
2053     
2054         if (nrecv.ge.1) then
2055             recvsize = xrecv( nrecv+1)-1
2056     
2057             allocate( crecvbuffer( recvsize*clen ) )
2058     
2059             if (lidebug.ge.1) then
2060                     call write_debug( name, 'recvsize, ubound(crecvbuffer,1) ', &
2061                             recvsize, ubound(crecvbuffer,1) )
2062     
2063                     call write_debug( name, 'ubound(xrecv,1) ', &
2064                                             ubound(xrecv,1) )
2065                     call write_debug( name, 'ubound(recvproc,1) ', &
2066                                             ubound(recvproc,1) )
2067                     call write_debug( name, 'ubound(recvtag,1) ', &
2068                                             ubound(recvtag,1) )
2069     
2070             endif
2071     
2072     !       -------------
2073     !       post receives
2074     !       -------------
2075             datatype = MPI_CHARACTER
2076             comm = communicator
2077     
2078             do ii=1,nrecv
2079                j1 = xrecv(ii)
2080                j2 = xrecv(ii+1)-1
2081     
2082                count = j2-j1+1
2083                count = count*clen
2084     
2085                source = recvproc( ii )
2086                tag = recvtag( ii )
2087     
2088     
2089                if (lidebug.ge.2) then
2090     
2091                     call write_debug(name, 'mpi_irecv: ii,j1,j2 ', ii,j1,j2 )
2092                     call write_debug(name, 'count, source, tag ', &
2093                                             count,source,tag )
2094                endif
2095     
2096                jpos = 1 + (j1-1)*clen
2097                call MPI_IRECV( crecvbuffer(jpos), count, datatype, source, tag, &
2098                             comm, request, ierror )
2099     
2100                call MPI_Check( 'sendrecv3_begin_1c:MPI_IRECV ', ierror )
2101     
2102                 recvrequest( ii ) = request
2103             enddo
2104     
2105        endif
2106     
2107     !  -----------------------
2108     !  post asynchronous sends
2109     !  -----------------------
2110     
2111        if (lidebug.ge.1) then
2112     
2113             call write_debug(name, 'post asynchronous sends ')
2114        endif
2115     
2116        if (nsend.ge.1) then
2117             sendsize = xsend( nsend+1)-1
2118     
2119             allocate( csendbuffer( sendsize*clen ) )
2120     
2121             if (lidebug.ge.1) then
2122     
2123                     call write_debug( name, 'sendsize, ubound(csendbuffer,1) ', &
2124                             sendsize, ubound(csendbuffer,1) )
2125     
2126                     call write_debug( name, 'ubound(xsend,1) ', &
2127                                             ubound(xsend,1) )
2128                     call write_debug( name, 'ubound(sendproc,1) ', &
2129                                             ubound(sendproc,1) )
2130                     call write_debug( name, 'ubound(sendtag,1) ', &
2131                                             ubound(sendtag,1) )
2132     
2133             endif
2134     
2135     
2136     
2137     
2138     !       -------------
2139     !       perform sends
2140     !       -------------
2141             datatype = MPI_CHARACTER
2142             comm = communicator
2143     
2144             do ii=1,nsend
2145     
2146     !       ----------------------------
2147     !       perform copy into sendbuffer
2148     !       ----------------------------
2149     
2150                 j1 = xsend(ii)
2151                 j2 = xsend(ii+1)-1
2152     
2153     
2154                 count = j2-j1+1
2155                 count = count*clen
2156     
2157                 do jj=j1,j2
2158                     ijk = sendijk( jj )
2159                     do ic=1,clen
2160                         jpos = (jj-1)*clen + ic
2161                         csendbuffer(jpos) = XX(ijk)(ic:ic)
2162                     enddo
2163                 enddo
2164     
2165                 dest = sendproc( ii )
2166                 tag = sendtag( ii )
2167     
2168                 if (lidebug.ge.2) then
2169     
2170                     call write_debug(name, 'mpi_isend: ii,j1,j2 ', ii,j1,j2)
2171                     call write_debug(name, 'count, dest, tag ', count,dest,tag )
2172                 endif
2173     
2174     
2175                 jpos = (j1-1)*clen + 1
2176                 call MPI_ISEND( csendbuffer(jpos), count, datatype, dest, tag, &
2177                             comm, request, ierror )
2178                call MPI_Check( 'sendrecv3_begin_1c:MPI_ISEND ', ierror )
2179     
2180                sendrequest( ii ) = request
2181             enddo
2182     
2183         endif
2184     
2185             return
2186             end subroutine sendrecv3_begin_1c
2187     
2188     
2189             subroutine sendrecv3_end_1d( XX, idebug )
2190     
2191               use functions
2192     
2193             implicit none
2194     
2195             double precision, intent(inout), dimension(:) :: XX
2196             integer, intent(in), optional :: idebug
2197     
2198             interface
2199     
2200             subroutine MPI_WAITANY(count, array_of_requests, jindex, &
2201                                     status, ierror)
2202             use mpi, only: MPI_STATUS_SIZE
2203     
2204             integer count
2205             integer array_of_requests(*)
2206             integer jindex
2207             integer status(MPI_STATUS_SIZE)
2208             integer ierror
2209             end subroutine MPI_WAITANY
2210     
2211             subroutine MPI_WAITALL( count, array_of_requests,  &
2212                     array_of_status, ierror )
2213             use mpi, only: MPI_STATUS_SIZE
2214     
2215             integer count
2216             integer array_of_requests(*)
2217             integer array_of_status( MPI_STATUS_SIZE,*)
2218             integer ierror
2219             end subroutine MPI_WAITALL
2220     
2221             end interface
2222     
2223     !       ---------------
2224     !       local variables
2225     !       ---------------
2226             character(len=80), parameter :: name = 'sendrecv3_end_1d'
2227     
2228             logical, parameter :: use_waitany = .false.
2229     
2230             integer :: lidebug
2231             integer :: jj,ijk,  jindex, ii,j1,j2, ierror
2232     
2233             integer, dimension(MPI_STATUS_SIZE) :: recv_status
2234             integer, dimension(:,:), pointer :: send_status
2235     
2236     !       --------------------------
2237     !       wait for sends to complete
2238     !       --------------------------
2239     
2240             lidebug = 0
2241             if (present(idebug)) then
2242                lidebug = idebug
2243             endif
2244     
2245             if (nsend .ge.1) then
2246     
2247             if (lidebug.ge.1) then
2248     
2249                call write_debug(name, &
2250                     'waiting for sends to complete, nsend  = ', nsend )
2251             endif
2252     
2253             allocate( send_status(MPI_STATUS_SIZE,nsend))
2254     
2255             if (use_persistent_message) then
2256                call MPI_WAITALL( nsend, send_persistent_request, &
2257                                         send_status, ierror )
2258             else
2259                call MPI_WAITALL( nsend, sendrequest, send_status, ierror )
2260             endif
2261     
2262             call MPI_Check( 'sendrecv3_end_1d:MPI_WAITALL ', ierror )
2263     
2264             deallocate( send_status )
2265             nullify( send_status )
2266     
2267     
2268             endif
2269     
2270     
2271     !       --------------------------
2272     !       wait for recvs to complete
2273     !       --------------------------
2274             if (nrecv.ge.1) then
2275     
2276             if (lidebug.ge.1) then
2277     
2278                call write_debug( name, &
2279                     'waiting for receives to complete, nrecv =  ', nrecv )
2280             endif
2281     
2282           if (use_waitany) then
2283             do ii=1,nrecv
2284     
2285                if (use_persistent_message) then
2286                   call MPI_WAITANY( nrecv, recv_persistent_request,  &
2287                             jindex, recv_status, ierror )
2288                else
2289                   call MPI_WAITANY( nrecv, recvrequest,   &
2290                             jindex, recv_status, ierror )
2291                endif
2292     
2293                call MPI_Check( 'sendrecv3_end_1d:MPI_WAITANY ', ierror )
2294     
2295                j1 = xrecv( jindex )
2296                j2 = xrecv( jindex + 1)-1
2297     
2298                if (lidebug.ge.2) then
2299                     call write_debug(name, 'jindex, j1,j2 ', jindex,j1,j2 )
2300                endif
2301     
2302                do jj=j1,j2
2303                  ijk = recvijk( jj )
2304                  XX(ijk) = drecvbuffer(jj)
2305                enddo
2306             enddo
2307           else
2308                if (use_persistent_message) then
2309                   call MPI_WAITALL( nrecv, recv_persistent_request, recv_status, ierror )
2310                else
2311                   call MPI_WAITALL( nrecv, recvrequest, recv_status, ierror )
2312                endif
2313                call MPI_Check( 'sendrecv3_end_1d:MPI_WAITALL recv ', ierror )
2314     
2315                j1 = xrecv(1)
2316                j2 = xrecv( nrecv +1)-1
2317                do jj=j1,j2
2318                    ijk = recvijk( jj )
2319                    XX(ijk) = drecvbuffer(jj)
2320                enddo
2321             endif
2322     
2323     
2324             endif
2325     
2326             return
2327             end subroutine sendrecv3_end_1d
2328     
2329     
2330             subroutine sendrecv3_end_1c( XX, idebug )
2331     
2332               use functions
2333     
2334             implicit none
2335     
2336             character(len=*), intent(inout), dimension(:) :: XX
2337             integer, intent(in), optional :: idebug
2338     
2339             interface
2340     
2341             subroutine MPI_WAITANY(count, array_of_requests, jindex, &
2342                                     status, ierror)
2343             use mpi, only: MPI_STATUS_SIZE
2344     
2345             integer count
2346             integer array_of_requests(*)
2347             integer jindex
2348             integer status(MPI_STATUS_SIZE)
2349             integer ierror
2350             end subroutine MPI_WAITANY
2351     
2352             subroutine MPI_WAITALL( count, array_of_requests,  &
2353                     array_of_status, ierror )
2354             use mpi, only: MPI_STATUS_SIZE
2355     
2356             integer count
2357             integer array_of_requests(*)
2358             integer array_of_status( MPI_STATUS_SIZE,*)
2359             integer ierror
2360             end subroutine MPI_WAITALL
2361     
2362             end interface
2363     
2364     !       ---------------
2365     !       local variables
2366     !       ---------------
2367             character(len=80), parameter :: name = 'sendrecv3_end_1c'
2368     
2369             integer :: ic, clen, jpos
2370     
2371             logical, parameter :: use_waitany = .false.
2372     
2373             integer :: lidebug
2374             integer :: jj,ijk,  jindex, ii,j1,j2, ierror
2375     
2376             integer, dimension(MPI_STATUS_SIZE) :: recv_status
2377             integer, dimension(:,:), pointer :: send_status
2378     
2379     !       --------------------------
2380     !       wait for sends to complete
2381     !       --------------------------
2382     
2383             lidebug = 0
2384             if (present(idebug)) then
2385                lidebug = idebug
2386             endif
2387     
2388             jpos = lbound(XX,1)
2389             clen = len(XX(jpos))
2390     
2391             if (nsend .ge.1) then
2392     
2393             if (lidebug.ge.1) then
2394     
2395                call write_debug(name, &
2396                     'waiting for sends to complete, nsend  = ', nsend )
2397             endif
2398     
2399     
2400             allocate( send_status(MPI_STATUS_SIZE,nsend))
2401     
2402             call MPI_WAITALL( nsend, sendrequest, send_status, ierror )
2403             call MPI_Check( 'sendrecv3_end_1c:MPI_WAITALL ', ierror )
2404     
2405             deallocate( send_status )
2406             nullify( send_status )
2407     
2408             deallocate( csendbuffer )
2409             nullify( csendbuffer )
2410     
2411             endif
2412     
2413     
2414     !       --------------------------
2415     !       wait for recvs to complete
2416     !       --------------------------
2417             if (nrecv.ge.1) then
2418     
2419             if (lidebug.ge.1) then
2420     
2421                call write_debug( name, &
2422                     'waiting for receives to complete, nrecv =  ', nrecv )
2423             endif
2424     
2425           if (use_waitany) then
2426             do ii=1,nrecv
2427                call MPI_WAITANY( nrecv, recvrequest, jindex, recv_status, ierror )
2428                call MPI_Check( 'sendrecv3_end_1c:MPI_WAITANY ', ierror )
2429     
2430                j1 = xrecv( jindex )
2431                j2 = xrecv( jindex + 1)-1
2432     
2433                if (lidebug.ge.2) then
2434                     call write_debug(name, 'jindex, j1,j2 ', jindex,j1,j2 )
2435                endif
2436     
2437                do jj=j1,j2
2438                  ijk = recvijk( jj )
2439     
2440                    do ic=1,clen
2441                     jpos = (jj-1)*clen + ic
2442                     XX(ijk)(ic:ic) = crecvbuffer(jpos)
2443                    enddo
2444     
2445     
2446                enddo
2447             enddo
2448           else
2449                call MPI_WAITALL( nrecv, recvrequest, recv_status, ierror )
2450                call MPI_Check( 'sendrecv3_end_1c:MPI_WAITALL recv ', ierror )
2451     
2452                j1 = xrecv(1)
2453                j2 = xrecv( nrecv +1)-1
2454                do jj=j1,j2
2455                    ijk = recvijk( jj )
2456     
2457                    do ic=1,clen
2458                     jpos = (jj-1)*clen + ic
2459                     XX(ijk)(ic:ic) = crecvbuffer(jpos)
2460                    enddo
2461     
2462     
2463     
2464                enddo
2465             endif
2466     
2467             deallocate( crecvbuffer )
2468             nullify( crecvbuffer )
2469     
2470             endif
2471     
2472             return
2473             end subroutine sendrecv3_end_1c
2474     
2475     
2476             subroutine sendrecv3_end_1i( XX, idebug )
2477     
2478               use functions
2479     
2480             implicit none
2481     
2482             integer, intent(inout), dimension(:) :: XX
2483             integer, intent(in), optional :: idebug
2484     
2485             interface
2486     
2487             subroutine MPI_WAITANY(count, array_of_requests, jindex, &
2488                                     status, ierror)
2489             use mpi, only: MPI_STATUS_SIZE
2490     
2491             integer count
2492             integer array_of_requests(*)
2493             integer jindex
2494             integer status(MPI_STATUS_SIZE)
2495             integer ierror
2496             end subroutine MPI_WAITANY
2497     
2498             subroutine MPI_WAITALL( count, array_of_requests,  &
2499                     array_of_status, ierror )
2500             use mpi, only: MPI_STATUS_SIZE
2501     
2502             integer count
2503             integer array_of_requests(*)
2504             integer array_of_status( MPI_STATUS_SIZE,*)
2505             integer ierror
2506             end subroutine MPI_WAITALL
2507     
2508             end interface
2509     
2510     !       ---------------
2511     !       local variables
2512     !       ---------------
2513             character(len=80), parameter :: name = 'sendrecv3_end_1i'
2514     
2515             logical, parameter :: use_waitany = .false.
2516     
2517             integer :: lidebug
2518             integer :: jj,ijk,  jindex, ii,j1,j2, ierror
2519     
2520             integer, dimension(MPI_STATUS_SIZE) :: recv_status
2521             integer, dimension(:,:), pointer :: send_status
2522     
2523     !       --------------------------
2524     !       wait for sends to complete
2525     !       --------------------------
2526     
2527             lidebug = 0
2528             if (present(idebug)) then
2529                lidebug = idebug
2530             endif
2531     
2532             if (nsend .ge.1) then
2533     
2534             if (lidebug.ge.1) then
2535     
2536                call write_debug(name, &
2537                     'waiting for sends to complete, nsend  = ', nsend )
2538             endif
2539     
2540             allocate( send_status(MPI_STATUS_SIZE,nsend))
2541     
2542             call MPI_WAITALL( nsend, sendrequest, send_status, ierror )
2543             call MPI_Check( 'sendrecv3_end_1i:MPI_WAITALL ', ierror )
2544     
2545             deallocate( send_status )
2546             nullify( send_status )
2547     
2548             deallocate( isendbuffer )
2549             nullify( isendbuffer )
2550     
2551             endif
2552     
2553     
2554     !       --------------------------
2555     !       wait for recvs to complete
2556     !       --------------------------
2557             if (nrecv.ge.1) then
2558     
2559             if (lidebug.ge.1) then
2560     
2561                call write_debug( name, &
2562                     'waiting for receives to complete, nrecv =  ', nrecv )
2563             endif
2564     
2565           if (use_waitany) then
2566             do ii=1,nrecv
2567                call MPI_WAITANY( nrecv, recvrequest, jindex, recv_status, ierror )
2568                call MPI_Check( 'sendrecv3_end_1i:MPI_WAITANY ', ierror )
2569     
2570                j1 = xrecv( jindex )
2571                j2 = xrecv( jindex + 1)-1
2572     
2573                if (lidebug.ge.2) then
2574                     call write_debug(name, 'jindex, j1,j2 ', jindex,j1,j2 )
2575                endif
2576     
2577                do jj=j1,j2
2578                  ijk = recvijk( jj )
2579                  XX(ijk) = irecvbuffer(jj)
2580                enddo
2581             enddo
2582           else
2583                call MPI_WAITALL( nrecv, recvrequest, recv_status, ierror )
2584                call MPI_Check( 'sendrecv3_end_1i:MPI_WAITALL recv ', ierror )
2585     
2586                j1 = xrecv(1)
2587                j2 = xrecv( nrecv +1)-1
2588                do jj=j1,j2
2589                    ijk = recvijk( jj )
2590                    XX(ijk) = irecvbuffer(jj)
2591                enddo
2592             endif
2593     
2594             deallocate( irecvbuffer )
2595             nullify( irecvbuffer )
2596     
2597             endif
2598     
2599             return
2600             end subroutine sendrecv3_end_1i
2601     
2602     
2603             subroutine send_recv3_1c( XX, ilayer, idebug )
2604             implicit none
2605     
2606             character(len=*),  dimension(:), intent(inout) :: XX
2607             integer, intent(in), optional :: ilayer,idebug
2608     
2609             integer :: lidebug, layer
2610     
2611             lidebug = 0
2612             if (present(idebug)) then
2613                lidebug = idebug
2614             endif
2615     
2616             layer = 1
2617             if (present(ilayer)) then
2618                layer = ilayer
2619             endif
2620     
2621             call sendrecv3_begin(XX,layer,lidebug)
2622             call sendrecv3_end( XX, lidebug )
2623     
2624             return
2625             end subroutine send_recv3_1c
2626     
2627             subroutine send_recv3_1d( XX, ilayer, idebug )
2628             implicit none
2629     
2630             double precision,  dimension(:), intent(inout) :: XX
2631             integer, intent(in), optional :: ilayer,idebug
2632     
2633             integer :: lidebug, layer
2634     
2635             lidebug = 0
2636             if (present(idebug)) then
2637                lidebug = idebug
2638             endif
2639     
2640             layer = 1
2641             if (present(ilayer)) then
2642                layer = ilayer
2643             endif
2644     
2645             call sendrecv3_begin(XX,layer,lidebug)
2646             call sendrecv3_end( XX, lidebug )
2647     
2648             return
2649             end subroutine send_recv3_1d
2650     
2651             subroutine send_recv3_2d( XX, ilayer, idebug )
2652             implicit none
2653     
2654             double precision,  dimension(:,:), intent(inout) :: XX
2655             integer, intent(in), optional :: ilayer,idebug
2656     
2657             integer :: lidebug, layer
2658             integer :: j
2659     
2660             lidebug = 0
2661             if (present(idebug)) then
2662                lidebug = idebug
2663             endif
2664     
2665             layer = 1
2666             if (present(ilayer)) then
2667                layer = ilayer
2668             endif
2669     
2670             do j=lbound(XX,2),ubound(XX,2)
2671               call sendrecv3_begin(XX(:,j),layer,lidebug)
2672               call sendrecv3_end( XX(:,j), lidebug )
2673             enddo
2674     
2675             return
2676             end subroutine send_recv3_2d
2677     
2678             subroutine send_recv3_3d( XX, ilayer, idebug )
2679             implicit none
2680     
2681             double precision,  dimension(:,:,:), intent(inout) :: XX
2682             integer, intent(in), optional :: ilayer,idebug
2683     
2684             integer :: lidebug, layer
2685             integer :: j,k
2686     
2687             lidebug = 0
2688             if (present(idebug)) then
2689                lidebug = idebug
2690             endif
2691     
2692             layer = 1
2693             if (present(ilayer)) then
2694                layer = ilayer
2695             endif
2696     
2697             do k=lbound(XX,3),ubound(XX,3)
2698             do j=lbound(XX,2),ubound(XX,2)
2699               call sendrecv3_begin(XX(:,j,k),layer,lidebug)
2700               call sendrecv3_end( XX(:,j,k), lidebug )
2701             enddo
2702             enddo
2703     
2704             return
2705             end subroutine send_recv3_3d
2706     
2707             subroutine send_recv3_1i( XX, ilayer, idebug )
2708             implicit none
2709     
2710             integer,  dimension(:), intent(inout) :: XX
2711             integer, intent(in), optional :: ilayer,idebug
2712     
2713             integer :: lidebug, layer
2714     
2715             lidebug = 0
2716             if (present(idebug)) then
2717                lidebug = idebug
2718             endif
2719     
2720             layer = 1
2721             if (present(ilayer)) then
2722                layer = ilayer
2723             endif
2724     
2725             call sendrecv3_begin(XX,layer,lidebug)
2726             call sendrecv3_end( XX, lidebug )
2727     
2728             return
2729             end subroutine send_recv3_1i
2730     
2731     
2732             end module sendrecv3
2733