All Classes Files Functions Variables Pages
module_init_mpibuffer.f90
Go to the documentation of this file.
1 !=====================================================================================================================================
2 ! EFISPEC3D !
3 ! (Elements FInis SPECtraux 3D) !
4 ! !
5 ! http://efispec.free.fr !
6 ! !
7 ! !
8 ! This file is part of EFISPEC3D !
9 ! Please refer to http://efispec.free.fr if you use it or part of it !
10 ! !
11 ! !
12 !1 ---> French License: CeCILL V2 !
13 ! !
14 ! Copyright BRGM 2009 contributeurs : Florent DE MARTIN !
15 ! David MICHEA !
16 ! Philippe THIERRY !
17 ! !
18 ! Contact: f.demartin at brgm.fr !
19 ! !
20 ! Ce logiciel est un programme informatique servant a resoudre l'equation du !
21 ! mouvement en trois dimensions via une methode des elements finis spectraux. !
22 ! !
23 ! Ce logiciel est regi par la licence CeCILL soumise au droit francais et !
24 ! respectant les principes de diffusion des logiciels libres. Vous pouvez !
25 ! utiliser, modifier et/ou redistribuer ce programme sous les conditions de la !
26 ! licence CeCILL telle que diffusee par le CEA, le CNRS et l'INRIA sur le site !
27 ! "http://www.cecill.info". !
28 ! !
29 ! En contrepartie de l'accessibilite au code source et des droits de copie, de !
30 ! modification et de redistribution accordes par cette licence, il n'est offert !
31 ! aux utilisateurs qu'une garantie limitee. Pour les memes raisons, seule une !
32 ! responsabilite restreinte pese sur l'auteur du programme, le titulaire des !
33 ! droits patrimoniaux et les concedants successifs. !
34 ! !
35 ! A cet egard l'attention de l'utilisateur est attiree sur les risques associes !
36 ! au chargement, a l'utilisation, a la modification et/ou au developpement et a !
37 ! la reproduction du logiciel par l'utilisateur etant donne sa specificite de !
38 ! logiciel libre, qui peut le rendre complexe a manipuler et qui le reserve donc !
39 ! a des developpeurs et des professionnels avertis possedant des connaissances !
40 ! informatiques approfondies. Les utilisateurs sont donc invites a charger et !
41 ! tester l'adequation du logiciel a leurs besoins dans des conditions permettant !
42 ! d'assurer la securite de leurs systemes et ou de leurs donnees et, plus !
43 ! generalement, a l'utiliser et l'exploiter dans les memes conditions de !
44 ! securite. !
45 ! !
46 ! Le fait que vous puissiez acceder a cet en-tete signifie que vous avez pris !
47 ! connaissance de la licence CeCILL et que vous en avez accepte les termes. !
48 ! !
49 ! !
50 !2 ---> International license: GNU GPL V3 !
51 ! !
52 ! EFISPEC3D is a computer program that solves the three-dimensional equations of !
53 ! motion using a finite spectral-element method. !
54 ! !
55 ! Copyright (C) 2009 Florent DE MARTIN !
56 ! !
57 ! Contact: f.demartin at brgm.fr !
58 ! !
59 ! This program is free software: you can redistribute it and/or modify it under !
60 ! the terms of the GNU General Public License as published by the Free Software !
61 ! Foundation, either version 3 of the License, or (at your option) any later !
62 ! version. !
63 ! !
64 ! This program is distributed in the hope that it will be useful, but WITHOUT ANY !
65 ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A !
66 ! PARTICULAR PURPOSE. See the GNU General Public License for more details. !
67 ! !
68 ! You should have received a copy of the GNU General Public License along with !
69 ! this program. If not, see http://www.gnu.org/licenses/. !
70 ! !
71 ! !
72 !3 ---> Third party libraries !
73 ! !
74 ! EFISPEC3D uses the following of third party libraries: !
75 ! !
76 ! --> METIS 5.1.0 !
77 ! see http://glaros.dtc.umn.edu/gkhome/metis/metis/overview !
78 ! !
79 ! --> Lib_VTK_IO !
80 ! see S. Zaghi's website: https://github.com/szaghi/Lib_VTK_IO !
81 ! !
82 ! --> INTERP_LINEAR !
83 ! see J. Burkardt website: http://people.sc.fsu.edu/~jburkardt/ !
84 ! !
85 ! --> SAC !
86 ! http://ds.iris.edu/files/sac-manual/ !
87 ! !
88 ! --> EXODUS II !
89 ! http://sourceforge.net/projects/exodusii/ !
90 ! !
91 ! --> NETCDF !
92 ! http://www.unidata.ucar.edu/software/netcdf/ !
93 ! !
94 ! --> HDF5 !
95 ! http://www.hdfgroup.org/HDF5/ !
96 ! !
97 ! Some of these libraries are located in directory lib and pre-compiled !
98 ! with intel compiler for x86_64 platform. !
99 ! !
100 ! !
101 !4 ---> Related Articles !
102 ! !
103 ! De Martin, F., Matsushima, M., Kawase, H. (BSSA, 2013) !
104 ! Impact of geometric effects on near-surface Green's functions !
105 ! doi:10.1785/0120130039 !
106 ! !
107 ! Aochi, H., Ducellier, A., Dupros, F., Delatre, M., Ulrich, T., De Martin, F., !
108 ! Yoshimi, M., (Pure Appl. Geophys. 2013) !
109 ! Finite Difference Simulations of Seismic Wave Propagation for the 2007 Mw 6.6 !
110 ! Niigata-ken Chuetsu-Oki Earthquake: Validity of Models and Reliable !
111 ! Input Ground Motion in the Near-Field !
112 ! doi:10.1007/s00024-011-0429-5 !
113 ! !
114 ! De Martin, F. (BSSA, 2011) !
115 ! Verification of a Spectral-Element Method Code for the Southern California !
116 ! Earthquake Center LOH.3 Viscoelastic Case !
117 ! doi:10.1785/0120100305 !
118 ! !
119 !=====================================================================================================================================
120 
123 
127 
128  implicit none
129 
130  private
131 
132  public :: init_mpi_buffers
133  private :: create_gll_buffer_recv
134  private :: remove_duplicate
135 
136  contains
137 
138 !
139 !
145 !***********************************************************************************************************************************************************************************
146  subroutine init_mpi_buffers()
147 !***********************************************************************************************************************************************************************************
148 
149  use mpi
150 
151  use mod_global_variables, only :&
152  ig_ngll&
153  ,ig_myrank&
154  ,ig_ncpu_neighbor&
155  ,ig_nhexa_outer&
156  ,ig_hexa_gll_glonum&
157  ,tg_cpu_neighbor&
158  ,ig_cpu_neighbor_info&
159  ,ig_nneighbor_all_kind&
160  ,rg_gll_coordinate&
161  ,ig_mpi_buffer_offset&
162  ,rg_mpi_buffer_send&
163  ,rg_mpi_buffer_recv&
164  ,ig_ndof&
165  ,ig_mpi_request_send&
166  ,ig_mpi_request_recv&
167  ,error_stop&
168  ,cg_prefix&
169  ,cg_myrank&
170  ,get_newunit&
171  ,lg_async_mpi_comm&
172  ,lg_output_debug_file&
173  ,ig_lst_unit&
174  ,ig_mpi_buffer_sizemax
175 
176  use mod_init_memory
177 
178  implicit none
179 
180  integer, parameter :: nface = 6
181  integer, parameter :: nedge = 12
182  integer, parameter :: nnode = 8
183 
184  real , dimension(:,:), allocatable :: gll_coord_send
185  real , dimension(:,:), allocatable :: gll_coord_recv
186 
187  integer, dimension(MPI_STATUS_SIZE) :: statut
188  integer :: i
189  integer :: j
190  integer :: ios
191  integer :: myunit_debug
192  integer :: icpu_neighbor
193  integer :: icpu
194  integer :: isurf
195  integer :: ngll_duplicate
196  integer :: ngll_send
197  integer :: ngll_recv
198  integer :: num_gll
199  integer :: surf_num
200  integer :: elt_num
201  integer :: ngll_unique
202  integer :: mpi_buffer_sizemax
203  integer :: buffer_sizemax
204  integer, dimension(:), allocatable :: buffer_gll_duplicate
205 
206  character(len= 6) :: cl_rank
207  character(len=255) :: info
208 
209  if (ig_myrank == 0) then
210  write(ig_lst_unit,'(" ",/,a)') "creating mpi buffers between cpu myrank and its neighbors..."
211  call flush(ig_lst_unit)
212  endif
213 
214  buffer_sizemax = (6*ig_ngll*ig_ngll + 12*ig_ngll + 8) * ig_nhexa_outer !DAVID: worst case, may be reduced
215 
216  ios = init_array_int(buffer_gll_duplicate,buffer_sizemax,"buffer_gll_duplicate")
217 
218  if (lg_async_mpi_comm) then
219 
220  ios = init_array_int(ig_mpi_request_send,ig_ncpu_neighbor,"ig_mpi_request_send")
221  ios = init_array_int(ig_mpi_request_recv,ig_ncpu_neighbor,"ig_mpi_request_recv")
222  ios = init_array_int(ig_mpi_buffer_offset,ig_ncpu_neighbor+1,"ig_mpi_buffer_offset")
223 
224  endif
225 
226  mpi_buffer_sizemax = 0
227 
228  do icpu_neighbor = 1,ig_ncpu_neighbor
229 
230  icpu = tg_cpu_neighbor(icpu_neighbor)%icpu
231  buffer_gll_duplicate(:) = 0
232  ngll_duplicate = 0
233 
234  do isurf = 1,ig_nneighbor_all_kind
235 
236  if (ig_cpu_neighbor_info(3,isurf) == icpu) then
237 
238  surf_num = ig_cpu_neighbor_info(2,isurf)
239  elt_num = ig_cpu_neighbor_info(1,isurf)
240 !
241 !----------------->put gll points from contact face in the buffer
242  if (surf_num <= nface) then
243  do i = 1,ig_ngll
244  do j = 1,ig_ngll
245  select case(surf_num)
246  case(1)
247  num_gll = ig_hexa_gll_glonum(j,i,1,elt_num)
248  case(2)
249  num_gll = ig_hexa_gll_glonum(j,1,i,elt_num)
250  case(3)
251  num_gll = ig_hexa_gll_glonum(ig_ngll,j,i,elt_num)
252  case(4)
253  num_gll = ig_hexa_gll_glonum(j,ig_ngll,i,elt_num)
254  case(5)
255  num_gll = ig_hexa_gll_glonum(1,j,i,elt_num)
256  case(6)
257  num_gll = ig_hexa_gll_glonum(j,i,ig_ngll,elt_num)
258  end select
259 
260  ngll_duplicate = ngll_duplicate + 1
261 
262  if (ngll_duplicate > buffer_sizemax) then
263  write(info,'(a)') "error in subroutine init_mpi_buffers: face in contact : size of buffer_gll_duplicate too small"
264  call error_stop(info)
265  endif
266 
267  buffer_gll_duplicate(ngll_duplicate) = num_gll
268 
269  enddo
270  enddo
271 !
272 !----------------->put gll points from contact edge in the buffer
273  else if(surf_num <= nface+nedge) then
274  do i = 1,ig_ngll
275  select case(surf_num - nface)
276  case(1)
277  num_gll = ig_hexa_gll_glonum(i,1,1,elt_num)
278  case(2)
279  num_gll = ig_hexa_gll_glonum(ig_ngll,i,1,elt_num)
280  case(3)
281  num_gll = ig_hexa_gll_glonum(i,ig_ngll,1,elt_num)
282  case(4)
283  num_gll = ig_hexa_gll_glonum(1,i,1,elt_num)
284  case(5)
285  num_gll = ig_hexa_gll_glonum(i,1,ig_ngll,elt_num)
286  case(6)
287  num_gll = ig_hexa_gll_glonum(ig_ngll,i,ig_ngll,elt_num)
288  case(7)
289  num_gll = ig_hexa_gll_glonum(i,ig_ngll,ig_ngll,elt_num)
290  case(8)
291  num_gll = ig_hexa_gll_glonum(1,i,ig_ngll,elt_num)
292  case(9)
293  num_gll = ig_hexa_gll_glonum(1,1,i,elt_num)
294  case(10)
295  num_gll = ig_hexa_gll_glonum(ig_ngll,1,i,elt_num)
296  case(11)
297  num_gll = ig_hexa_gll_glonum(ig_ngll,ig_ngll,i,elt_num)
298  case(12)
299  num_gll = ig_hexa_gll_glonum(1,ig_ngll,i,elt_num)
300  end select
301 
302  ngll_duplicate = ngll_duplicate + 1
303 
304  if (ngll_duplicate > buffer_sizemax) then
305  write(info,'(a)') "error in subroutine init_mpi_buffers: edge in contact : size of buffer_gll_duplicate too small"
306  call error_stop(info)
307  endif
308 
309  buffer_gll_duplicate(ngll_duplicate) = num_gll
310 
311  enddo
312 !
313 !----------------->put gll points from contact corner nodes in the buffer
314  else
315  select case(surf_num - (nface+nedge))
316  case(1)
317  num_gll = ig_hexa_gll_glonum(1,1,1,elt_num)
318  case(2)
319  num_gll = ig_hexa_gll_glonum(ig_ngll,1,1,elt_num)
320  case(3)
321  num_gll = ig_hexa_gll_glonum(ig_ngll,ig_ngll,1,elt_num)
322  case(4)
323  num_gll = ig_hexa_gll_glonum(1,ig_ngll,1,elt_num)
324  case(5)
325  num_gll = ig_hexa_gll_glonum(1,1,ig_ngll,elt_num)
326  case(6)
327  num_gll = ig_hexa_gll_glonum(ig_ngll,1,ig_ngll,elt_num)
328  case(7)
329  num_gll = ig_hexa_gll_glonum(ig_ngll,ig_ngll,ig_ngll,elt_num)
330  case(8)
331  num_gll = ig_hexa_gll_glonum(1,ig_ngll,ig_ngll,elt_num)
332  end select
333 
334  ngll_duplicate = ngll_duplicate + 1
335 
336  if (ngll_duplicate > buffer_sizemax) then
337  write(info,'(a)') "error in subroutine init_mpi_buffers: node in contact : size of buffer_gll_duplicate too small"
338  call error_stop(info)
339  endif
340 
341  buffer_gll_duplicate(ngll_duplicate) = num_gll
342 
343  endif
344  endif
345  enddo
346 !
347 !--------->remove duplicates
348  call remove_duplicate(buffer_gll_duplicate,ngll_duplicate,tg_cpu_neighbor(icpu_neighbor)%gll_send,ngll_unique)
349 
350  if (mpi_buffer_sizemax < ngll_unique) mpi_buffer_sizemax = ngll_unique
351  tg_cpu_neighbor(icpu_neighbor)%ngll = ngll_unique
352 !
353 !--------->ig_mpi_buffer_offset(x) -> give offset in rg_mpi_buffer for xth connected proc (start from 0 !!!)
354  if (lg_async_mpi_comm) then
355  ig_mpi_buffer_offset(icpu_neighbor+1) = ig_mpi_buffer_offset(icpu_neighbor) + tg_cpu_neighbor(icpu_neighbor)%ngll*ig_ndof
356  endif
357 
358  if (lg_output_debug_file) then
359  open(unit=get_newunit(myunit_debug),file="debug."//trim(cg_prefix)//".mpi.buffer.cpu."//trim(cg_myrank))
360  write(unit=myunit_debug,fmt='(3(a,i6),a)') "from proc ",ig_myrank," to proc ",icpu," : ",ngll_unique," gll points in MPI buffer"
361  do i = 1,ngll_unique
362  write(unit=myunit_debug,fmt='(i10,3(e14.7,1x))') tg_cpu_neighbor(icpu_neighbor)%gll_send(i)&
363  ,rg_gll_coordinate(1,tg_cpu_neighbor(icpu_neighbor)%gll_send(i))&
364  ,rg_gll_coordinate(2,tg_cpu_neighbor(icpu_neighbor)%gll_send(i))&
365  ,rg_gll_coordinate(3,tg_cpu_neighbor(icpu_neighbor)%gll_send(i))
366  enddo
367  close(myunit_debug)
368  endif
369 
370  enddo
371 
372 !
373 !----->check if cpus connected together send and receive the same number of gll
374  do icpu_neighbor = 1,ig_ncpu_neighbor
375 
376  icpu = tg_cpu_neighbor(icpu_neighbor)%icpu
377  ngll_send = tg_cpu_neighbor(icpu_neighbor)%ngll
378 
379  call mpi_sendrecv(ngll_send,1,mpi_integer,icpu,89,ngll_recv,1,mpi_integer,icpu,89,mpi_comm_world,statut,ios)
380 
381  if (ngll_send /= ngll_recv) then
382  write(cl_rank,'(i6.6)') icpu
383  write(info,'(a)') "error in subroutine init_mpi_buffers: different number ngll_send and ngll_recv betwen cpu"//trim(cg_myrank)//" and "//trim(cl_rank)
384  call error_stop(info)
385  endif
386 
387  enddo
388 
389 !
390 !----->set maximum size of MPI buffers between cpu myrank and its neighbors
391  ig_mpi_buffer_sizemax = mpi_buffer_sizemax
392 
393 !
394 !----->initialize memory of GLL nodes coordinates
395 
396  ios = init_array_real(gll_coord_send,ig_ncpu_neighbor,ig_ndof*mpi_buffer_sizemax,"gll_coord_send")
397 
398  ios = init_array_real(gll_coord_recv,ig_ncpu_neighbor,ig_ndof*mpi_buffer_sizemax,"gll_coord_recv")
399 
400 !
401 !----->send/recv 2D buffer's glls coordinates for linked procs
402  do icpu_neighbor = 1,ig_ncpu_neighbor
403 
404  icpu = tg_cpu_neighbor(icpu_neighbor)%icpu
405 
406  !
407  !fill coords to send
408  do i = 1, tg_cpu_neighbor(icpu_neighbor)%ngll
409  gll_coord_send((i-1)*ig_ndof+1,icpu_neighbor) = rg_gll_coordinate(1,tg_cpu_neighbor(icpu_neighbor)%gll_send(i))
410  gll_coord_send((i-1)*ig_ndof+2,icpu_neighbor) = rg_gll_coordinate(2,tg_cpu_neighbor(icpu_neighbor)%gll_send(i))
411  gll_coord_send((i-1)*ig_ndof+3,icpu_neighbor) = rg_gll_coordinate(3,tg_cpu_neighbor(icpu_neighbor)%gll_send(i))
412  enddo
413 
414  call mpi_sendrecv(gll_coord_send(1,icpu_neighbor),ig_ndof*tg_cpu_neighbor(icpu_neighbor)%ngll,mpi_real,icpu,90,gll_coord_recv(1,icpu_neighbor),ig_ndof*tg_cpu_neighbor(icpu_neighbor)%ngll,mpi_real,icpu,90,mpi_comm_world,statut,ios)
415 
416  if (ios /= 0) then
417  write(info,'(a)') "error in subroutine init_mpi_buffers while sending gll coordinates"
418  call error_stop(info)
419  endif
420 
421  enddo
422 
423 !
424 !----->build eqres
425  do icpu_neighbor = 1,ig_ncpu_neighbor
426  call create_gll_buffer_recv(gll_coord_send(1,icpu_neighbor),gll_coord_recv(1,icpu_neighbor),tg_cpu_neighbor(icpu_neighbor)%ngll,tg_cpu_neighbor(icpu_neighbor)%gll_recv,tg_cpu_neighbor(icpu_neighbor)%gll_send,mpi_buffer_sizemax)
427  enddo
428 
429 !
430 !----->deallocate arrays
431  deallocate(gll_coord_send)
432  deallocate(gll_coord_recv)
433  deallocate(ig_cpu_neighbor_info)
434 
435 !
436 !----->allocate array for mpi buffers
437  if (lg_async_mpi_comm) then
438 
439  ios = init_array_real(rg_mpi_buffer_send,ig_mpi_buffer_offset(ig_ncpu_neighbor+1),"rg_mpi_buffer_send")
440 
441  ios = init_array_real(rg_mpi_buffer_recv,ig_mpi_buffer_offset(ig_ncpu_neighbor+1),"rg_mpi_buffer_recv")
442 
443  endif
444 
445  if (ig_myrank == 0) then
446  write(ig_lst_unit,'(a)') "done"
447  call flush(ig_lst_unit)
448  endif
449 
450 !***********************************************************************************************************************************************************************************
451  end subroutine init_mpi_buffers
452 !***********************************************************************************************************************************************************************************
453 
454 !
455 !
466 !***********************************************************************************************************************************************************************************
467  subroutine create_gll_buffer_recv(gll_coord_send, gll_coord_recv, ngll, gll_recv, gll_send, max_size)
468 !***********************************************************************************************************************************************************************************
469 
470  use mpi
471 
472  use mod_global_variables, only :&
473  ig_ncpu_neighbor&
474  ,ig_myrank&
475  ,ig_ndof&
476  ,error_stop
477 
478  use mod_init_memory
479 
480  implicit none
481 
482  integer, intent(in) :: ngll
483  integer, intent(in) :: max_size
484  integer, intent(in) , dimension(ngll) :: gll_send
485  integer, intent(out), allocatable, dimension(:) :: gll_recv
486 
487  real , intent(in), dimension(IG_NDOF*max_size) :: gll_coord_send
488  real , intent(in), dimension(IG_NDOF*max_size) :: gll_coord_recv
489 
490  integer :: i
491  integer :: j
492  integer :: myindex
493  integer :: ios
494  integer :: imin
495 
496  real :: xs
497  real :: ys
498  real :: zs
499  real :: xr
500  real :: yr
501  real :: zr
502  real :: dist
503  real :: mindist
504 
505 
506  ios = init_array_int(gll_recv,ngll,"gll_recv")
507 
508 ext: do i = 1,ngll
509 
510  mindist = huge(dist)
511  imin = 0
512  xr = gll_coord_recv((i-1)*ig_ndof+1)
513  yr = gll_coord_recv((i-1)*ig_ndof+2)
514  zr = gll_coord_recv((i-1)*ig_ndof+3)
515 
516 int: do j = 1,ngll
517 
518  xs = gll_coord_send((j-1)*ig_ndof+1)
519  ys = gll_coord_send((j-1)*ig_ndof+2)
520  zs = gll_coord_send((j-1)*ig_ndof+3)
521 
522  dist = sqrt((xs-xr)**2 + (ys-yr)**2 + (zs-zr)**2)
523 
524  if ( dist == 0.0 ) then
525  gll_recv(i) = gll_send(j)
526  cycle ext
527  elseif (dist < mindist) then
528  mindist = dist
529  myindex = j
530  endif
531 
532  enddo int
533 
534  gll_recv(i) = gll_send(myindex)
535 
536  enddo ext
537 
538  return
539 
540 !***********************************************************************************************************************************************************************************
541  end subroutine create_gll_buffer_recv
542 !***********************************************************************************************************************************************************************************
543 
544 !
545 !
552 !***********************************************************************************************************************************************************************************
553  subroutine remove_duplicate(x1,n,x3,m)
554 !***********************************************************************************************************************************************************************************
555 
556  implicit none
557 
558  integer, intent( in) :: n
559  integer, intent( in) , dimension(n) :: x1
560  integer, intent(out), allocatable, dimension(:) :: x3
561  integer, intent(out) :: m
562 
563  integer , dimension(n) :: x2 !array without duplicate (size n)
564  integer :: i
565  integer :: j
566 
567  m = 1
568  x2(1) = x1(1)
569 
570 outer: do i=2,n
571  do j=1,m
572  if (x2(j) == x1(i)) then
573 !
574 !----------------->Found a match so start looking again
575  cycle outer
576  endif
577  enddo
578 !
579 !--------->No match found so add it to array x2
580  m = m + 1
581  x2(m) = x1(i)
582  enddo outer
583 
584  allocate(x3(m))
585  x3(1:m) = x2(1:m)
586 
587  return
588 
589 !***********************************************************************************************************************************************************************************
590  end subroutine remove_duplicate
591 !***********************************************************************************************************************************************************************************
592 
593 end module mod_init_mpibuffer
This module contains subroutines to allocate arrays and to compute an approximation of the total RAM ...
subroutine, public init_mpi_buffers()
This subroutine searches for common GLL nodes between cpu myrank and its neighbor cpus...
This module defines all global variables of EFISPEC3D. Scalar variables are initialized directly in t...
Interface init_array_real to redirect allocation to n-rank arrays.
subroutine, private remove_duplicate(x1, n, x3, m)
This subroutine removes duplicated GLLs between cpu myrank and its connected cpus.
subroutine, private create_gll_buffer_recv(gll_coord_send, gll_coord_recv, ngll, gll_recv, gll_send, max_size)
This subroutine creates GLLs buffers received by cpu myrank from other cpus. Identical GLLs between c...
This module contains subroutines to initialize MPI buffers between cpu myrank and its neighbor cpus...
Interface init_array_int to redirect allocation to n-rank arrays.