All Classes Files Functions Variables Pages
module_global_variables.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 
129 
130  use mpi
131 
132  implicit none
133 
134 !
135 !
136 !***********************************************************************************************************************************************************************************
137 !parameters
138 !***********************************************************************************************************************************************************************************
139 
141  logical, parameter :: LG_VISCO = .false.
142 
144  logical, parameter :: LG_SNAPSHOT_VTK = .false.
145 
147  logical, parameter :: LG_SNAPSHOT_GMT = .true.
148 
150  logical, parameter :: LG_OUTPUT_MEDIUM_VTK = .false.
151 
153  integer, parameter :: IG_LAGRANGE_ORDER = 4
154 
156  integer, parameter :: IG_NGLL = IG_LAGRANGE_ORDER + 1
157 
159  integer, parameter :: IG_LST_UNIT = 10
160 
162  integer, parameter :: IG_NDOF = 3
163 
165  integer, parameter :: IG_NRELAX = 8
166 
168  real , parameter :: RG_NEWMARK_GAMMA = 0.5
169 
171  real , parameter :: RG_PI = 3.141592654
172 
174  real , parameter :: EPSILON_MACHINE = epsilon(EPSILON_MACHINE)
175 
177  real , parameter :: TINY_REAL = tiny(TINY_REAL)
178 
180  logical, parameter :: LG_ASYNC_MPI_COMM = .false.
181 
183  logical, parameter :: LG_OUTPUT_CPUTIME = .false.
184 
186  logical, parameter :: LG_OUTPUT_DEBUG_FILE = .false.
187 
188 !
189 !
190 !***********************************************************************************************************************************************************************************
191 !variables (scalars, or n-order tensors)
192 !***********************************************************************************************************************************************************************************
193 
195  real, dimension(:,:,:,:,:), allocatable :: rg_dcsource_gll_force
196 
198  real, dimension(:,:,:,:), allocatable :: rg_hexa_gll_jacobian_det
199 
201  real, dimension(:,:,:,:), allocatable :: rg_hexa_gll_dxidx
202 
204  real, dimension(:,:,:,:), allocatable :: rg_hexa_gll_dxidy
205 
207  real, dimension(:,:,:,:), allocatable :: rg_hexa_gll_dxidz
208 
210  real, dimension(:,:,:,:), allocatable :: rg_hexa_gll_detdx
211 
213  real, dimension(:,:,:,:), allocatable :: rg_hexa_gll_detdy
214 
216  real, dimension(:,:,:,:), allocatable :: rg_hexa_gll_detdz
217 
219  real, dimension(:,:,:,:), allocatable :: rg_hexa_gll_dzedx
220 
222  real, dimension(:,:,:,:), allocatable :: rg_hexa_gll_dzedy
223 
225  real, dimension(:,:,:,:), allocatable :: rg_hexa_gll_dzedz
226 
228  real, dimension(:,:,:,:), allocatable :: rg_hexa_gll_rho
229 
231  real, dimension(:,:,:,:), allocatable :: rg_hexa_gll_rhovs2
232 
234  real, dimension(:,:,:,:), allocatable :: rg_hexa_gll_rhovp2
235 
237  real, dimension(:,:,:), allocatable :: rg_quadp_gll_rhovs
238 
240  real, dimension(:,:,:), allocatable :: rg_quadp_gll_rhovp
241 
243  real, dimension(:,:,:), allocatable :: rg_quadp_gll_jaco_det
244 
246  real, dimension(:,:,:,:), allocatable :: rg_quadp_gll_normal
247 
249  real, dimension(:,:,:,:,:), allocatable :: rg_hexa_gll_wkqs
250 
252  real, dimension(:,:,:,:,:), allocatable :: rg_hexa_gll_wkqp
253 
255  real, dimension(:,:,:,:,:), allocatable :: rg_hexa_gll_ksixx
256 
258  real, dimension(:,:,:,:,:), allocatable :: rg_hexa_gll_ksiyy
259 
261  real, dimension(:,:,:,:,:), allocatable :: rg_hexa_gll_ksizz
262 
264  real, dimension(:,:,:,:,:), allocatable :: rg_hexa_gll_ksixy
265 
267  real, dimension(:,:,:,:,:), allocatable :: rg_hexa_gll_ksixz
268 
270  real, dimension(:,:,:,:,:), allocatable :: rg_hexa_gll_ksiyz
271 
273  real, dimension(IG_NRELAX,3), parameter :: RG_RELAX_COEFF = [ &
274  1.72333e-3, &
275  1.80701e-3, &
276  5.38887e-3, &
277  1.99322e-2, &
278  8.49833e-2, &
279  4.09335e-1, &
280  2.05951e+0, &
281  13.26290e+0, &
282  1.66958e-2, &
283  3.81644e-2, &
284  9.84666e-3, &
285  -1.36803e-2, &
286  -2.85125e-2, &
287  -5.37309e-2, &
288  -6.65035e-2, &
289  -1.33696e-1, &
290  8.98758e-2, &
291  6.84635e-2, &
292  9.67052e-2, &
293  1.20172e-1, &
294  1.30728e-1, &
295  1.38746e-1, &
296  1.40705e-1, &
297  2.14647e-1 &
298  ]
299 
301  real, dimension(IG_NRELAX) :: rg_mem_var_exp
302 
304  real, dimension(:,:), allocatable :: rg_gll_displacement
305 
307  real, dimension(:,:), allocatable :: rg_gll_velocity
308 
310  real, dimension(:,:), allocatable :: rg_gll_acceleration
311 
313  real, dimension(:), allocatable :: rg_gnode_abscissa
314 
316  real, dimension(:,:), allocatable :: rg_gnode_abscissa_dist
317 
319  real, dimension(:), allocatable :: rg_gll_mass_matrix
320 
322  real, dimension(:), allocatable :: rg_gnode_x
323 
325  real, dimension(:), allocatable :: rg_gnode_y
326 
328  real, dimension(:), allocatable :: rg_gnode_z
329 
331  real, dimension(:), allocatable :: rg_receiver_snapshot_z
332 
334  integer, dimension(:), allocatable :: ig_receiver_snapshot_locnum
335 
337  real, dimension(:), allocatable :: rg_dcsource_user_func
338 
340  real, dimension(:), allocatable :: rg_sfsource_user_func
341 
343  real, dimension(:), allocatable :: rg_mpi_buffer_send
344 
346  real, dimension(:), allocatable :: rg_mpi_buffer_recv
347 
349  real :: rg_dt = 0.0
350 
352  real :: rg_dt2 = 0.0
353 
355  real :: rg_simu_current_time = 0.0
356 
358  real :: rg_simu_total_time = 0.0
359 
361  real :: rg_mesh_xmax = 0.0
362 
364  real :: rg_mesh_xmin = 0.0
365 
367  real :: rg_mesh_ymax = 0.0
368 
370  real :: rg_mesh_ymin = 0.0
371 
373  real :: rg_mesh_zmax = 0.0
374 
376  real :: rg_mesh_zmin = 0.0
377 
379  real :: rg_receiver_snapshot_dxdy = 0.0
380 
382  real :: rg_receiver_snapshot_dx = 0.0
383 
385  real :: rg_receiver_snapshot_dy = 0.0
386 
388  integer :: ig_receiver_snapshot_nx = 0
389 
391  integer :: ig_receiver_snapshot_ny = 0
392 
394  integer, target, dimension(:,:,:,:), allocatable :: ig_hexa_gll_glonum
395 
397  integer, dimension(:,:,:), allocatable :: ig_quadp_gll_glonum
398 
400  integer, dimension(:,:,:), allocatable :: ig_quadf_gll_glonum
401 
403  integer, dimension(:,:), allocatable :: ig_hexa_gnode_glonum
404 
406  integer, dimension(:,:), allocatable :: ig_quadp_gnode_glonum
407 
409  integer, dimension(:,:), allocatable :: ig_quadf_gnode_glonum
410 
412  integer, dimension(:), allocatable :: ig_hexa_gnode_xiloc
413 
415  integer, dimension(:), allocatable :: ig_hexa_gnode_etloc
416 
418  integer, dimension(:), allocatable :: ig_hexa_gnode_zeloc
419 
421  integer, dimension(:), allocatable :: ig_quad_gnode_xiloc
422 
424  integer, dimension(:), allocatable :: ig_quad_gnode_etloc
425 
427  integer, dimension(:), allocatable :: ig_hexa_receiver_unit
428 
430  integer, dimension(:), allocatable :: ig_quad_receiver_unit
431 
433  integer :: ig_mpi_buffer_sizemax = 0
434 
436  integer, dimension(:), allocatable :: ig_mpi_request_send
437 
439  integer, dimension(:), allocatable :: ig_mpi_request_recv
440 
442  integer, dimension(:), allocatable :: ig_mpi_buffer_offset
443 
445  integer, dimension(:), allocatable :: ig_hexa_material_number
446 
448  integer(kind=1) , dimension(:), allocatable :: ig_material_type
449 
451  integer, dimension(:), allocatable :: ig_receiver_snapshot_glonum
452 
454  integer, dimension(:), allocatable :: ig_receiver_snapshot_mpi_shift
455 
457  integer, dimension(:), allocatable :: ig_receiver_snapshot_total_number
458 
460  integer :: ig_ncpu = 0
461 
463  integer :: ig_myrank = 0
464 
466  integer :: ig_ncpu_neighbor = 0
467 
469  integer, dimension(:,:), allocatable :: ig_cpu_neighbor_info
470 
472  integer :: ig_nhexa = 0
473 
475  integer :: ig_nhexa_outer = 0
476 
478  integer :: ig_nhexa_inner = 0
479 
481  integer :: ig_nquad_parax = 0
482 
484  integer :: ig_nquad_fsurf = 0
485 
487  integer :: ig_mesh_nnode = 0
488 
490  integer :: ig_hexa_nnode = 0
491 
493  integer :: ig_quad_nnode = 0
494 
496  integer :: ig_line_nnode = 0
497 
499  integer, dimension(3,8) :: ig_hexa_node2gll
500 
502  integer :: ig_ndt = 0
503 
505  integer :: ig_idt = 0
506 
508  integer :: ig_receiver_saving_incr = 1
509 
511  integer :: ig_snapshot_saving_incr = 0
512 
514  integer :: ig_snapshot_volume_saving_incr = 0
515 
517  logical :: lg_boundary_absorption = .true.
518 
520  integer :: ig_ndcsource = 0
521 
523  integer :: ig_nsfsource = 0
524 
526  integer :: ig_nreceiver_hexa = 0
527 
529  integer :: ig_nreceiver_quad = 0
530 
532  integer :: ig_nreceiver_snapshot = 0
533 
535  integer :: ig_medium_type = 0
536 
538  integer :: ig_nmaterial = 0
539 
541  integer, dimension(:), allocatable :: ig_quadp_neighbor_hexa
542 
544  integer, dimension(:), allocatable :: ig_quadp_neighbor_hexaface
545 
547  integer :: ig_nneighbor_all_kind = 0
548 
550  integer :: ig_cpu_name_len = 0
551 
553  character(len=MPI_MAX_PROCESSOR_NAME) :: cg_cpu_name
554 
556  character(len=92) :: cg_prefix
557 
559  character(len= 6) :: cg_myrank
560 
562  character(len= 6) :: cg_ncpu
563 
565  logical :: lg_snapshot = .false.
566 
568  logical :: lg_snapshot_displacement = .false.
569 
571  logical :: lg_snapshot_velocity = .false.
572 
574  logical :: lg_snapshot_acceleration = .false.
575 
577  logical :: lg_snapshot_volume = .false.
578 
580  logical :: lg_snapshot_volume_displacement = .false.
581 
583  logical :: lg_snapshot_volume_velocity = .false.
584 
586  logical :: lg_snapshot_volume_acceleration = .false.
587 
589  real, dimension(:,:), allocatable :: rg_gll_coordinate
590 
592  real, dimension(IG_NGLL) :: rg_gll_weight
593 
595  real, dimension(IG_NGLL) :: rg_gll_abscissa
596 
598  real, dimension(IG_NGLL,IG_NGLL) :: rg_gll_lagrange_deriv
599 
601  real, dimension(IG_NGLL,IG_NGLL) :: rg_gll_abscissa_dist
602 
604  integer :: ig_ngll_total = 0
605 
607  integer, dimension(:), allocatable :: ig_hexa_snapshot
608 
610  integer, dimension(:), allocatable :: ig_vtk_hexa_conn_snapshot
611 
613  integer, dimension(:), allocatable :: ig_vtk_node_gll_glonum_snapshot
614 
616  real, dimension(:), allocatable :: rg_vtk_node_x_snapshot
617 
619  real, dimension(:), allocatable :: rg_vtk_node_y_snapshot
620 
622  real, dimension(:), allocatable :: rg_vtk_node_z_snapshot
623 
625  integer(kind=1), dimension(:), allocatable :: ig_vtk_cell_type_snapshot
626 
628  integer, dimension(:), allocatable :: ig_vtk_offset_snapshot
629 
631  integer :: ig_vtk_nhexa_snapshot
632 
634  integer :: ig_vtk_nnode_snapshot
635 
637  real :: rg_snapshot_volume_xmin
638 
640  real :: rg_snapshot_volume_xmax
641 
643  real :: rg_snapshot_volume_ymin
644 
646  real :: rg_snapshot_volume_ymax
647 
649  real :: rg_snapshot_volume_zmin
650 
652  real :: rg_snapshot_volume_zmax
653 
654 !
655 !
656 !***********************************************************************************************************************************************************************************
657 !single force point sources
658 !***********************************************************************************************************************************************************************************
659 
662  real :: x
663  real :: y
664  real :: z
665  real :: fac
666  real :: rise_time
667  real :: shift_time
668  real :: var1
669  real :: dmin
670  integer :: icur
671  integer :: idir
672  integer :: cpu
673  integer :: iel
674  integer :: iequ
675  integer :: kgll
676  integer :: lgll
677  integer :: mgll
678  end type
679 
680  type(type_single_force_source), allocatable, dimension(:) :: tg_sfsource
681 
682 !
683 !
684 !***********************************************************************************************************************************************************************************
685 !type for double couple sources
686 !***********************************************************************************************************************************************************************************
687 
690  real :: x
691  real :: y
692  real :: z
693  real :: mxx
694  real :: myy
695  real :: mzz
696  real :: mxy
697  real :: mxz
698  real :: myz
699  real :: shift_time
700  real :: rise_time
701  real :: xi
702  real :: et
703  real :: ze
704  real :: dmin
705  real :: str
706  real :: dip
707  real :: rak
708  real :: mw
709  integer :: icur
710  integer :: cpu
711  integer :: iel
712  integer :: kgll
713  integer :: lgll
714  integer :: mgll
715  end type
716 
717  type(type_double_couple_source), allocatable, dimension(:) :: tg_dcsource
718 
719 !
720 !
721 !***********************************************************************************************************************************************************************************
722 !type for receivers inside hexahedron
723 !***********************************************************************************************************************************************************************************
724 
727  real :: x
728  real :: y
729  real :: z
730  real :: xi
731  real :: eta
732  real :: zeta
733  real :: dmin
734  real , dimension(IG_NGLL,IG_NGLL,IG_NGLL) :: lag
735  integer, dimension(IG_NGLL,IG_NGLL,IG_NGLL) :: gll
736  integer :: cpu
737  integer :: iel
738  integer :: kgll
739  integer :: lgll
740  integer :: mgll
741  integer :: rglo
742  end type
743 
744  type(type_receiver_hexa), allocatable, dimension(:) :: tg_receiver_hexa
745 
746 !
747 !
748 !***********************************************************************************************************************************************************************************
749 !type for receivers inside quad
750 !***********************************************************************************************************************************************************************************
751 
754  real :: x
755  real :: y
756  real :: z
757  real :: xi
758  real :: eta
759  real :: dmin
760  real , dimension(IG_NGLL,IG_NGLL) :: lag
761  real :: pgd_x
762  real :: pgd_y
763  real :: pgd_z
764  real :: pgv_x
765  real :: pgv_y
766  real :: pgv_z
767  real :: pgv_xyz
768  real :: pga_x
769  real :: pga_y
770  real :: pga_z
771  integer, dimension(IG_NGLL,IG_NGLL) :: gll
772  integer :: cpu
773  integer :: iel
774  integer :: lgll
775  integer :: mgll
776  integer :: rglo
777  end type
778 
779  type(type_receiver_quad), allocatable, dimension(:) :: tg_receiver_snapshot_quad
780 
781  type(type_receiver_quad), allocatable, dimension(:) :: tg_receiver_quad
782 
783 !
784 !
785 !***********************************************************************************************************************************************************************************
786 !type for cpus connected to cpu myrank
787 !***********************************************************************************************************************************************************************************
788 
791  integer :: ngll
792  integer :: icpu
793  integer, allocatable, dimension(:) :: gll_send
794  integer, allocatable, dimension(:) :: gll_recv
795  endtype
796 
797  type(type_cpu_neighbor), allocatable, dimension(:) :: tg_cpu_neighbor
798 
799 !
800 !
801 !***********************************************************************************************************************************************************************************
802 !type for linear elastic properties of materials
803 !***********************************************************************************************************************************************************************************
804 
807  real :: dens
808  real :: svel
809  real :: pvel
810  real :: rvel
811  real :: pois
812  real :: lam1
813  real :: lam2
814  real :: youn
815  real :: bulk
816  real :: pwmo
817  end type
818 
819  type(type_elastic_material), allocatable, dimension(:) :: tg_elastic_material
820 
821 !
822 !
823 !***********************************************************************************************************************************************************************************
824 !type for viscoelastic properties of materials
825 !***********************************************************************************************************************************************************************************
826 
829  real :: freq
830  real :: qs
831  real :: qp
832  real :: uswm
833  real :: upwm
834  real, dimension(:), allocatable :: wkqs
835  real, dimension(:), allocatable :: wkqp
836  end type
837 
838  type(type_viscoelastic_material), allocatable, dimension(:) :: tg_viscoelastic_material
839 
840 
841 !***********************************************************************************************************************************************************************************
842 !***********************************************************************************************************************************************************************************
843 !***********************************************************************************************************************************************************************************
844 !***********************************************************************************************************************************************************************************
845 
846  public :: get_newunit
847  public :: error_stop
848  public :: info_all_cpu
849  public :: sweep_blanks
850  public :: strupcase
851  public :: strlowcase
852 
853  contains
855 !
856 !
864 !***********************************************************************************************************************************************************************************
865  integer function get_newunit(myunit)
866 !***********************************************************************************************************************************************************************************
867 
868  implicit none
869 
870  integer, intent(out), optional :: myunit
871  integer, parameter :: lun_min=101
872  integer, parameter :: lun_max=20101
873  logical :: is_opened
874  integer :: lun
875 
876  get_newunit=-1
877  do lun=lun_min,lun_max
878  inquire(unit=lun,opened=is_opened)
879  if (.not. is_opened) then
880  get_newunit=lun
881  exit
882  endif
883  enddo
884  if (present(myunit)) myunit=get_newunit
885 
886 !***********************************************************************************************************************************************************************************
887  end function get_newunit
888 !***********************************************************************************************************************************************************************************
889 
890 !
891 !
894 !***********************************************************************************************************************************************************************************
895  subroutine error_stop(info,r)
896 !***********************************************************************************************************************************************************************************
897 
898  use mpi
899 
900  implicit none
901 
902  character(len=*), intent(in) :: info
903  real , intent(in), optional :: r
904  integer :: ios
905 
906  if (present(r)) then
907  write(*,'(a,i5,1x,a,E14.7)') "cpu ",ig_myrank,trim(adjustl(info))//" ",r
908  else
909  write(*,'(a,i5,1x,a)') "cpu ",ig_myrank,trim(adjustl(info))
910  endif
911 
912  call mpi_abort(mpi_comm_world,100,ios)
913  stop
914 
915 !***********************************************************************************************************************************************************************************
916  end subroutine error_stop
917 !***********************************************************************************************************************************************************************************
918 
919 !
920 !
924 !***********************************************************************************************************************************************************************************
925  subroutine info_all_cpu(i,info)
926 !***********************************************************************************************************************************************************************************
927 
928  use mpi
929 
930  implicit none
931 
932  integer , intent(in) :: i
933  character(len=255), intent(in) :: info
934 
935  integer, dimension(ig_ncpu) :: buffer
936  integer :: ios
937  integer :: icpu
938  integer(kind=8) :: isum
939 
940  isum = 0
941 
942  call mpi_gather(i,1,mpi_integer,buffer,1,mpi_integer,0,mpi_comm_world,ios)
943 
944  if (ig_myrank == 0) then
945 
946  write(ig_lst_unit,'(a)') " "
947 
948  do icpu = 1,ig_ncpu
949  write(ig_lst_unit,'(a,i6,a,i12,1x,a)') "cpu ",icpu-1," computes ",buffer(icpu),trim(info)
950  isum = isum + buffer(icpu)
951  enddo
952 
953  write(ig_lst_unit,'(" -----------------------",/,a,i29)') "total = ",isum
954 
955  endif
956 
957  return
958 !***********************************************************************************************************************************************************************************
959  end subroutine info_all_cpu
960 !***********************************************************************************************************************************************************************************
961 
962 !
963 !
966 !***********************************************************************************************************************************************************************************
967  character(len=100) function sweep_blanks(in_str)
968 !***********************************************************************************************************************************************************************************
969 
970  implicit none
971 
972  character(len=100), intent(in) :: in_str
973  character(len=100) :: out_str
974  character(len=1) :: ch
975  integer :: j
976 
977  out_str = " "
978  do j=1, len_trim(in_str)
979 !
980 !----->get j-th char
981  ch = in_str(j:j)
982 
983  if (ch .ne. " ") then
984  out_str = trim(out_str) // ch
985  endif
986  sweep_blanks = out_str
987  enddo
988 !***********************************************************************************************************************************************************************************
989  end function sweep_blanks
990 !***********************************************************************************************************************************************************************************
991 
992 !
993 !
996 !***********************************************************************************************************************************************************************************
997  function strupcase (input_string) result (output_string)
998 !***********************************************************************************************************************************************************************************
999 
1000  implicit none
1001 
1002  character(len=*), parameter :: lower_case = 'abcdefghijklmnopqrstuvwxyz'
1003  character(len=*), parameter :: upper_case = 'ABCDEFGHIJKLMNOPQRSTUVWXYZ'
1004 
1005  character(len=*), intent( in ) :: input_string
1006  character(len=len(input_string)) :: output_string
1007 
1008  integer :: i
1009  integer :: n
1010 
1011  output_string = input_string
1012 
1013  do i = 1,len(output_string)
1014  n = index(lower_case,output_string(i:i))
1015  if (n /= 0) output_string(i:i) = upper_case(n:n)
1016  enddo
1017 
1018 !***********************************************************************************************************************************************************************************
1019  end function strupcase
1020 !***********************************************************************************************************************************************************************************
1021 
1022 !
1023 !
1026 !***********************************************************************************************************************************************************************************
1027  function strlowcase (input_string) result (output_string)
1028 !***********************************************************************************************************************************************************************************
1029 
1030  implicit none
1031 
1032  character(len=*), parameter :: lower_case = 'abcdefghijklmnopqrstuvwxyz'
1033  character(len=*), parameter :: upper_case = 'ABCDEFGHIJKLMNOPQRSTUVWXYZ'
1034 
1035  character(len=*), intent(in) :: input_string
1036  character(len=len(input_string)) :: output_string
1037 
1038  integer :: i
1039  integer :: n
1040 
1041  output_string = input_string
1042 
1043  do i = 1,len(output_string)
1044  n = index(upper_case,output_string(i:i))
1045  if (n /= 0) output_string(i:i) = lower_case(n:n)
1046  enddo
1047 
1048 !***********************************************************************************************************************************************************************************
1049  end function strlowcase
1050 !***********************************************************************************************************************************************************************************
1052 
1053 end module mod_global_variables
type that gathers information about cpus connected to cpu myrank (cpus not connected to cpu myrank ar...
type for receivers (i.e., stations) located inside quadrangle elements
This module defines all global variables of EFISPEC3D. Scalar variables are initialized directly in t...
type for viscoelastic properties of materials
type for receivers (i.e., stations) located inside hexahedron elements
type for linear elastic properties of materials