All Classes Files Functions Variables Pages
module_mpi_assemble.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  use mpi
129 
130  implicit none
131 
132  public :: assemble_force_sync_comm
135 
136  contains
137 
138 !
139 !
142 !***********************************************************************************************************************************************************************************
144 !***********************************************************************************************************************************************************************************
145 
146  use mpi
147 
148  use mod_global_variables, only :&
149  tg_cpu_neighbor&
150  ,rg_gll_acceleration&
151  ,ig_ndof&
152  ,ig_ngll_total&
153  ,ig_myrank&
154  ,ig_mpi_buffer_sizemax&
155  ,ig_ncpu_neighbor
156 
157  implicit none
158 
159  real , dimension(IG_NDOF,ig_mpi_buffer_sizemax) :: dummyr
160  real , dimension(IG_NDOF,ig_mpi_buffer_sizemax,ig_ncpu_neighbor) :: dlacce
161 
162  integer, dimension(MPI_STATUS_SIZE) :: statut
163  integer :: igll
164  integer :: icon
165  integer :: ngll
166  integer :: ios
167  !
168  !
169  !******************************************
170  !fill buffer dlacce with gll to be send
171  !******************************************
172  do icon = 1,ig_ncpu_neighbor
173  do igll = 1,tg_cpu_neighbor(icon)%ngll
174 
175  dlacce(1,igll,icon) = rg_gll_acceleration(1,tg_cpu_neighbor(icon)%gll_send(igll))
176  dlacce(2,igll,icon) = rg_gll_acceleration(2,tg_cpu_neighbor(icon)%gll_send(igll))
177  dlacce(3,igll,icon) = rg_gll_acceleration(3,tg_cpu_neighbor(icon)%gll_send(igll))
178 
179  enddo
180  enddo
181  !
182  !
183  !**********************************************************
184  !synchrone comm: exchange buffers between cpu in contact
185  !**********************************************************
186  do icon = 1,ig_ncpu_neighbor
187 
188  ngll = tg_cpu_neighbor(icon)%ngll
189 
190  call mpi_sendrecv(dlacce(1,1,icon),ig_ndof*ngll,mpi_real,tg_cpu_neighbor(icon)%icpu,300,dummyr(1,1),ig_ndof*ngll,mpi_real,tg_cpu_neighbor(icon)%icpu,300,mpi_comm_world,statut,ios)
191 
192  do igll = 1,ngll
193 
194  rg_gll_acceleration(1,tg_cpu_neighbor(icon)%gll_recv(igll)) = rg_gll_acceleration(1,tg_cpu_neighbor(icon)%gll_recv(igll)) + dummyr(1,igll)
195  rg_gll_acceleration(2,tg_cpu_neighbor(icon)%gll_recv(igll)) = rg_gll_acceleration(2,tg_cpu_neighbor(icon)%gll_recv(igll)) + dummyr(2,igll)
196  rg_gll_acceleration(3,tg_cpu_neighbor(icon)%gll_recv(igll)) = rg_gll_acceleration(3,tg_cpu_neighbor(icon)%gll_recv(igll)) + dummyr(3,igll)
197 
198  enddo
199 
200  enddo
201 
202  return
203 !***********************************************************************************************************************************************************************************
204  end subroutine assemble_force_sync_comm
205 !***********************************************************************************************************************************************************************************
206 
207 !
208 !
210 !***********************************************************************************************************************************************************************************
212 !***********************************************************************************************************************************************************************************
213 
214  use mpi
215 
216  use mod_global_variables, only :&
217  tg_cpu_neighbor&
218  ,rg_gll_acceleration&
219  ,ig_ndof&
220  ,ig_ngll_total&
221  ,ig_ncpu&
222  ,ig_myrank&
223  ,ig_mpi_buffer_offset&
224  ,rg_mpi_buffer_send&
225  ,rg_mpi_buffer_recv&
226  ,ig_mpi_request_send&
227  ,ig_mpi_request_recv&
228  ,ig_mpi_buffer_sizemax&
229  ,ig_ncpu_neighbor
230 
231  implicit none
232 
233  real, dimension(IG_NDOF,ig_mpi_buffer_sizemax,ig_ncpu_neighbor) :: dlacce
234  integer :: igll
235  integer :: icon
236  integer :: ngll
237  integer :: ios
238  integer :: dof_offset
239  integer :: proc_offset
240  !
241  !
242  !*****************************************************************************************************************************************
243  !fill buffer dlacce with gll to be send
244  !*****************************************************************************************************************************************
245  do icon = 1,ig_ncpu_neighbor
246  do igll =1,tg_cpu_neighbor(icon)%ngll
247 
248  dlacce(1,igll,icon) = rg_gll_acceleration(1,tg_cpu_neighbor(icon)%gll_send(igll))
249  dlacce(2,igll,icon) = rg_gll_acceleration(2,tg_cpu_neighbor(icon)%gll_send(igll))
250  dlacce(3,igll,icon) = rg_gll_acceleration(3,tg_cpu_neighbor(icon)%gll_send(igll))
251 
252  enddo
253  enddo
254  !
255  !
256  !*****************************************************************************************************************************************
257  !wait until all send buffers are released
258  !*****************************************************************************************************************************************
259  if ( (ig_ncpu > 1) .and. (ig_mpi_request_send(1) /= 0) ) call mpi_waitall(ig_ncpu_neighbor, ig_mpi_request_send, mpi_statuses_ignore, ios)
260  !
261  !
262  !*****************************************************************************************************************************************
263  !asynchrone comm: exchange buffers between cpu in contact
264  !*****************************************************************************************************************************************
265  do icon = 1,ig_ncpu_neighbor
266 
267  ngll = tg_cpu_neighbor(icon)%ngll
268  proc_offset = ig_mpi_buffer_offset(icon)
269 
270  !
271  !fill buffer for cpu myrank
272  do igll = 1,ngll
273 
274  !beginning of x buffer = 0
275  dof_offset = 0
276  rg_mpi_buffer_send(proc_offset+dof_offset+igll) = dlacce(1,igll,icon)
277 
278  !beginning of y buffer = ngll
279  dof_offset = ngll
280  rg_mpi_buffer_send(proc_offset+dof_offset+igll) = dlacce(2,igll,icon)
281 
282  !beginning of z buffer = 2*ngll
283  dof_offset = 2*ngll
284  rg_mpi_buffer_send(proc_offset+dof_offset+igll) = dlacce(3,igll,icon)
285 
286  enddo
287 
288  call mpi_irecv(rg_mpi_buffer_recv(proc_offset+1),ig_ndof*ngll,mpi_real,tg_cpu_neighbor(icon)%icpu,300,mpi_comm_world,ig_mpi_request_recv(icon),ios)
289  call mpi_isend(rg_mpi_buffer_send(proc_offset+1),ig_ndof*ngll,mpi_real,tg_cpu_neighbor(icon)%icpu,300,mpi_comm_world,ig_mpi_request_send(icon),ios)
290 
291  enddo
292 
293 
294  return
295 !***********************************************************************************************************************************************************************************
296  end subroutine assemble_force_async_comm_init
297 !***********************************************************************************************************************************************************************************
298 
299 !
300 !
303 !***********************************************************************************************************************************************************************************
305 !***********************************************************************************************************************************************************************************
306 
307  use mpi
308 
309  use mod_global_variables, only :&
310  tg_cpu_neighbor&
311  ,rg_gll_acceleration&
312  ,ig_ndof&
313  ,ig_ngll_total&
314  ,ig_myrank&
315  ,ig_mpi_buffer_offset&
316  ,rg_mpi_buffer_send&
317  ,rg_mpi_buffer_recv&
318  ,ig_mpi_request_send&
319  ,ig_mpi_request_recv&
320  ,ig_mpi_buffer_sizemax&
321  ,ig_ncpu_neighbor
322 
323  implicit none
324 
325  integer :: igll
326  integer :: ngll
327  integer :: ios
328  integer :: dof_offset
329  integer :: proc_offset
330  integer :: nb_msg
331  integer :: r_index
332  integer, dimension(MPI_STATUS_SIZE) :: statut
333 
334  nb_msg = 0
335 
336  do while(nb_msg < ig_ncpu_neighbor)
337  !
338  !---->wait for any received message
339  call mpi_waitany(ig_ncpu_neighbor, ig_mpi_request_recv, r_index, statut, ios)
340 
341  nb_msg = nb_msg + 1
342  ngll = tg_cpu_neighbor(r_index)%ngll
343  proc_offset = ig_mpi_buffer_offset(r_index)
344  !
345  !---->sum contributions for this proc
346 
347  !beginning of x buffer = 0
348  dof_offset = 0
349  do igll = 1,ngll
350  rg_gll_acceleration(1,tg_cpu_neighbor(r_index)%gll_recv(igll)) = rg_gll_acceleration(1,tg_cpu_neighbor(r_index)%gll_recv(igll)) + rg_mpi_buffer_recv(proc_offset+dof_offset+igll)
351  enddo
352 
353  !beginning of y buffer = ngll
354  dof_offset = ngll
355  do igll = 1,ngll
356  rg_gll_acceleration(2,tg_cpu_neighbor(r_index)%gll_recv(igll)) = rg_gll_acceleration(2,tg_cpu_neighbor(r_index)%gll_recv(igll)) + rg_mpi_buffer_recv(proc_offset+dof_offset+igll)
357  enddo
358 
359  !beginning of z buffer = 2*ngll
360  dof_offset = 2*ngll
361  do igll = 1,ngll
362  rg_gll_acceleration(3,tg_cpu_neighbor(r_index)%gll_recv(igll)) = rg_gll_acceleration(3,tg_cpu_neighbor(r_index)%gll_recv(igll)) + rg_mpi_buffer_recv(proc_offset+dof_offset+igll)
363  enddo
364 
365  enddo
366 
367  return
368 !***********************************************************************************************************************************************************************************
369  end subroutine assemble_force_async_comm_end
370 !***********************************************************************************************************************************************************************************
371 
372 end module mod_mpi_assemble
subroutine, public assemble_force_async_comm_end()
subroutine to finalize forces assembly between connected cpus by MPI asynchrone communications.
subroutine, public assemble_force_sync_comm()
subroutine to assemble forces between connected cpus by MPI synchrone communications.
This module defines all global variables of EFISPEC3D. Scalar variables are initialized directly in t...
This module contains subroutines to assemble forces between cpu myrank and its neighbor cpus by synch...
subroutine, public assemble_force_async_comm_init()
subroutine to initialize forces assembly between connected cpus by MPI asynchrone communications...