-
Notifications
You must be signed in to change notification settings - Fork 146
/
adj_seismogram_Tromp2005.f90
240 lines (192 loc) · 7.83 KB
/
adj_seismogram_Tromp2005.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
!========================================================================
!
! S P E C F E M 2 D
! -----------------
!
! Main historical authors: Dimitri Komatitsch and Jeroen Tromp
! CNRS, France
! and Princeton University, USA
! (there are currently many more authors!)
! (c) October 2017
!
! This software is a computer program whose purpose is to solve
! the two-dimensional viscoelastic anisotropic or poroelastic wave equation
! using a spectral-element method (SEM).
!
! This program is free software; you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation; either version 3 of the License, or
! (at your option) any later version.
!
! This program is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! You should have received a copy of the GNU General Public License along
! with this program; if not, write to the Free Software Foundation, Inc.,
! 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
!
! The full text of the license is available in file "LICENSE".
!
!========================================================================
program adj_seismogram
! This program cuts a certain portion of the seismograms and convert it
! into the adjoint source for generating banana-dougnut kernels
! rm -rf xadj_seismogram ; gfortran adj_seismogram.f90 -o xadj_seismogram
implicit none
!--------------------------------------------------------------------------------
! USER PARAMETERS
!--------------------------------------------------------------------------------
! please check these parameters with your setting in Par_file
integer, parameter :: NSTEP = 3000
integer, parameter :: nrec = 1
double precision, parameter :: t0 = 8.0 ! labeled as 'time delay'
double precision, parameter :: deltat = 0.02
!--------------------------------------------------------------------------------
double precision, parameter :: EPS = 1.d-40
integer :: itime,icomp,istart,iend,nlen,irec,NDIM,NDIMr,adj_comp
integer :: ier
double precision :: time,tstart(nrec),tend(nrec),dummy,time0
double precision, dimension(NSTEP) :: time_window
double precision :: seism(NSTEP,3),Nnorm,seism_win(NSTEP)
double precision :: seism_veloc(NSTEP),seism_accel(NSTEP),ft_bar(NSTEP)
character(len=150), dimension(nrec) :: station_name
character(len=3) :: compr(2),comp(3)
character(len=150) :: filename
NDIM = 3
comp = (/"BXX","BXY","BXZ"/)
! number of components
NDIMr = 2 ! P-SV
!NDIMr=1 ! SH (membrane)
compr = (/"BXX","BXZ"/) ! P-SV
!compr = (/"BXY","tmp"/) ! SH (membrane)
! list of stations
station_name(1) = 'S0001'
! KEY: 'absolute' time interval for window
tstart(1) = 27.d0 + t0
tend(1) = 32.d0 + t0
! chose the component for the adjoint source (adj_comp = 1:X, 2:Y, 3:Z)
adj_comp = 1
! user output
print *,'adjoint source - seismogram:'
print *
print *,'parameters:'
print *,' NSTEP = ',NSTEP
print *,' deltat = ',deltat
print *,' nrec = ',nrec
print *,' t0 = ',t0
print *
print *,'setup:'
print *,' seismogram components = ',NDIMr
print *,' seismogram labels = ',compr(1),' / ',compr(2)
print *,' station name = ',trim(station_name(1))
print *
print *,' time window start/end = ',tstart(1),tend(1)
print *,' adjoint source trace component (1 == X/2 == Y/3 == Z) = ',adj_comp
print *
do irec = 1,nrec
do icomp = 1,NDIMr
filename = 'OUTPUT_FILES/'//'AA.'//trim(station_name(irec))//'.'// compr(icomp) // '.semd'
open(unit = 10, file = trim(filename),status='old',iostat=ier)
if (ier /= 0) then
print *,'Error opening file: ',trim(filename)
stop 'Error opening trace file'
endif
! counts lines
nlen = 0
do while (ier == 0)
read(10,*,iostat=ier) time , dummy
if (ier == 0) nlen = nlen + 1
enddo
rewind(10)
! checks with settings
if (nlen /= NSTEP) then
print *,'Error: number of lines ',nlen,' in file ',trim(filename),' do not match setting NSTEP = ',NSTEP
stop 'Error mismatch NSTEP'
endif
! reads in trace
do itime = 1,NSTEP
read(10,*) time , seism(itime,icomp)
! stores start time for time step size
if (itime == 1) time0 = time
! checks time step size
if (itime == 2) then
if (abs((time - time0)-deltat) > 1.e-5) then
print *,'Error: time step size ',(time-time0),' in file ',trim(filename),' does not math deltat = ',deltat
stop 'Error mismatch deltat'
endif
endif
enddo
enddo
if (NDIMr == 2) then
seism(:,3) = seism(:,2)
seism(:,2) = 0.d0
else
seism(:,2) = seism(:,1)
seism(:,1) = 0.d0
seism(:,3) = 0.d0
endif
close(10)
istart = max(floor(tstart(irec)/deltat),1)
iend = min(floor(tend(irec)/deltat),NSTEP)
print *,'istart =',istart, 'iend =', iend
print *,'tstart =',istart*deltat, 'tend =', iend*deltat
if (istart >= iend) stop 'check istart,iend'
nlen = iend - istart +1
do icomp = 1, NDIM
print *,comp(icomp)
filename = 'SEM/'//'AA.'//trim(station_name(irec))//'.'// comp(icomp) // '.adj'
open(unit = 11, file = trim(filename),status='unknown',iostat=ier)
if (ier /= 0) then
print *,'Error opening file: ',trim(filename)
stop 'Error opening SEM adjoint source file'
endif
time_window(:) = 0.d0
seism_win(:) = seism(:,icomp)
seism_veloc(:) = 0.d0
seism_accel(:) = 0.d0
do itime = istart,iend
! cosine window
!time_window(itime) = 1.d0 - cos(pi*(itime-1)/NSTEP+1)**10
! Welch window
time_window(itime) = 1.d0 - (2* (dble(itime) - istart)/(iend-istart) -1.d0)**2
enddo
! first time derivative (by finite-differences)
do itime = 2,NSTEP-1
seism_veloc(itime) = (seism_win(itime+1) - seism_win(itime-1))/(2*deltat)
enddo
seism_veloc(1) = (seism_win(2) - seism_win(1))/deltat
seism_veloc(NSTEP) = (seism_win(NSTEP) - seism_win(NSTEP-1))/deltat
! second time derivative
do itime = 2,NSTEP-1
seism_accel(itime) = (seism_veloc(itime+1) - seism_veloc(itime-1))/(2*deltat)
enddo
seism_accel(1) = (seism_veloc(2) - seism_veloc(1))/deltat
seism_accel(NSTEP) = (seism_veloc(NSTEP) - seism_veloc(NSTEP-1))/deltat
! cross-correlation traveltime adjoint source
Nnorm = deltat * sum(time_window(:) * seism_win(:) * seism_accel(:))
!Nnorm = deltat * sum(time_window(:) * seism_veloc(:) * seism_veloc(:))
if (abs(Nnorm) > EPS) then
!ft_bar(:) = -seism_veloc(:) * time_window(:) / Nnorm
ft_bar(:) = seism_veloc(:) * time_window(:) / Nnorm
print *,'Norm =', Nnorm
else
print *,'Norm < EPS for file, zeroeing out trace'
print *,'Norm =', Nnorm
ft_bar(:) = 0.d0
endif
do itime = 1,NSTEP
if (icomp == adj_comp) then
write(11,*) (itime-1)*deltat - t0, ft_bar(itime)
else
write(11,*) (itime-1)*deltat - t0, 0.d0
endif
enddo
enddo
close(11)
enddo
print *,'*************************'
print *,'The input files (AA.S****.BXX/BXY/BXZ.adj) needed to run the adjoint simulation are in folder SEM/'
print *,'*************************'
end program adj_seismogram