1# http://pyrocko.org - GPLv3 

2# 

3# The Pyrocko Developers, 21st Century 

4# ---|P------/S----------~Lg---------- 

5from __future__ import absolute_import, division 

6 

7import numpy as num 

8import logging 

9import os 

10import shutil 

11import math 

12import copy 

13import signal 

14 

15from tempfile import mkdtemp 

16from subprocess import Popen, PIPE 

17from os.path import join as pjoin 

18 

19from pyrocko import gf 

20from pyrocko import trace, util, cake 

21from pyrocko.moment_tensor import MomentTensor, symmat6 

22from pyrocko.guts import Float, Int, Tuple, List, Complex, Bool, Object, String 

23 

24km = 1e3 

25 

26guts_prefix = 'pf' 

27 

28Timing = gf.meta.Timing 

29 

30logger = logging.getLogger('pyrocko.fomosto.qseis') 

31 

32# how to call the programs 

33program_bins = { 

34 'qseis.2006': 'fomosto_qseis2006', 

35 'qseis.2006a': 'fomosto_qseis2006a', 

36} 

37 

38 

39def have_backend(): 

40 have_any = False 

41 for cmd in [[exe] for exe in program_bins.values()]: 

42 try: 

43 p = Popen(cmd, stdout=PIPE, stderr=PIPE, stdin=PIPE) 

44 (stdout, stderr) = p.communicate() 

45 have_any = True 

46 

47 except OSError: 

48 pass 

49 

50 return have_any 

51 

52 

53qseis_components = 'r t z v'.split() 

54qseis_greenf_names = ('ex', 'ss', 'ds', 'cl', 'fz', 'fh') 

55 

56 

57def nextpow2(i): 

58 return 2**int(math.ceil(math.log(i)/math.log(2.))) 

59 

60 

61def str_float_vals(vals): 

62 return ' '.join('%e' % val for val in vals) 

63 

64 

65def str_int_vals(vals): 

66 return ' '.join('%i' % val for val in vals) 

67 

68 

69def str_str_vals(vals): 

70 return ' '.join("'%s'" % val for val in vals) 

71 

72 

73def scl(cs): 

74 if not cs: 

75 return '\n#' 

76 

77 return '\n'+' '.join('(%e,%e)' % (c.real, c.imag) for c in cs) 

78 

79 

80def cake_model_to_config(mod): 

81 k = 1000. 

82 srows = [] 

83 ref_depth = 0. 

84 for i, row in enumerate(mod.to_scanlines()): 

85 depth, vp, vs, rho, qp, qs = row 

86 if i == 0: 

87 ref_depth = depth 

88 

89 row = [(depth-ref_depth)/k, vp/k, vs/k, rho/k, qp, qs] 

90 srows.append('%i %s' % (i+1, str_float_vals(row))) 

91 

92 return '\n'.join(srows), len(srows), ref_depth 

93 

94 

95class QSeisSourceMech(Object): 

96 pass 

97 

98 

99class QSeisSourceMechMT(QSeisSourceMech): 

100 mnn = Float.T(default=1.0) 

101 mee = Float.T(default=1.0) 

102 mdd = Float.T(default=1.0) 

103 mne = Float.T(default=0.0) 

104 mnd = Float.T(default=0.0) 

105 med = Float.T(default=0.0) 

106 

107 def string_for_config(self): 

108 return '1 %(mnn)15e %(mee)15e %(mdd)15e ' \ 

109 '%(mne)15e %(med)15e %(mnd)15e ' % self.__dict__ 

110 

111 

112class QSeisSourceMechSDR(QSeisSourceMech): 

113 m_iso = Float.T(default=0.0) 

114 m_clvd = Float.T(default=0.0) 

115 m_dc = Float.T(default=1.0e9) 

116 strike = Float.T(default=0.0) 

117 dip = Float.T(default=90.0) 

118 rake = Float.T(default=0.0) 

119 

120 def string_for_config(self): 

121 return '2 %(m_iso)15e %(m_clvd)15e %(m_dc)15e ' \ 

122 '%(strike)15e %(dip)15e %(rake)15e ' % self.__dict__ 

123 

124 

125class QSeisPropagationFilter(Object): 

126 min_depth = Float.T(default=0.0) 

127 max_depth = Float.T(default=0.0) 

128 filtered_phase = Int.T(default=0) 

129 

130 def string_for_config(self): 

131 return '%(min_depth)15e %(max_depth)15e ' \ 

132 '%(filtered_phase)i' % self.__dict__ 

133 

134 

135class QSeisPoleZeroFilter(Object): 

136 constant = Complex.T(default=(1+0j)) 

137 poles = List.T(Complex.T()) 

138 zeros = List.T(Complex.T()) 

139 

140 def string_for_config(self, version=None): 

141 if version == '2006a': 

142 return '(%e,%e)\n%i%s\n%i%s' % ( 

143 self.constant.real, self.constant.imag, 

144 len(self.zeros), scl(self.zeros), 

145 len(self.poles), scl(self.poles)) 

146 elif version == '2006': 

147 return '%e\n%i%s\n%i%s' % ( 

148 abs(self.constant), 

149 len(self.zeros), scl(self.zeros), 

150 len(self.poles), scl(self.poles)) 

151 

152 

153class QSeisConfig(Object): 

154 

155 qseis_version = String.T(default='2006') 

156 time_region = Tuple.T(2, Timing.T(), default=( 

157 Timing('-10'), Timing('+890'))) 

158 

159 cut = Tuple.T(2, Timing.T(), optional=True) 

160 fade = Tuple.T(4, Timing.T(), optional=True) 

161 relevel_with_fade_in = Bool.T(default=False) 

162 

163 sw_algorithm = Int.T(default=0) 

164 slowness_window = Tuple.T(4, Float.T(default=0.0)) 

165 wavenumber_sampling = Float.T(default=2.5) 

166 aliasing_suppression_factor = Float.T(default=0.1) 

167 

168 filter_surface_effects = Int.T(default=0) 

169 filter_shallow_paths = Int.T(default=0) 

170 filter_shallow_paths_depth = Float.T(default=0.0) 

171 propagation_filters = List.T(QSeisPropagationFilter.T()) 

172 receiver_filter = QSeisPoleZeroFilter.T(optional=True) 

173 

174 sw_flat_earth_transform = Int.T(default=0) 

175 

176 gradient_resolution_vp = Float.T(default=0.0) 

177 gradient_resolution_vs = Float.T(default=0.0) 

178 gradient_resolution_density = Float.T(default=0.0) 

179 

180 wavelet_duration_samples = Float.T(default=0.0) 

181 wavelet_type = Int.T(default=2) 

182 user_wavelet_samples = List.T(Float.T()) 

183 

184 def items(self): 

185 return dict(self.T.inamevals(self)) 

186 

187 

188class QSeisConfigFull(QSeisConfig): 

189 

190 time_start = Float.T(default=0.0) 

191 time_reduction_velocity = Float.T(default=0.0) 

192 time_window = Float.T(default=900.0) 

193 

194 source_depth = Float.T(default=10.0) 

195 source_mech = QSeisSourceMech.T( 

196 optional=True, 

197 default=QSeisSourceMechMT.D()) 

198 

199 receiver_depth = Float.T(default=0.0) 

200 receiver_distances = List.T(Float.T()) 

201 nsamples = Int.T(default=256) 

202 

203 gf_sw_source_types = Tuple.T(6, Int.T(), default=(1, 1, 1, 1, 0, 0)) 

204 

205 gf_filenames = Tuple.T(6, String.T(), default=qseis_greenf_names) 

206 

207 seismogram_filename = String.T(default='seis') 

208 

209 receiver_azimuths = List.T(Float.T()) 

210 

211 earthmodel_1d = gf.meta.Earthmodel1D.T(optional=True) 

212 earthmodel_receiver_1d = gf.meta.Earthmodel1D.T(optional=True) 

213 

214 @staticmethod 

215 def example(): 

216 conf = QSeisConfigFull() 

217 conf.receiver_distances = [2000.] 

218 conf.receiver_azimuths = [0.] 

219 conf.time_start = -10.0 

220 conf.time_reduction_velocity = 15.0 

221 conf.earthmodel_1d = cake.load_model().extract(depth_max='cmb') 

222 conf.earthmodel_receiver_1d = None 

223 conf.sw_flat_earth_transform = 1 

224 return conf 

225 

226 def get_output_filenames(self, rundir): 

227 return [pjoin(rundir, self.seismogram_filename+'.t'+c) 

228 for c in qseis_components] 

229 

230 def get_output_filenames_gf(self, rundir): 

231 return [pjoin(rundir, fn+'.t'+c) 

232 for fn in self.gf_filenames for c in qseis_components] 

233 

234 def string_for_config(self): 

235 

236 def aggregate(xx): 

237 return len(xx), '\n'.join( 

238 [''] + [x.string_for_config() for x in xx]) 

239 

240 assert len(self.receiver_distances) > 0 

241 assert len(self.receiver_distances) == len(self.receiver_azimuths) 

242 assert self.earthmodel_1d is not None 

243 

244 d = self.__dict__.copy() 

245 

246 # fixing these switches here to reduce the amount of wrapper code 

247 d['sw_distance_unit'] = 1 # always give distances in [km] 

248 d['sw_t_reduce'] = 1 # time reduction always as velocity [km/s] 

249 d['sw_equidistant'] = 0 # always give all distances and azimuths 

250 d['sw_irregular_azimuths'] = 1 

251 

252 d['n_distances'] = len(self.receiver_distances) 

253 d['str_distances'] = str_float_vals(self.receiver_distances) 

254 d['str_azimuths'] = str_float_vals(self.receiver_azimuths) 

255 

256 model_str, nlines, ref_depth = cake_model_to_config(self.earthmodel_1d) 

257 d['n_model_lines'] = nlines 

258 d['model_lines'] = model_str 

259 

260 if self.earthmodel_receiver_1d: 

261 model_str, nlines, ref_depth2 = cake_model_to_config( 

262 self.earthmodel_receiver_1d) 

263 assert ref_depth == ref_depth2 

264 else: 

265 model_str = "# no receiver side model" 

266 nlines = 0 

267 

268 d['n_model_receiver_lines'] = nlines 

269 d['model_receiver_lines'] = model_str 

270 

271 d['str_slowness_window'] = str_float_vals(self.slowness_window) 

272 if self.propagation_filters and ref_depth != 0.0: 

273 raise QSeisError( 

274 'Earth model must start with zero depth if ' 

275 'propagation_filters are set.') 

276 

277 d['n_depth_ranges'], d['str_depth_ranges'] = \ 

278 aggregate(self.propagation_filters) 

279 

280 if self.wavelet_type == 0: # user wavelet 

281 d['str_w_samples'] = '\n' \ 

282 + '%i\n' % len(self.user_wavelet_samples) \ 

283 + str_float_vals(self.user_wavelet_samples) 

284 else: 

285 d['str_w_samples'] = '' 

286 

287 if self.receiver_filter: 

288 d['str_receiver_filter'] = self.receiver_filter.string_for_config( 

289 self.qseis_version) 

290 else: 

291 if self.qseis_version == '2006a': 

292 d['str_receiver_filter'] = '(1.0,0.0)\n0\n#\n0' 

293 else: 

294 d['str_receiver_filter'] = '1.0\n0\n#\n0' 

295 

296 d['str_gf_sw_source_types'] = str_int_vals(self.gf_sw_source_types) 

297 d['str_gf_filenames'] = str_str_vals(self.gf_filenames) 

298 

299 if self.source_mech: 

300 d['str_source'] = '%s \'%s\'' % ( 

301 self.source_mech.string_for_config(), 

302 self.seismogram_filename) 

303 else: 

304 d['str_source'] = '0' 

305 

306 d['source_depth_rel'] = d['source_depth'] - ref_depth / km 

307 d['receiver_depth_rel'] = d['receiver_depth'] - ref_depth / km 

308 

309 template = '''# autogenerated QSEIS input by qseis.py 

310# 

311# This is the input file of FORTRAN77 program "qseis06" for calculation of 

312# synthetic seismograms based on a layered halfspace earth model. 

313# 

314# by 

315# Rongjiang Wang <wang@gfz-potsdam.de> 

316# GeoForschungsZentrum Potsdam 

317# Telegrafenberg, D-14473 Potsdam, Germany 

318# 

319# Last modified: Potsdam, Nov., 2006 

320# 

321# = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = 

322# If not specified, SI Unit System is used overall! 

323# 

324# Coordinate systems: 

325# cylindrical (z,r,t) with z = downward, 

326# r = from source outward, 

327# t = azmuth angle from north to east; 

328# cartesian (x,y,z) with x = north, 

329# y = east, 

330# z = downward; 

331# = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = = 

332# 

333# SOURCE PARAMETERS 

334# ================= 

335# 1. source depth [km] 

336#------------------------------------------------------------------------------ 

337 %(source_depth_rel)e |dble: source_depth; 

338#------------------------------------------------------------------------------ 

339# 

340# RECEIVER PARAMETERS 

341# =================== 

342# 1. receiver depth [km] 

343# 2. switch for distance sampling role (1/0 = equidistant/irregular); switch 

344# for unit used (1/0 = km/deg) 

345# 3. number of distance samples 

346# 4. if equidistant, then start and end trace distance (> 0); else distance 

347# list (please order the receiver distances from small to large) 

348# 5. (reduced) time begin [sec] & length of time window [sec], number of time 

349# samples (<= 2*nfmax in qsglobal.h) 

350# 6. switch for unit of the following time reduction parameter: 1 = velocity 

351# [km/sec], 0 = slowness [sec/deg]; time reduction parameter 

352#------------------------------------------------------------------------------ 

353 %(receiver_depth_rel)e |dble: receiver_depth; 

354 %(sw_equidistant)i %(sw_distance_unit)i |int: sw_equidistant, sw_d_unit; 

355 %(n_distances)i |int: no_distances; 

356 %(str_distances)s |dble: d_1,d_n; or d_1,d_2, ...(no comments in between!); 

357 %(time_start)e %(time_window)e %(nsamples)i |dble: t_start,t_window; int: no_t_samples; 

358 %(sw_t_reduce)i %(time_reduction_velocity)e |int: sw_t_reduce; dble: t_reduce; 

359#------------------------------------------------------------------------------ 

360# 

361# WAVENUMBER INTEGRATION PARAMETERS 

362# ================================= 

363# 1. select slowness integration algorithm (0 = suggested for full wave-field 

364# modelling; 1 or 2 = suggested when using a slowness window with narrow 

365# taper range - a technique for suppressing space-domain aliasing); 

366# 2. 4 parameters for low and high slowness (Note 1) cut-offs [s/km] with 

367# tapering: 0 < slw1 < slw2 defining cosine taper at the lower end, and 0 < 

368# slw3 < slw4 defining the cosine taper at the higher end. default values 

369# will be used in case of inconsistent input of the cut-offs (possibly with 

370# much more computational effort); 

371# 3. parameter for sampling rate of the wavenumber integration (1 = sampled 

372# with the spatial Nyquist frequency, 2 = sampled with twice higher than 

373# the Nyquist, and so on: the larger this parameter, the smaller the space- 

374# domain aliasing effect, but also the more computation effort); 

375# 4. the factor for suppressing time domain aliasing (> 0 and <= 1) (Note 2). 

376#------------------------------------------------------------------------------ 

377 %(sw_algorithm)i |int: sw_algorithm; 

378 %(str_slowness_window)s |dble: slw(1-4); 

379 %(wavenumber_sampling)e |dble: sample_rate; 

380 %(aliasing_suppression_factor)e |dble: supp_factor; 

381#------------------------------------------------------------------------------ 

382# 

383# OPTIONS FOR PARTIAL SOLUTIONS 

384# (only applied to the source-site structure) 

385# =========================================== 

386# 

387# 1. switch for filtering free surface effects (0 = with free surface, i.e., 

388# do not select this filter; 1 = without free surface; 2 = without free 

389# surface but with correction on amplitude and wave form. Note switch 2 

390# can only be used for receivers at the surface) 

391# 2. switch for filtering waves with a shallow penetration depth (concerning 

392# their whole trace from source to receiver), penetration depth limit [km] 

393# 

394# if this option is selected, waves whose travel path never exceeds the 

395# given depth limit will be filtered ("seismic nuting"). the condition for 

396# selecting this filter is that the given shallow path depth limit should 

397# be larger than both source and receiver depth. 

398# 

399# 3. number of depth ranges where the following selected up/down-sp2oing P or 

400# SV waves should be filtered 

401# 4. the 1. depth range: upper and lower depth [km], switch for filtering P 

402# or SV wave in this depth range: 

403# 

404# switch no: 1 2 3 4 other 

405# filtered phase: P(up) P(down) SV(up) SV(down) Error 

406# 

407# 5. the 2. ... 

408# 

409# The partial solution options are useful tools to increase the numerical 

410# significance of desired wave phases. Especially when the desired phases 

411# are smaller than the undesired phases, these options should be selected 

412# and carefully combined. 

413#------------------------------------------------------------------------------ 

414 %(filter_surface_effects)i |int: isurf; 

415 %(filter_shallow_paths)i %(filter_shallow_paths_depth)e |int: sw_path_filter; dble:shallow_depth_limit; 

416 %(n_depth_ranges)i %(str_depth_ranges)s 

417#------------------------------------------------------------------------------ 

418# 

419# SOURCE TIME FUNCTION (WAVELET) PARAMETERS (Note 3) 

420# ================================================== 

421# 1. wavelet duration [unit = time sample rather than sec!], that is about 

422# equal to the half-amplitude cut-off period of the wavelet (> 0. if <= 0, 

423# then default value = 2 time samples will be used), and switch for the 

424# wavelet form (0 = user's own wavelet; 1 = default wavelet: normalized 

425# square half-sinusoid for simulating a physical delta impulse; 2 = tapered 

426# Heaviside wavelet, i.e. integral of wavelet 1) 

427# 2. IF user's own wavelet is selected, then number of the wavelet time samples 

428# (<= 1024), and followed by 

429# 3. equidistant wavelet time samples 

430# 4 ...(continue) (! no comment lines allowed between the time sample list!) 

431# IF default, delete line 2, 3, 4 ... or comment them out! 

432#------------------------------------------------------------------------------ 

433 %(wavelet_duration_samples)e %(wavelet_type)i%(str_w_samples)s 

434#------------------------------------------------------------------------------ 

435# 

436# FILTER PARAMETERS OF RECEIVERS (SEISMOMETERS OR HYDROPHONES) 

437# ============================================================ 

438# 1. constant coefficient (normalization factor) 

439# 2. number of roots (<= nrootmax in qsglobal.h) 

440# 3. list of the root positions in the complex format (Re,Im). If no roots, 

441# comment out this line 

442# 4. number of poles (<= npolemax in qsglobal.h) 

443# 5. list of the pole positions in the complex format (Re,Im). If no poles, 

444# comment out this line 

445#------------------------------------------------------------------------------ 

446 %(str_receiver_filter)s 

447#------------------------------------------------------------------------------ 

448# 

449# OUTPUT FILES FOR GREEN'S FUNCTIONS (Note 4) 

450# =========================================== 

451# 1. selections of source types (yes/no = 1/0) 

452# 2. file names of Green's functions (please give the names without extensions, 

453# which will be appended by the program automatically: *.tz, *.tr, *.tt 

454# and *.tv are for the vertical, radial, tangential, and volume change (for 

455# hydrophones) components, respectively) 

456#------------------------------------------------------------------------------ 

457# explosion strike-slip dip-slip clvd single_f_v single_f_h 

458#------------------------------------------------------------------------------ 

459 %(str_gf_sw_source_types)s 

460 %(str_gf_filenames)s 

461#------------------------------------------------------------------------------ 

462# OUTPUT FILES FOR AN ARBITRARY POINT DISLOCATION SOURCE 

463# (for applications to earthquakes) 

464# ====================================================== 

465# 1. selection (0 = not selected; 1 or 2 = selected), if (selection = 1), then 

466# the 6 moment tensor elements [N*m]: Mxx, Myy, Mzz, Mxy, Myz, Mzx (x is 

467# northward, y is eastward and z is downard); else if (selection = 2), then 

468# Mis [N*m] = isotropic moment part = (MT+MN+MP)/3, Mcl = CLVD moment part 

469# = (2/3)(MT+MP-2*MN), Mdc = double-couple moment part = MT-MN, Strike [deg], 

470# Dip [deg] and Rake [deg]. 

471# 

472# Note: to use this option, the Green's functions above should be computed 

473# (selection = 1) if they do not exist already. 

474# 

475# north(x) 

476# / 

477# /\ strike 

478# *-----------------------> east(y) 

479# |\ \ 

480# |-\ \ 

481# | \ fault plane \ 

482# |90 \ \ 

483# |-dip\ \ 

484# | \ \ 

485# | \ \ 

486# downward(z) \-----------------------\\ 

487# 

488# 2. switch for azimuth distribution of the stations (0 = uniform azimuth, 

489# else = irregular azimuth angles) 

490# 3. list of the azimuth angles [deg] for all stations given above (if the 

491# uniform azimuth is selected, then only one azimuth angle is required) 

492# 

493#------------------------------------------------------------------------------ 

494# Mis Mcl Mdc Strike Dip Rake File 

495#------------------------------------------------------------------------------ 

496# 2 0.00 1.00 6.0E+19 120.0 30.0 25.0 'seis' 

497#------------------------------------------------------------------------------ 

498# Mxx Myy Mzz Mxy Myz Mzx File 

499#------------------------------------------------------------------------------ 

500%(str_source)s 

501%(sw_irregular_azimuths)i 

502%(str_azimuths)s 

503#------------------------------------------------------------------------------ 

504# 

505# GLOBAL MODEL PARAMETERS (Note 5) 

506# ================================ 

507# 1. switch for flat-earth-transform 

508# 2. gradient resolution [%%] of vp, vs, and ro (density), if <= 0, then default 

509# values (depending on wave length at cut-off frequency) will be used 

510#------------------------------------------------------------------------------ 

511 %(sw_flat_earth_transform)i |int: sw_flat_earth_transform; 

512 %(gradient_resolution_vp)e %(gradient_resolution_vs)e %(gradient_resolution_density)e |dble: vp_res, vs_res, ro_res; 

513#------------------------------------------------------------------------------ 

514# 

515# LAYERED EARTH MODEL 

516# (SHALLOW SOURCE + UNIFORM DEEP SOURCE/RECEIVER STRUCTURE) 

517# ========================================================= 

518# 1. number of data lines of the layered model (source site) 

519#------------------------------------------------------------------------------ 

520 %(n_model_lines)i |int: no_model_lines; 

521#------------------------------------------------------------------------------ 

522# 

523# MULTILAYERED MODEL PARAMETERS (source site) 

524# =========================================== 

525# no depth[km] vp[km/s] vs[km/s] ro[g/cm^3] qp qs 

526#------------------------------------------------------------------------------ 

527%(model_lines)s 

528#------------------------------------------------------------------------------ 

529# 

530# LAYERED EARTH MODEL 

531# (ONLY THE SHALLOW RECEIVER STRUCTURE) 

532# ===================================== 

533# 1. number of data lines of the layered model 

534# 

535# Note: if the number = 0, then the receiver site is the same as the 

536# source site, else different receiver-site structure is considered. 

537# please be sure that the lowest interface of the receiver-site 

538# structure given given below can be found within the source-site 

539# structure, too. 

540# 

541#------------------------------------------------------------------------------ 

542 %(n_model_receiver_lines)i |int: no_model_lines; 

543#------------------------------------------------------------------------------ 

544# 

545# MULTILAYERED MODEL PARAMETERS (shallow receiver-site structure) 

546# =============================================================== 

547# no depth[km] vp[km/s] vs[km/s] ro[g/cm^3] qp qs 

548#------------------------------------------------------------------------------ 

549%(model_receiver_lines)s 

550#---------------------------------end of all inputs---------------------------- 

551 

552 

553Note 1: 

554 

555The slowness is defined by inverse value of apparent wave velocity = sin(i)/v 

556with i = incident angle and v = true wave velocity. 

557 

558Note 2: 

559 

560The suppression of the time domain aliasing is achieved by using the complex 

561frequency technique. The suppression factor should be a value between 0 and 1. 

562If this factor is set to 0.1, for example, the aliasing phase at the reduced 

563time begin is suppressed to 10%%. 

564 

565Note 3: 

566 

567The default basic wavelet function (option 1) is (2/tau)*sin^2(pi*t/tau), 

568for 0 < t < tau, simulating physical delta impuls. Its half-amplitude cut-off 

569frequency is 1/tau. To avoid high-frequency noise, tau should not be smaller 

570than 4-5 time samples. 

571 

572Note 4: 

573 

574 Double-Couple m11/ m22/ m33/ m12/ m23/ m31 Azimuth_Factor_(tz,tr,tv)/(tt) 

575 ============================================================================ 

576 explosion 1.0/ 1.0/ 1.0/ -- / -- / -- 1.0 / 0.0 

577 strike-slip -- / -- / -- / 1.0/ -- / -- sin(2*azi) / cos(2*azi) 

578 1.0/-1.0/ -- / -- / -- / -- cos(2*azi) / -sin(2*azi) 

579 dip-slip -- / -- / -- / -- / -- / 1.0 cos(azi) / sin(azi) 

580 -- / -- / -- / -- / 1.0/ -- sin(azi) / -cos(azi) 

581 clvd -0.5/-0.5/ 1.0/ -- / -- / -- 1.0 / 0.0 

582 ============================================================================ 

583 Single-Force fx / fy / fz Azimuth_Factor_(tz,tr,tv)/(tt) 

584 ============================================================================ 

585 fz -- / -- / 1.0 1.0 / 0.0 

586 fx 1.0/ -- / -- cos(azi) / sin(azi) 

587 fy -- / 1.0/ -- sin(azi) / -cos(azi) 

588 ============================================================================ 

589 

590Note 5: 

591 

592Layers with a constant gradient will be discretized with a number of homogeneous 

593sublayers. The gradient resolutions are then used to determine the maximum 

594allowed thickness of the sublayers. If the resolutions of Vp, Vs and Rho 

595(density) require different thicknesses, the smallest is first chosen. If this 

596is even smaller than 1%% of the characteristic wavelength, then the latter is 

597taken finally for the sublayer thickness. 

598''' # noqa 

599 

600 return (template % d).encode('ascii') 

601 

602 

603class QSeisError(gf.store.StoreError): 

604 pass 

605 

606 

607class Interrupted(gf.store.StoreError): 

608 def __str__(self): 

609 return 'Interrupted.' 

610 

611 

612class QSeisRunner(object): 

613 

614 def __init__(self, tmp=None, keep_tmp=False): 

615 self.tempdir = mkdtemp(prefix='qseisrun-', dir=tmp) 

616 self.keep_tmp = keep_tmp 

617 self.config = None 

618 

619 def run(self, config): 

620 self.config = config 

621 

622 input_fn = pjoin(self.tempdir, 'input') 

623 

624 with open(input_fn, 'wb') as f: 

625 input_str = config.string_for_config() 

626 logger.debug('===== begin qseis input =====\n' 

627 '%s===== end qseis input =====' % input_str.decode()) 

628 f.write(input_str) 

629 

630 program = program_bins['qseis.%s' % config.qseis_version] 

631 

632 old_wd = os.getcwd() 

633 

634 os.chdir(self.tempdir) 

635 

636 interrupted = [] 

637 

638 def signal_handler(signum, frame): 

639 os.kill(proc.pid, signal.SIGTERM) 

640 interrupted.append(True) 

641 

642 original = signal.signal(signal.SIGINT, signal_handler) 

643 try: 

644 try: 

645 proc = Popen(program, stdin=PIPE, stdout=PIPE, stderr=PIPE) 

646 

647 except OSError: 

648 os.chdir(old_wd) 

649 raise QSeisError( 

650 '''could not start qseis executable: "%s" 

651Available fomosto backends and download links to the modelling codes are listed 

652on 

653 

654 https://pyrocko.org/docs/current/apps/fomosto/backends.html 

655 

656''' % program) 

657 

658 (output_str, error_str) = proc.communicate(b'input\n') 

659 

660 finally: 

661 signal.signal(signal.SIGINT, original) 

662 

663 if interrupted: 

664 raise KeyboardInterrupt() 

665 

666 logger.debug('===== begin qseis output =====\n' 

667 '%s===== end qseis output =====' % output_str.decode()) 

668 

669 errmess = [] 

670 if proc.returncode != 0: 

671 errmess.append( 

672 'qseis had a non-zero exit state: %i' % proc.returncode) 

673 

674 if error_str: 

675 logger.warning( 

676 'qseis emitted something via stderr:\n\n%s' 

677 % error_str.decode()) 

678 

679 # errmess.append('qseis emitted something via stderr') 

680 

681 if output_str.lower().find(b'error') != -1: 

682 errmess.append("the string 'error' appeared in qseis output") 

683 

684 if errmess: 

685 self.keep_tmp = True 

686 

687 os.chdir(old_wd) 

688 raise QSeisError(''' 

689===== begin qseis input ===== 

690%s===== end qseis input ===== 

691===== begin qseis output ===== 

692%s===== end qseis output ===== 

693===== begin qseis error ===== 

694%s===== end qseis error ===== 

695%s 

696qseis has been invoked as "%s" 

697in the directory %s'''.lstrip() % ( 

698 input_str.decode(), output_str.decode(), error_str.decode(), 

699 '\n'.join(errmess), program, self.tempdir)) 

700 

701 self.qseis_output = output_str 

702 self.qseis_error = error_str 

703 

704 os.chdir(old_wd) 

705 

706 def get_traces(self, which='seis'): 

707 

708 if which == 'seis': 

709 fns = self.config.get_output_filenames(self.tempdir) 

710 components = qseis_components 

711 

712 elif which == 'gf': 

713 fns = self.config.get_output_filenames_gf(self.tempdir) 

714 components = [ 

715 fn+'.t'+c 

716 for fn in self.config.gf_filenames for c in qseis_components] 

717 else: 

718 raise Exception( 

719 'get_traces: which argument should be "seis" or "gf"') 

720 

721 traces = [] 

722 distances = self.config.receiver_distances 

723 azimuths = self.config.receiver_azimuths 

724 for comp, fn in zip(components, fns): 

725 if not os.path.exists(fn): 

726 continue 

727 

728 data = num.loadtxt(fn, skiprows=1, dtype=float) 

729 nsamples, ntraces = data.shape 

730 ntraces -= 1 

731 vred = self.config.time_reduction_velocity 

732 deltat = (data[-1, 0] - data[0, 0])/(nsamples-1) 

733 

734 for itrace, distance, azimuth in zip( 

735 range(ntraces), distances, azimuths): 

736 

737 tmin = self.config.time_start 

738 if vred != 0.0: 

739 tmin += distance / vred 

740 

741 tmin += deltat 

742 tr = trace.Trace( 

743 '', '%04i' % itrace, '', comp, 

744 tmin=tmin, deltat=deltat, ydata=data[:, itrace+1], 

745 meta=dict( 

746 distance=distance*km, 

747 azimuth=azimuth)) 

748 

749 traces.append(tr) 

750 

751 return traces 

752 

753 def __del__(self): 

754 if self.tempdir: 

755 if not self.keep_tmp: 

756 shutil.rmtree(self.tempdir) 

757 self.tempdir = None 

758 else: 

759 logger.warning( 

760 'not removing temporary directory: %s' % self.tempdir) 

761 

762 

763class QSeisGFBuilder(gf.builder.Builder): 

764 def __init__(self, store_dir, step, shared, block_size=None, tmp=None, 

765 force=False): 

766 

767 self.store = gf.store.Store(store_dir, 'w') 

768 

769 if block_size is None: 

770 block_size = (1, 1, 100) 

771 

772 if len(self.store.config.ns) == 2: 

773 block_size = block_size[1:] 

774 

775 gf.builder.Builder.__init__( 

776 self, self.store.config, step, block_size=block_size, force=force) 

777 

778 baseconf = self.store.get_extra('qseis') 

779 

780 conf = QSeisConfigFull(**baseconf.items()) 

781 conf.earthmodel_1d = self.store.config.earthmodel_1d 

782 conf.earthmodel_receiver_1d = self.store.config.earthmodel_receiver_1d 

783 

784 deltat = 1.0/self.gf_config.sample_rate 

785 

786 if 'time_window_min' not in shared: 

787 d = self.store.make_timing_params( 

788 conf.time_region[0], conf.time_region[1], 

789 force=force) 

790 

791 shared['time_window_min'] = d['tlenmax_vred'] 

792 shared['time_start'] = d['tmin_vred'] 

793 shared['time_reduction_velocity'] = d['vred'] / km 

794 

795 time_window_min = shared['time_window_min'] 

796 conf.time_start = shared['time_start'] 

797 

798 conf.time_reduction_velocity = shared['time_reduction_velocity'] 

799 

800 conf.nsamples = nextpow2(int(round(time_window_min / deltat)) + 1) 

801 conf.time_window = (conf.nsamples-1)*deltat 

802 

803 self.qseis_config = conf 

804 

805 self.tmp = tmp 

806 if self.tmp is not None: 

807 util.ensuredir(self.tmp) 

808 

809 def cleanup(self): 

810 self.store.close() 

811 

812 def work_block(self, index): 

813 if len(self.store.config.ns) == 2: 

814 (sz, firstx), (sz, lastx), (ns, nx) = \ 

815 self.get_block_extents(index) 

816 

817 rz = self.store.config.receiver_depth 

818 else: 

819 (rz, sz, firstx), (rz, sz, lastx), (nr, ns, nx) = \ 

820 self.get_block_extents(index) 

821 

822 conf = copy.deepcopy(self.qseis_config) 

823 

824 logger.info('Starting block %i / %i' % 

825 (index+1, self.nblocks)) 

826 

827 conf.source_depth = float(sz/km) 

828 conf.receiver_depth = float(rz/km) 

829 

830 runner = QSeisRunner(tmp=self.tmp) 

831 

832 dx = self.gf_config.distance_delta 

833 

834 distances = num.linspace(firstx, firstx + (nx-1)*dx, nx).tolist() 

835 

836 if distances[-1] < self.gf_config.distance_max: 

837 # add global max distance, because qseis does some adjustments with 

838 # this value 

839 distances.append(self.gf_config.distance_max) 

840 

841 mex = (MomentTensor(m=symmat6(1, 1, 1, 0, 0, 0)), 

842 {'r': (0, +1), 'z': (1, +1)}) 

843 

844 mmt1 = (MomentTensor(m=symmat6(1, 0, 0, 1, 0, 0)), 

845 {'r': (0, +1), 't': (3, +1), 'z': (5, +1)}) 

846 mmt2 = (MomentTensor(m=symmat6(0, 0, 0, 0, 1, 1)), 

847 {'r': (1, +1), 't': (4, +1), 'z': (6, +1)}) 

848 mmt3 = (MomentTensor(m=symmat6(0, 0, 1, 0, 0, 0)), 

849 {'r': (2, +1), 'z': (7, +1)}) 

850 mmt4 = (MomentTensor(m=symmat6(0, 1, 0, 0, 0, 0)), 

851 {'r': (8, +1), 'z': (9, +1)}) 

852 

853 component_scheme = self.store.config.component_scheme 

854 off = 0 

855 if component_scheme == 'elastic8': 

856 off = 8 

857 elif component_scheme == 'elastic10': 

858 off = 10 

859 

860 msf = (None, { 

861 'fz.tr': (off+0, +1), 

862 'fh.tr': (off+1, +1), 

863 'fh.tt': (off+2, -1), 

864 'fz.tz': (off+3, +1), 

865 'fh.tz': (off+4, +1)}) 

866 if component_scheme == 'elastic2': 

867 gfsneeded = (1, 0, 0, 0, 0, 0) 

868 gfmapping = [mex] 

869 

870 if component_scheme == 'elastic5': 

871 gfsneeded = (0, 0, 0, 0, 1, 1) 

872 gfmapping = [msf] 

873 

874 elif component_scheme == 'elastic8': 

875 gfsneeded = (1, 1, 1, 1, 0, 0) 

876 gfmapping = [mmt1, mmt2, mmt3] 

877 

878 elif component_scheme == 'elastic10': 

879 gfsneeded = (1, 1, 1, 1, 0, 0) 

880 gfmapping = [mmt1, mmt2, mmt3, mmt4] 

881 

882 elif component_scheme == 'elastic13': 

883 gfsneeded = (1, 1, 1, 1, 1, 1) 

884 gfmapping = [mmt1, mmt2, mmt3, msf] 

885 

886 elif component_scheme == 'elastic15': 

887 gfsneeded = (1, 1, 1, 1, 1, 1) 

888 gfmapping = [mmt1, mmt2, mmt3, mmt4, msf] 

889 

890 conf.gf_sw_source_types = gfsneeded 

891 conf.receiver_distances = [d/km for d in distances] 

892 conf.receiver_azimuths = [0.0] * len(distances) 

893 

894 for mt, gfmap in gfmapping: 

895 if mt: 

896 m = mt.m() 

897 f = float 

898 conf.source_mech = QSeisSourceMechMT( 

899 mnn=f(m[0, 0]), mee=f(m[1, 1]), mdd=f(m[2, 2]), 

900 mne=f(m[0, 1]), mnd=f(m[0, 2]), med=f(m[1, 2])) 

901 else: 

902 conf.source_mech = None 

903 

904 if any(conf.gf_sw_source_types) or conf.source_mech is not None: 

905 runner.run(conf) 

906 

907 if any(c in gfmap for c in qseis_components): 

908 rawtraces = runner.get_traces('seis') 

909 else: 

910 rawtraces = runner.get_traces('gf') 

911 

912 interrupted = [] 

913 

914 def signal_handler(signum, frame): 

915 interrupted.append(True) 

916 

917 original = signal.signal(signal.SIGINT, signal_handler) 

918 self.store.lock() 

919 try: 

920 for itr, tr in enumerate(rawtraces): 

921 if tr.channel not in gfmap: 

922 continue 

923 

924 x = tr.meta['distance'] 

925 if x > firstx + (nx-1)*dx: 

926 continue 

927 

928 ig, factor = gfmap[tr.channel] 

929 

930 if len(self.store.config.ns) == 2: 

931 args = (sz, x, ig) 

932 else: 

933 args = (rz, sz, x, ig) 

934 

935 if conf.cut: 

936 tmin = self.store.t(conf.cut[0], args[:-1]) 

937 tmax = self.store.t(conf.cut[1], args[:-1]) 

938 

939 if None in (tmin, tmax): 

940 self.warn( 

941 'Failed cutting {} traces. ' + 

942 'Failed to determine time window') 

943 continue 

944 

945 tr.chop(tmin, tmax) 

946 

947 tmin = tr.tmin 

948 tmax = tr.tmax 

949 

950 if conf.fade: 

951 ta, tb, tc, td = [ 

952 self.store.t(v, args[:-1]) for v in conf.fade] 

953 

954 if None in (ta, tb, tc, td): 

955 continue 

956 

957 if not (ta <= tb and tb <= tc and tc <= td): 

958 raise QSeisError( 

959 'invalid fade configuration ' 

960 '(it should be (ta <= tb <= tc <= td) but ' 

961 'ta=%g, tb=%g, tc=%g, td=%g)' % ( 

962 ta, tb, tc, td)) 

963 

964 t = tr.get_xdata() 

965 fin = num.interp(t, [ta, tb], [0., 1.]) 

966 fout = num.interp(t, [tc, td], [1., 0.]) 

967 anti_fin = 1. - fin 

968 anti_fout = 1. - fout 

969 

970 y = tr.ydata 

971 

972 sum_anti_fin = num.sum(anti_fin) 

973 sum_anti_fout = num.sum(anti_fout) 

974 

975 if sum_anti_fin != 0.0: 

976 yin = num.sum(anti_fin*y) / sum_anti_fin 

977 else: 

978 yin = 0.0 

979 

980 if sum_anti_fout != 0.0: 

981 yout = num.sum(anti_fout*y) / sum_anti_fout 

982 else: 

983 yout = 0.0 

984 

985 y2 = anti_fin*yin + fin*fout*y + anti_fout*yout 

986 

987 if conf.relevel_with_fade_in: 

988 y2 -= yin 

989 

990 tr.set_ydata(y2) 

991 

992 gf_tr = gf.store.GFTrace.from_trace(tr) 

993 gf_tr.data *= factor 

994 

995 try: 

996 self.store.put(args, gf_tr) 

997 except gf.store.DuplicateInsert: 

998 self.warn('{} insertions_skipped (duplicates)') 

999 

1000 finally: 

1001 self.log_warnings(index+1, logger) 

1002 self.store.unlock() 

1003 signal.signal(signal.SIGINT, original) 

1004 

1005 if interrupted: 

1006 raise KeyboardInterrupt() 

1007 

1008 conf.gf_sw_source_types = (0, 0, 0, 0, 0, 0) 

1009 

1010 logger.info('Done with block %i / %i' % 

1011 (index+1, self.nblocks)) 

1012 

1013 

1014def init(store_dir, variant, config_params=None): 

1015 if variant is None: 

1016 variant = '2006' 

1017 

1018 if variant not in ('2006', '2006a'): 

1019 raise gf.store.StoreError('unsupported variant: %s' % variant) 

1020 

1021 modelling_code_id = 'qseis.%s' % variant 

1022 

1023 qseis = QSeisConfig(qseis_version=variant) 

1024 qseis.time_region = ( 

1025 gf.meta.Timing('begin-50'), 

1026 gf.meta.Timing('end+100')) 

1027 

1028 qseis.cut = ( 

1029 gf.meta.Timing('begin-50'), 

1030 gf.meta.Timing('end+100')) 

1031 

1032 qseis.wavelet_duration_samples = 0.001 

1033 qseis.sw_flat_earth_transform = 1 

1034 

1035 store_id = os.path.basename(os.path.realpath(store_dir)) 

1036 

1037 d = dict( 

1038 id=store_id, 

1039 ncomponents=10, 

1040 sample_rate=0.2, 

1041 receiver_depth=0*km, 

1042 source_depth_min=10*km, 

1043 source_depth_max=20*km, 

1044 source_depth_delta=10*km, 

1045 distance_min=100*km, 

1046 distance_max=1000*km, 

1047 distance_delta=10*km, 

1048 earthmodel_1d=cake.load_model().extract(depth_max='cmb'), 

1049 modelling_code_id=modelling_code_id, 

1050 tabulated_phases=[ 

1051 gf.meta.TPDef( 

1052 id='begin', 

1053 definition='p,P,p\\,P\\,Pv_(cmb)p'), 

1054 gf.meta.TPDef( 

1055 id='end', 

1056 definition='2.5'), 

1057 gf.meta.TPDef( 

1058 id='P', 

1059 definition='!P'), 

1060 gf.meta.TPDef( 

1061 id='S', 

1062 definition='!S'), 

1063 gf.meta.TPDef( 

1064 id='p', 

1065 definition='!p'), 

1066 gf.meta.TPDef( 

1067 id='s', 

1068 definition='!s')]) 

1069 

1070 if config_params is not None: 

1071 d.update(config_params) 

1072 

1073 config = gf.meta.ConfigTypeA(**d) 

1074 

1075 config.validate() 

1076 return gf.store.Store.create_editables( 

1077 store_dir, config=config, extra={'qseis': qseis}) 

1078 

1079 

1080def build(store_dir, force=False, nworkers=None, continue_=False, step=None, 

1081 iblock=None): 

1082 

1083 return QSeisGFBuilder.build( 

1084 store_dir, force=force, nworkers=nworkers, continue_=continue_, 

1085 step=step, iblock=iblock)