All Classes Files Functions Variables Pages
efispec.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 
129 
131 program efispec3d
132 
133  use mpi
134 
135  use mod_global_variables, only :&
136  ig_ndt&
137  ,rg_simu_current_time&
138  ,ig_idt&
139  ,rg_dt&
140  ,rg_dt2&
141  ,ig_lst_unit&
142  ,lg_output_cputime&
143  ,lg_snapshot_vtk&
144  ,ig_lagrange_order&
145  ,lg_snapshot&
146  ,lg_snapshot_volume&
147  ,ig_nreceiver_snapshot&
148  ,ig_receiver_saving_incr&
149  ,ig_ncpu&
150  ,ig_nhexa&
151  ,ig_myrank&
152  ,cg_myrank&
153  ,cg_cpu_name&
154  ,ig_cpu_name_len&
155  ,cg_ncpu&
156  ,cg_prefix&
157  ,get_newunit&
158  ,error_stop&
159  ,lg_boundary_absorption&
160  ,ig_medium_type&
161  ,ig_nhexa_outer&
162  ,lg_async_mpi_comm
163 
164  use mod_init_efi
165 
166  use mod_init_mesh
167 
169 
170  use mod_solver
171 
172  use mod_init_medium
173 
175 
176  use mod_snapshot_surface , only :&
181 
182  use mod_snapshot_volume , only:&
186 
187  use mod_init_memory
188 
189  use mod_mpi_assemble
190 
191  use mod_receiver , only :&
195 
196  use mod_source , only :&
200 
201  implicit none
202 
203  double precision, allocatable, dimension(:) :: buffer_double
204  double precision, allocatable, dimension(:) :: buffer_double_all_cpu
205 
206  double precision :: start_time
207  double precision :: end_time
208 
209  double precision :: dltim1
210  double precision :: dltim2
211 
212  double precision :: time_init_mesh
213  double precision :: time_init_gll
214  double precision :: time_init_mpi_buffers
215  double precision :: time_init_medium
216  double precision :: time_init_time_step
217  double precision :: time_init_source
218  double precision :: time_init_receiver
219  double precision :: time_init_jacobian
220  double precision :: time_init_mass_matrix
221  double precision :: time_compute_dcsource
222  double precision :: time_init_snapshot
223  double precision :: time_memory_consumption
224 
225  double precision :: t_tloop
226  double precision :: t_newmark
227  double precision :: t_forext
228  double precision :: t_forabs
229  double precision :: t_forint_inner
230  double precision :: t_forint_outer
231  double precision :: t_resol
232  double precision :: t_snafsu
233  double precision :: t_outavd
234  double precision :: t_linkfo
235  double precision :: t_link_init
236 
237  integer, parameter :: ntime_init=12
238  integer, allocatable, dimension(:) :: buffer_integer
239  integer :: ios
240  integer :: icpu
241  integer :: itime
242  integer :: unit_time
243 
244  character(len= 95) :: myfmt
245  character(len=255) :: info
246 
247 !
248 !
249 !*************************************************************************************************************************
250 !*************************************************************************************************************************
251 !phase 1: initialization
252 !*************************************************************************************************************************
253 !*************************************************************************************************************************
254 
255 
256  !
257  !
258  !**********************************************************************************************************************
259  call init_mpi()
260  start_time = mpi_wtime()
261  !**********************************************************************************************************************
262 
263 
264  !
265  !
266  !**********************************************************************************************************************
267  call init_input_variables()
268  !**********************************************************************************************************************
269 
270 
271  !
272  !
273  !**********************************************************************************************************************
274  time_init_mesh = mpi_wtime()
275 
276  call init_mesh()
277 
278  time_init_mesh = mpi_wtime() - time_init_mesh
279  !**********************************************************************************************************************
280 
281 
282  !
283  !
284  !**********************************************************************************************************************
285  time_init_gll = mpi_wtime()
286 
287  call init_gll_nodes()
289 
290  time_init_gll = mpi_wtime() - time_init_gll
291  !**********************************************************************************************************************
292 
293 
294  !
295  !
296  !**********************************************************************************************************************
297  time_init_mpi_buffers = mpi_wtime()
298 
299  call init_mpi_buffers()
300 
301  time_init_mpi_buffers = mpi_wtime() - time_init_mpi_buffers
302  !**********************************************************************************************************************
303 
304 
305  !
306  !
307  !**********************************************************************************************************************
308  time_init_medium = mpi_wtime()
309 
310  if (ig_medium_type == 0) then
311 
312  if (ig_myrank == 0) write(ig_lst_unit,'(" ",/,a)') "generating medium with constant mechanical properties per hexa"
313  call init_hexa_medium()
314 
315  else
316 
317  write(info,'(a)') "error in efispec: illegal medium type"
318  call error_stop(info)
319 
320  endif
321 
322  call init_quadp_medium()
323 
324  time_init_medium = mpi_wtime() - time_init_medium
325  !**********************************************************************************************************************
326 
327  !
328  !
329  !**********************************************************************************************************************
330  time_init_time_step = mpi_wtime()
331 
332  call init_time_step()
333 
334  time_init_time_step = mpi_wtime() - time_init_time_step
335  !**********************************************************************************************************************
336 
337 
338  !
339  !
340  !**********************************************************************************************************************
341  time_init_source = mpi_wtime()
342 
345 
346  time_init_source = mpi_wtime() - time_init_source
347  !**********************************************************************************************************************
348 
349 
350  !
351  !
352  !**********************************************************************************************************************
353  time_init_receiver = mpi_wtime()
354 
355  call init_hexa_receiver()
356  call init_quad_receiver()
357 
358  time_init_receiver = mpi_wtime() - time_init_receiver
359  !**********************************************************************************************************************
360 
361 
362  !
363  !
364  !**********************************************************************************************************************
365  time_init_jacobian = mpi_wtime()
366 
369 
370  time_init_jacobian = mpi_wtime() - time_init_jacobian
371  !**********************************************************************************************************************
372 
373 
374  !
375  !
376  !**********************************************************************************************************************
377  time_init_mass_matrix = mpi_wtime()
378 
379  call init_mass_matrix()
380 
381  time_init_mass_matrix = mpi_wtime() - time_init_mass_matrix
382  !**********************************************************************************************************************
383 
384 
385  !
386  !
387  !**********************************************************************************************************************
388  time_compute_dcsource = mpi_wtime()
389 
391 
392  time_compute_dcsource = mpi_wtime() - time_compute_dcsource
393  !**********************************************************************************************************************
394 
395 
396  !
397  !
398  !**********************************************************************************************************************
399  time_init_snapshot = mpi_wtime()
400 
401  if (lg_snapshot) call init_snapshot_surface()
402 
403  if (lg_snapshot_volume) call init_snapshot_volume()
404 
405  time_init_snapshot = mpi_wtime() - time_init_snapshot
406  !**********************************************************************************************************************
407 
408 
409  !
410  !
411  !**********************************************************************************************************************
412  time_memory_consumption = mpi_wtime()
413 
414  call memory_consumption()
415 
416  time_memory_consumption = mpi_wtime() - time_memory_consumption
417  !**********************************************************************************************************************
418 
419 
420  end_time = (mpi_wtime() - start_time)
421 
422  call mpi_reduce(end_time,dltim1,1,mpi_double_precision,mpi_max,0,mpi_comm_world,ios)
423 
424 
425  !
426  !
427  !**********************************************************************************************************************
428  !summarize elapsed time in subroutines for initialization
429  !**********************************************************************************************************************
430  allocate(buffer_double(ntime_init))
431 
432  allocate(buffer_double_all_cpu(ntime_init*ig_ncpu))
433 
434  buffer_double( 1) = time_init_mesh
435  buffer_double( 2) = time_init_gll
436  buffer_double( 3) = time_init_mpi_buffers
437  buffer_double( 4) = time_init_medium
438  buffer_double( 5) = time_init_source
439  buffer_double( 6) = time_init_receiver
440  buffer_double( 7) = time_init_jacobian
441  buffer_double( 8) = time_init_mass_matrix
442  buffer_double( 9) = time_compute_dcsource
443  buffer_double(10) = time_init_snapshot
444  buffer_double(11) = time_memory_consumption
445  buffer_double(12) = time_init_time_step
446 
447  call mpi_gather(buffer_double,ntime_init,mpi_double_precision,buffer_double_all_cpu,ntime_init,mpi_double_precision,0,mpi_comm_world,ios)
448 
449  if (ig_myrank == 0) then
450 
451  write(myfmt,'(i)') ntime_init
452  myfmt = "(a,i6,1x,"//trim(adjustl(myfmt))//"(e14.7,1x))"
453 
454  write(unit=ig_lst_unit,fmt='(" ",/,a)') "elapsed time for initialization"
455  write(unit=ig_lst_unit,fmt='( a)') " init_mesh init_gll_nodes init_mpi_buff init_medium init_source init_receiver init_jaco init_mass init_dcs init_snapshot init_memory init_time_step"
456 
457  do icpu = 1,ig_ncpu
458  write(unit=ig_lst_unit,fmt=trim(myfmt)) "cpu ",icpu,(buffer_double_all_cpu((icpu-1)*ntime_init+itime),itime=1,ntime_init)
459  enddo
460 
461  endif
462 
463  deallocate(buffer_double)
464  deallocate(buffer_double_all_cpu)
465 
466 !
467 !
468 !***********************************************************************************************************************
469 !***********************************************************************************************************************
470 !phase 2: time loop
471 !***********************************************************************************************************************
472 !***********************************************************************************************************************
473  call mpi_barrier(mpi_comm_world,ios)
474 
475  start_time = mpi_wtime()
476 
477  if (lg_output_cputime) then
478  open(unit=get_newunit(unit_time),file=trim(cg_prefix)//".time.cpu."//trim(cg_myrank),status='replace')
479  endif
480 
481  if (ig_myrank == 0) then
482  write(ig_lst_unit,'(" ",/,a)') "starting time loop"
483  write(ig_lst_unit,'(" ",/,a)') " -->time of simulation"
484  endif
485 
486  do ig_idt = 1,ig_ndt
487 
488  !
489  !
490  !*****************************************************************************************************************
491  !->time of the simulation
492  !*****************************************************************************************************************
493  rg_simu_current_time = (ig_idt-1)*rg_dt
494 
495  !
496  !
497  !*****************************************************************************************************************
498  !->output time step in file *.lst
499  !*****************************************************************************************************************
500  if ( (ig_myrank == 0) .and. mod(ig_idt-1,1000) == 0 ) then
501  write(ig_lst_unit,'(7X,E14.7)') rg_simu_current_time
502  call flush(ig_lst_unit)
503  endif
504 
505  !
506  !
507  !*****************************************************************************************************************
508  !->compute cpu time for each cpu for each time steps
509  !*****************************************************************************************************************
510  if (lg_output_cputime) then
511  t_tloop = mpi_wtime()
512  endif
513 
514  !
515  !
516  !*****************************************************************************************************************
517  !->compute displacements at step n+1
518  !*****************************************************************************************************************
519  if (lg_output_cputime) then
520  t_newmark = mpi_wtime()
521  endif
522 
523  call newmark_ini()
524 
525  if (lg_output_cputime) then
526  t_newmark = mpi_wtime() - t_newmark
527  endif
528 
529  !
530  !
531  !*****************************************************************************************************************
532  !->compute external forces
533  !*****************************************************************************************************************
534  if (lg_output_cputime) then
535  t_forext = mpi_wtime()
536  endif
537 
539 
540  if (lg_output_cputime) then
541  t_forext = mpi_wtime() - t_forext
542  endif
543 
544  !
545  !
546  !*****************************************************************************************************************
547  !->compute boundary absorption forces
548  !*****************************************************************************************************************
549  if (lg_output_cputime) then
550  t_forabs = mpi_wtime()
551  endif
552 
553  if (lg_boundary_absorption) call compute_absorption_forces()
554 
555  if (lg_output_cputime) then
556  t_forabs = mpi_wtime() - t_forabs
557  endif
558 
559  !
560  !
561  !*****************************************************************************************************************
562  !->compute internal forces for inner and outer elements
563  !*****************************************************************************************************************
564  if (lg_output_cputime) then
565  t_forint_outer = mpi_wtime()
566  endif
567 
568  if (ig_lagrange_order == 4) call compute_internal_forces_order4(1,ig_nhexa_outer)
569  if (ig_lagrange_order == 5) call compute_internal_forces_order5(1,ig_nhexa_outer)
570  if (ig_lagrange_order == 6) call compute_internal_forces_order6(1,ig_nhexa_outer)
571 
572  if (lg_output_cputime) then
573  t_forint_outer = mpi_wtime() - t_forint_outer
574  endif
575 
576  if (lg_output_cputime) then
577  t_link_init = mpi_wtime()
578  endif
579 
580  if (lg_async_mpi_comm) call assemble_force_async_comm_init()
581 
582  if (lg_output_cputime) then
583  t_link_init = mpi_wtime() - t_link_init
584  endif
585 
586  if (lg_output_cputime) then
587  t_forint_inner = mpi_wtime()
588  endif
589 
590  if (ig_lagrange_order == 4) call compute_internal_forces_order4(ig_nhexa_outer+1,ig_nhexa)
591  if (ig_lagrange_order == 5) call compute_internal_forces_order5(ig_nhexa_outer+1,ig_nhexa)
592  if (ig_lagrange_order == 6) call compute_internal_forces_order6(ig_nhexa_outer+1,ig_nhexa)
593 
594  if (lg_output_cputime) then
595  t_forint_inner = mpi_wtime() - t_forint_inner
596  endif
597 
598  !
599  !
600  !*****************************************************************************************************************
601  !->assemble forces for linked cpu
602  !*****************************************************************************************************************
603  if (lg_output_cputime) then
604  t_linkfo = mpi_wtime()
605  endif
606 
607  if (lg_async_mpi_comm) then
609  else
611  endif
612 
613  if (lg_output_cputime) then
614  t_linkfo = mpi_wtime() - t_linkfo
615  endif
616 
617  !
618  !
619  !*****************************************************************************************************************
620  !->solve for acceleration and velocity at step n+1
621  !*****************************************************************************************************************
622  if (lg_output_cputime) then
623  t_resol = mpi_wtime()
624  endif
625 
626  call newmark_end()
627 
628  if (lg_output_cputime) then
629  t_resol = mpi_wtime() - t_resol
630  endif
631 
632  !
633  !
634  !*****************************************************************************************************************
635  !->output snapshot of the surface of the domain + compute PGD, PGV and PGA
636  !*****************************************************************************************************************
637  if (lg_output_cputime) then
638  t_snafsu = mpi_wtime()
639  endif
640 
641  if (lg_snapshot) call write_snapshot_surface()
642 
643  if (lg_snapshot_volume) call write_snapshot_volume()
644 
645  if (lg_output_cputime) then
646  t_snafsu = mpi_wtime() - t_snafsu
647  endif
648 
649  !
650  !
651  !*****************************************************************************************************************
652  !->output acceleration, velocity and displacement time history of receivers
653  !*****************************************************************************************************************
654  if (lg_output_cputime) then
655  t_outavd = mpi_wtime()
656  endif
657 
658  if ( (ig_idt == 1) .or. mod((ig_idt-1),ig_receiver_saving_incr) == 0 ) then
659 
660  call write_receiver_output()
661 
662  endif
663 
664  if (lg_output_cputime) then
665  t_outavd = mpi_wtime() - t_outavd
666  endif
667 
668  !
669  !
670  !*****************************************************************************************************************
671  !->output cpu time for each cpu for each time step + detail of different subroutines
672  !*****************************************************************************************************************
673  if (lg_output_cputime) then
674 
675  t_tloop = mpi_wtime() - t_tloop
676 
677  write(unit_time,'(i10,1x,11(e22.15,1x))') ig_idt & !1 time step
678  ,t_tloop & !2 computation of entire time loop
679  ,t_newmark & !3 computation of displacement (step n+1)
680  ,t_forext & !4 computation of external forces
681  ,t_forabs & !5 computation of boundary absorption forces
682  ,t_forint_outer & !6 computation of internal forces of outer elements
683  ,t_forint_inner & !7 computation of internal forces of inner elements
684  ,t_linkfo & !8 computation of assemble forces
685  ,t_resol & !9 computation of acceleration and velocity (step n+1)
686  ,t_snafsu & !10 computation of free surface snapshot
687  ,t_outavd & !11 computation of receiver time history
688  ,t_link_init !12 assemble_force_async_comm_init
689 
690  call flush(unit_time)
691 
692  endif
693 
694  enddo!loop on time step
695 
696 !
697 !
698 !*********************************************************************************************************************
699 !*********************************************************************************************************************
700 !phase 3: output elapsed time info and finalize computation
701 !*********************************************************************************************************************
702 !*********************************************************************************************************************
703 
704  if (lg_snapshot) then
705 
706  !***************************************************************************************************************
707  !write peak ground motion at the free surface
708  !***************************************************************************************************************
710 
711  !***************************************************************************************************************
712  !write VTK collection file if surface snapshots are saved in VTK format
713  !***************************************************************************************************************
714  if (lg_snapshot_vtk) call write_collection_vtk_surf()
715 
716  endif
717 
718  !
719  !
720  !***************************************************************************************************************
721  !write VTK collection file if volume snapshot is activated
722  !***************************************************************************************************************
723  if (lg_snapshot_volume) then
724 
726 
727  endif
728 
729  if (lg_output_cputime) then
730  close(unit_time)
731  endif
732 
733  !
734  !
735  !***************************************************************************************************************
736  !compute elapsed time for time loop for all cpu
737  !***************************************************************************************************************
738  end_time = (mpi_wtime() - start_time)
739  call mpi_reduce(end_time,dltim2,1,mpi_double_precision,mpi_max,0,mpi_comm_world,ios)
740 
741  !
742  !
743  !***************************************************************************************************************
744  !cpu 0 gather elapsed time for time loop of each cpu
745  !***************************************************************************************************************
746  allocate(buffer_double_all_cpu(ig_ncpu))
747  call mpi_gather(end_time,1,mpi_double_precision,buffer_double_all_cpu,1,mpi_double_precision,0,mpi_comm_world,ios)
748 
749  !
750  !
751  !***************************************************************************************************************
752  !cpu 0 gather number of hexa of each cpu
753  !***************************************************************************************************************
754  allocate(buffer_integer(ig_ncpu))
755  call mpi_gather(ig_nhexa,1,mpi_integer,buffer_integer,1,mpi_integer,0,mpi_comm_world,ios)
756 
757  !
758  !
759  !***************************************************************************************************************
760  !output results in file *.lst
761  !***************************************************************************************************************
762  if (ig_myrank == 0) then
763  write(ig_lst_unit,'(" ",/,a,e15.7,a)') "elapsed time for initialization = ",dltim1," s"
764  write(ig_lst_unit,'(a,e15.7,a)') "elapsed time for time loop computation = ",dltim2," s"
765  write(ig_lst_unit,'(a,e15.7,a)') "total elapsed time for computation = ",dltim1+dltim2," s"
766 
767  write(ig_lst_unit,'("",/,a,i0,a)' ) "average time per time step and per hexa (order ",ig_lagrange_order,") for the simulation"
768  do icpu = 1,ig_ncpu
769  write(ig_lst_unit,'(a,i8,1x,e15.7,a)') " -->cpu ",icpu-1,buffer_double_all_cpu(icpu)/(dble(ig_ndt)*dble(buffer_integer(icpu)))," s"
770  enddo
771  endif
772 
773  call mpi_finalize(ios)
774  stop
775 
776 end program efispec3d
subroutine, public init_input_variables()
This subroutine initializes the simulation by writing header of listing file *.lst and reading variab...
This module contains subroutines to read mesh files and creates GLL nodes global numbering in cpu myr...
This module contains subroutines to allocate arrays and to compute an approximation of the total RAM ...
This module contains subroutines to compute information about receivers and to write receivers' time ...
subroutine, public compute_external_force()
This subroutine sets external forces of the system for double couple and single force point sources...
subroutine, public compute_internal_forces_order5(elt_start, elt_end)
This subroutine computes internal forces for spectral-elements of order 5. Stress-strain relationshi...
subroutine, public init_gll_nodes()
This subroutine computes GLL nodes abscissa and weight in the reference domain [-1:1] as well as deri...
subroutine, public init_mpi()
This subroutine initializes MPI.
This module contains subroutines to initialize medium for hexahedron and quadrangle elements...
subroutine, public init_time_step()
This subroutine returns the time step mod_global_variables::rg_dt that makes the simulation stable...
This module contains subroutines to compute Newmark explicit time marching scheme, external forces , internal forces and boundary traction forces of the system .
subroutine, public init_quadp_medium()
This subroutine fills medium properties at GLL nodes of paraxial quadrangle elements (i...
program efispec3d
Definition: efispec.f90:131
This module contains subroutines to compute and write snapshots of a mesh composed by hexahedron elem...
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.
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...
subroutine, public init_jacobian_matrix_quad()
This subroutine computes the determinant of Jacobian matrix and normal unit vector of quadrangle elem...
subroutine, public init_snapshot_surface()
This subroutine generates a structured grid of receivers on the free surface.
subroutine, public init_quad_receiver()
This subroutine reads all receivers in file *.fsr; determines to which cpu they belong; computes thei...
subroutine, public compute_internal_forces_order4(elt_start, elt_end)
This subroutine computes internal forces for spectral-elements of order 4. Stress-strain relationshi...
subroutine, public memory_consumption()
Subroutine to compute an approximation of total RAM used by global variables of EFISPEC3D. See module mod_global_variables.
subroutine, public compute_double_couple_source()
This subroutine pre-computes all double couple point sources defined by type mod_global_variables::ty...
This module contains subroutines to initialize some global variable vectors and matrices.
subroutine, public init_double_couple_source()
This subroutine reads all double couple point sources in file *.dcs; sets double couple point sources...
subroutine, public init_gll_nodes_coordinates()
This subroutine computes GLL nodes -coordinates in the physical domain, cartesian right-handed coordi...
This module contains subroutines to compute information about sources.
subroutine, public init_single_force_source()
This subroutine reads all single force point sources in file *.sfs; determines to which cpu belong si...
subroutine, public init_mesh()
This subroutine reads mesh files *.inp for cpu myrank and creates GLL numbering of hexahedron and qua...
subroutine, public init_hexa_receiver()
This subroutine reads all receivers in file *.vor; determines to which cpu they belong; computes thei...
subroutine, public write_peak_ground_motion()
This subroutine writes peak ground motions files in GMT or VTK formats.
subroutine, public write_collection_vtk_vol()
This subroutine selects which ParaView collection files should be written (displacement, velocity and/or acceleration collections) depending on value of variables mod_global_variables::lg_snapshot_volume_displacement, mod_global_variables::lg_snapshot_volume_velocity and mod_global_variables::lg_snapshot_volume_acceleration.
subroutine, public write_receiver_output()
This subroutine writes -displacements, velocities and accelerations at receivers. ...
subroutine, public write_collection_vtk_surf()
This subroutine selects which ParaView collection files should be written (displacement, velocity and/or acceleration collections) depending on value of variables mod_global_variables::lg_snapshot_displacement, mod_global_variables::lg_snapshot_velocity and mod_global_variables::lg_snapshot_acceleration.
This module contains subroutines to initialize MPI buffers between cpu myrank and its neighbor cpus...
subroutine, public compute_absorption_forces()
This subroutine computes absorption forces for any spectral-elements order. A so-called 'P1' explici...
subroutine, public init_snapshot_volume()
This subroutine intializes the VTK mesh used for volume snapshots by calling 'get_efi_hexa' and 'init...
This module contains subroutines to initialize the time step of the simulation.
This module contains subroutines to assemble forces between cpu myrank and its neighbor cpus by synch...
subroutine, public init_hexa_medium()
This subroutine fills medium properties (elastic or viscoelastic) at GLL nodes of hexahedron elements...
This module contains subroutines to compute and write snapshots of the free surface (either in GMT or...
subroutine, public compute_internal_forces_order6(elt_start, elt_end)
This subroutine computes internal forces for spectral-elements of order 6. Stress-strain relationshi...
subroutine, public init_mass_matrix()
This subroutine computes and assembles the mass matrix of the system .
subroutine, public newmark_end()
This subroutine finalizes Newmark time marching scheme at step n+1.
subroutine, public assemble_force_async_comm_init()
subroutine to initialize forces assembly between connected cpus by MPI asynchrone communications...
subroutine, public init_jacobian_matrix_hexa()
This subroutine computes Jacobian matrix and its determinant for hexahedron elements.
subroutine, public write_snapshot_surface()
This subroutine computes and writes -displacements, velocities and accelerations at receivers used fo...
subroutine, public write_snapshot_volume()
This subroutine manages output format for volume snapshot. For now, only VTK format is available...
subroutine, public newmark_ini()
This subroutine initializes Newmark time marching scheme at step n+1.