[cig-commits] r12651 - in seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta: setup src

dkomati1 at geodynamics.org dkomati1 at geodynamics.org
Fri Aug 15 19:44:41 PDT 2008


Author: dkomati1
Date: 2008-08-15 19:44:40 -0700 (Fri, 15 Aug 2008)
New Revision: 12651

Modified:
   seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/setup/constants.h
   seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/memory_eval.f90
   seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_1066a.f90
   seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_ak135.f90
   seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_iasp91.f90
   seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_prem.f90
   seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_ref.f90
   seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/read_compute_parameters.F90
   seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/specfem3D.F90
   seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/write_seismograms.F90
Log:
suppressed flag PATCH_TIKIR_PARTLY_RESTORE which is now useless because we will
not cut the mesh into multiples of 8 yet.
added flag PUT_SOURCE_IN_EACH_REGION to put a source in the outer core and in the inner
core as well for the stability tests performed with the serial code.


Modified: seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/setup/constants.h
===================================================================
--- seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/setup/constants.h	2008-08-16 01:07:06 UTC (rev 12650)
+++ seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/setup/constants.h	2008-08-16 02:44:40 UTC (rev 12651)
@@ -181,6 +181,9 @@
 ! flag to only create the mesh but not start the solver (for instance to check the mesh obtained)
   logical, parameter :: MESHER_ONLY = .false.
 
+! flag to put a fictitious source in each region in the case of a serial test
+  logical, parameter :: PUT_SOURCE_IN_EACH_REGION = .true.
+
 !------------------------------------------------------
 !----------- do not modify anything below -------------
 !------------------------------------------------------
@@ -454,8 +457,6 @@
 
 ! to suppress the crustal layers (replaced by an extension of the mantle: R_EARTH is not modified, but no more crustal doubling)
   logical, parameter :: SUPPRESS_CRUSTAL_MESH = .false.
-!! DK DK temporary patch to cut in multiples of 8 instead of 16; will do better later
-  logical, parameter :: PATCH_TIKIR_PARTLY_RESTORE = .false.
 
 ! to add a fourth doubling at the bottom of the outer core
   logical, parameter :: ADD_4TH_DOUBLING = .false.

Modified: seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/memory_eval.f90
===================================================================
--- seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/memory_eval.f90	2008-08-16 01:07:06 UTC (rev 12650)
+++ seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/memory_eval.f90	2008-08-16 02:44:40 UTC (rev 12651)
@@ -102,8 +102,7 @@
 !! DK DK the above code by David is incorrect when SUPPRESS_CRUSTAL_MESH = .true.
 !! DK DK therefore apply a temporary patch for now
 !! DK DK we should fix the above code later
-!! DK DK temporary patch to cut in multiples of 8 instead of 16; will do better later
-  if(SUPPRESS_CRUSTAL_MESH .and. PATCH_TIKIR_PARTLY_RESTORE) nspec_tiso = NSPEC(IREGION_CRUST_MANTLE)
+  if(SUPPRESS_CRUSTAL_MESH) nspec_tiso = NSPEC(IREGION_CRUST_MANTLE)
 
 ! define static size of the arrays whose size depends on logical tests
 

Modified: seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_1066a.f90
===================================================================
--- seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_1066a.f90	2008-08-16 01:07:06 UTC (rev 12650)
+++ seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_1066a.f90	2008-08-16 02:44:40 UTC (rev 12651)
@@ -788,7 +788,7 @@
   M1066a_V%vs_1066a(159) =   2.58140000000000
   M1066a_V%vs_1066a(160) =   2.58220000000000
 
-  if (SUPPRESS_CRUSTAL_MESH .and. .not. PATCH_TIKIR_PARTLY_RESTORE) then
+  if (SUPPRESS_CRUSTAL_MESH) then
     M1066a_V%vp_1066a(158:160) = M1066a_V%vp_1066a(157)
     M1066a_V%vs_1066a(158:160) = M1066a_V%vs_1066a(157)
     M1066a_V%density_1066a(158:160) = M1066a_V%density_1066a(157)

Modified: seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_ak135.f90
===================================================================
--- seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_ak135.f90	2008-08-16 01:07:06 UTC (rev 12650)
+++ seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_ak135.f90	2008-08-16 02:44:40 UTC (rev 12651)
@@ -727,7 +727,7 @@
   Mak135_V%vs_ak135(143) =   3.20000000000000
   Mak135_V%vs_ak135(144) =   3.20000000000000
 
-  if (SUPPRESS_CRUSTAL_MESH .and. .not. PATCH_TIKIR_PARTLY_RESTORE) then
+  if (SUPPRESS_CRUSTAL_MESH) then
     Mak135_V%vp_ak135(137:144) = Mak135_V%vp_ak135(136)
     Mak135_V%vs_ak135(137:144) = Mak135_V%vs_ak135(136)
     Mak135_V%density_ak135(137:144) = Mak135_V%density_ak135(136)

Modified: seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_iasp91.f90
===================================================================
--- seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_iasp91.f90	2008-08-16 01:07:06 UTC (rev 12650)
+++ seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_iasp91.f90	2008-08-16 02:44:40 UTC (rev 12651)
@@ -186,7 +186,7 @@
       Qmu=600.0d0
       Qkappa=57827.0d0
 
-  else if (SUPPRESS_CRUSTAL_MESH .and. .not. PATCH_TIKIR_PARTLY_RESTORE) then
+  else if (SUPPRESS_CRUSTAL_MESH) then
 !! DK DK extend the Moho up to the surface instead of the crust
           vp = 8.78541d0-0.74953d0*(RMOHO / R_EARTH)
           vs = 6.706231d0-2.248585d0*(RMOHO / R_EARTH)

Modified: seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_prem.f90
===================================================================
--- seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_prem.f90	2008-08-16 01:07:06 UTC (rev 12650)
+++ seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_prem.f90	2008-08-16 02:44:40 UTC (rev 12651)
@@ -180,7 +180,7 @@
     Qmu=80.0d0
     Qkappa=57827.0d0
   else
-  if(CRUSTAL .and. (.not. SUPPRESS_CRUSTAL_MESH .or. PATCH_TIKIR_PARTLY_RESTORE)) then
+  if(CRUSTAL .and. .not. SUPPRESS_CRUSTAL_MESH) then
 ! fill with PREM mantle and later add CRUST2.0
     if(r > R80) then
       drhodr=0.6924d0
@@ -201,7 +201,7 @@
       Qkappa=57827.0d0
 
 
-    else if (SUPPRESS_CRUSTAL_MESH .and. .not. PATCH_TIKIR_PARTLY_RESTORE) then
+    else if (SUPPRESS_CRUSTAL_MESH) then
 !! DK DK extend the Moho up to the surface instead of the crust
       drhodr=0.6924d0
       rho = 2.6910d0+0.6924d0*(RMOHO / R_EARTH)

Modified: seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_ref.f90
===================================================================
--- seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_ref.f90	2008-08-16 01:07:06 UTC (rev 12650)
+++ seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/model_ref.f90	2008-08-16 02:44:40 UTC (rev 12651)
@@ -7361,7 +7361,7 @@
  1.00000000000000 , &
  1.00000000000000 /)
 
-  if (SUPPRESS_CRUSTAL_MESH .and. .not. PATCH_TIKIR_PARTLY_RESTORE) then
+  if (SUPPRESS_CRUSTAL_MESH) then
     Mref_V%density_ref(718:750) = Mref_V%density_ref(717)
     Mref_V%vpv_ref(718:750) = Mref_V%vpv_ref(717)
     Mref_V%vph_ref(718:750) = Mref_V%vph_ref(717)

Modified: seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/read_compute_parameters.F90
===================================================================
--- seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/read_compute_parameters.F90	2008-08-16 01:07:06 UTC (rev 12650)
+++ seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/read_compute_parameters.F90	2008-08-16 02:44:40 UTC (rev 12651)
@@ -1141,24 +1141,14 @@
   CUT_SUPERBRICK_XI = .false.
   CUT_SUPERBRICK_ETA = .false.
 
-!! DK DK temporary patch
-  if(PATCH_TIKIR_PARTLY_RESTORE .and. .not. SUPPRESS_CRUSTAL_MESH) &
-    stop 'temporary patch: must set SUPPRESS_CRUSTAL_MESH = .true. if PATCH_TIKIR_PARTLY_RESTORE = .true.'
-
   if (SUPPRESS_CRUSTAL_MESH .and. .not. ADD_4TH_DOUBLING) then
 
     if(mod(NEX_XI,8) /= 0) stop 'NEX_XI must be a multiple of 8'
     if(mod(NEX_ETA,8) /= 0) stop 'NEX_ETA must be a multiple of 8'
-!! DK DK temporary patch
- if(PATCH_TIKIR_PARTLY_RESTORE) then
-    if(mod(NEX_XI/8,NPROC_XI) /= 0) stop 'NEX_XI must be a multiple of 8*NPROC_XI'
-    if(mod(NEX_ETA/8,NPROC_ETA) /= 0) stop 'NEX_ETA must be a multiple of 8*NPROC_ETA'
- else
     if(mod(NEX_XI/4,NPROC_XI) /= 0) stop 'NEX_XI must be a multiple of 4*NPROC_XI'
     if(mod(NEX_ETA/4,NPROC_ETA) /= 0) stop 'NEX_ETA must be a multiple of 4*NPROC_ETA'
     if(mod(NEX_XI/8,NPROC_XI) /= 0) CUT_SUPERBRICK_XI = .true.
     if(mod(NEX_ETA/8,NPROC_ETA) /= 0) CUT_SUPERBRICK_ETA = .true.
- endif
 
   else if (SUPPRESS_CRUSTAL_MESH .or. .not. ADD_4TH_DOUBLING) then
 

Modified: seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/specfem3D.F90
===================================================================
--- seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/specfem3D.F90	2008-08-16 01:07:06 UTC (rev 12650)
+++ seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/specfem3D.F90	2008-08-16 02:44:40 UTC (rev 12651)
@@ -2108,6 +2108,21 @@
             NGLOB2DMAX_XMIN_XMAX(IREGION_OUTER_CORE),NGLOB2DMAX_YMIN_YMAX(IREGION_OUTER_CORE),NGLOB2DMAX_XY_VAL_OC,NCHUNKS)
 #endif
 
+#ifndef USE_MPI
+!! DK DK put a fictitious source in each region in the case of a serial test if needed
+  if(PUT_SOURCE_IN_EACH_REGION) then
+    stf = 1.d-6 * comp_source_time_function(dble(it-1)*DT-t0,10.d0)
+! distinguish between single and double precision for reals
+    if(CUSTOM_REAL == SIZE_REAL) then
+      stf_used = sngl(stf)
+    else
+      stf_used = stf
+    endif
+    iglob = ibool_outer_core(2,2,2,2)
+    accel_outer_core(iglob) = accel_outer_core(iglob) + stf_used
+  endif
+#endif
+
 ! multiply by the inverse of the mass matrix and update velocity
   do i=1,NGLOB_OUTER_CORE
     accel_outer_core(i) = accel_outer_core(i)*rmass_outer_core(i)
@@ -2393,6 +2408,21 @@
     veloc_crust_mantle(:,i) = veloc_crust_mantle(:,i) + deltatover2*accel_crust_mantle(:,i)
   enddo
 
+#ifndef USE_MPI
+!! DK DK put a fictitious source in each region in the case of a serial test if needed
+  if(PUT_SOURCE_IN_EACH_REGION) then
+    stf = 1.d-6 * comp_source_time_function(dble(it-1)*DT-t0,10.d0)
+! distinguish between single and double precision for reals
+    if(CUSTOM_REAL == SIZE_REAL) then
+      stf_used = sngl(stf)
+    else
+      stf_used = stf
+    endif
+    iglob = ibool_inner_core(2,2,2,2)
+    accel_inner_core(3,iglob) = accel_inner_core(3,iglob) + stf_used
+  endif
+#endif
+
   do i=1,NGLOB_INNER_CORE
     accel_inner_core(1,i) = accel_inner_core(1,i)*rmass_inner_core(i)
     accel_inner_core(2,i) = accel_inner_core(2,i)*rmass_inner_core(i)
@@ -2459,6 +2489,7 @@
 
 ! write the current or final seismograms
   if(seismo_current == NTSTEP_BETWEEN_OUTPUT_SEISMOS .or. it == it_end) then
+
       call write_seismograms(myrank,seismograms,number_receiver_global,station_name, &
             network_name,stlat,stlon,stele,nrec,nrec_local,DT,t0,it_end, &
             yr_SAC,jda_SAC,ho_SAC,mi_SAC,sec_SAC,t_cmt_SAC, &
@@ -2468,13 +2499,20 @@
             OUTPUT_SEISMOS_SAC_BINARY,ROTATE_SEISMOGRAMS_RT,NTSTEP_BETWEEN_OUTPUT_SEISMOS, &
             seismo_offset,seismo_current,WRITE_SEISMOGRAMS_BY_MASTER, &
             SAVE_ALL_SEISMOS_IN_ONE_FILE,USE_BINARY_FOR_LARGE_FILE,one_seismogram)
+
       if(myrank==0) then
         write(IMAIN,*)
-        write(IMAIN,*) ' Total number of time steps written: ', it-it_begin+1
+        write(IMAIN,*) ' Total number of time steps written: ',it-it_begin+1
         write(IMAIN,*)
       endif
+
+! prepare to shift to the next interval to store seismograms
     seismo_offset = seismo_offset + seismo_current
     seismo_current = 0
+
+! clean seismogram array
+    seismograms(:,:,:) = 0._CUSTOM_REAL
+
   endif
 #endif
 

Modified: seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/write_seismograms.F90
===================================================================
--- seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/write_seismograms.F90	2008-08-16 01:07:06 UTC (rev 12650)
+++ seismo/3D/SPECFEM3D_GLOBE/trunk/version41_beta/src/write_seismograms.F90	2008-08-16 02:44:40 UTC (rev 12651)
@@ -88,6 +88,7 @@
 ! new flags to decide on seismogram type BS BS 06/2007
   logical OUTPUT_SEISMOS_ASCII_TEXT, OUTPUT_SEISMOS_SAC_ALPHANUM, &
           OUTPUT_SEISMOS_SAC_BINARY
+
 ! flag whether seismograms are ouput for North-East-Z component or Radial-Transverse-Z
   logical ROTATE_SEISMOGRAMS_RT
 



More information about the cig-commits mailing list