-
Notifications
You must be signed in to change notification settings - Fork 33
/
EDKSmp.py
828 lines (664 loc) · 27.4 KB
/
EDKSmp.py
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
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
'''
A bunch of routines to handle EDKS
Written by F. Ortega in 2010.
Modified by R. Jolivet in 2014.
Modified by R. Jolivet in 2017 (multiprocessing added for point dropping)
'''
# Externals
import os
import struct
import sys
import numpy as np
import copy
import multiprocessing as mp
# Scipy
from scipy.io import FortranFile
import scipy.interpolate as sciint
# Initialize a class to allow multiprocessing for EDKS interpolation in Python
class interpolator(mp.Process):
'''
Multiprocessing class runing the edks interpolation.
This class requires one to build the interpolator in advance.
Args:
* interpolators : List of interpolators
* queue : Instance of mp.Queue
* depths : depths (first dimnesion of the interpolators)
* distas : distances (second dimension of the interpolators)
* istart : starting point
* iend : ending point
Returns:
* None
'''
# ----------------------------------------------------------------------
# Initialize
def __init__(self, interpolators, queue, depths, distas, istart, iend):
# Save things
self.interpolators = interpolators
self.depths = depths
self.distas = distas
self.istart = istart
self.iend = iend
# Save the queue
self.queue = queue
# Initialize the process
super(interpolator, self).__init__()
# All done
return
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
# Run method
def run(self):
'''
Run the interpolation
'''
# Interpolate
values = []
for inter in self.interpolators:
values.append(inter(np.vstack((self.depths[self.istart:self.iend],
self.distas[self.istart:self.iend])).T))
# Save start/end
values.append((self.istart, self.iend))
# Store output
self.queue.put(values)
# All done
return
# ----------------------------------------------------------------------
# Initialize a class to allow multiprocessing to drop points
class pointdropper(mp.Process):
'''
Initialize the multiprocessing class to run the point dropper.
This class drops point sources in the triangular or rectangular mesh.
Args:
* fault : Instance of Fault.py
* queue : Instance of mp.Queue
* charArea : Characteristic area of the subfaults
* istart : Index of the first patch to deal with
* iend : Index of the last pacth to deal with
Returns:
* None
'''
# ----------------------------------------------------------------------
# Initialize
def __init__(self, fault, queue, charArea, istart, iend):
# Save the fault
self.fault = copy.deepcopy(fault)
self.charArea = charArea
self.istart = istart
self.iend = iend
#print(istart, iend)
# Save the queue
self.queue = queue
# Initialize the Process
super(pointdropper, self).__init__()
# All done
return
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
# Run routine needed by multiprocessing
def run(self):
'''
Run the subpatch construction
'''
# Create lists
Ids, Xs, Ys, Zs, Strike, Dip, Area = [], [], [], [], [], [], []
allSplitted = []
# Iterate overthe patches
for i in range(self.istart, self.iend):
# Get patch
patch = self.fault.patch[i]
# Check if the Area is bigger than the target
area = self.fault.patchArea(patch)
if area>self.charArea[i]:
keepGoing = True
tobeSplitted = [patch]
splittedPatches = []
else:
keepGoing = False
print('Be carefull, patch {} has not been refined into point sources'.format(self.fault.getindex(patch)))
print('Possible causes: Area = {}, Nodes = {}'.format(area, patch))
tobeSplitted = []
splittedPatches = [patch]
# Iterate
while keepGoing:
# Take a patch
p = tobeSplitted.pop()
# Split into 4 patches
Splitted = self.fault.splitPatch(p)
# Check the area
for splitted in Splitted:
# get area
area = self.fault.patchArea(splitted)
# check
if area<self.charArea[i]:
splittedPatches.append(splitted)
else:
tobeSplitted.append(splitted)
# Do we continue?
if len(tobeSplitted)==0:
keepGoing = False
# Do we have a limit
if hasattr(self.fault, 'maximumSources'):
if len(splittedPatches)>=self.fault.maximumSources:
keepGoing = False
# When all done get their centers
geometry = [self.fault.getpatchgeometry(p, center=True)[:3] for p in splittedPatches]
x, y, z = zip(*geometry)
strike, dip = self.fault.getpatchgeometry(patch)[5:7]
strike = np.ones((len(x),))*strike
strike = strike.tolist()
dip = np.ones((len(x),))*dip
dip = dip.tolist()
areas = [self.fault.patchArea(p) for p in splittedPatches]
ids = np.ones((len(x),))*(i)
ids = ids.astype(int).tolist()
# Save
Ids += ids
Xs += x
Ys += y
Zs += z
Strike += strike
Dip += dip
Area += areas
allSplitted += splittedPatches
# Put in the Queue
self.queue.put([Ids, Xs, Ys, Zs, Strike, Dip, Area, allSplitted])
# all done
return
# ----------------------------------------------------------------------
# end of pointdropper
# ----------------------------------------------------------------------
def dropSourcesInPatches(fault, verbose=False, returnSplittedPatches=False):
'''
From a fault object, returns sources to be given to sum_layered_sub.
The number of sources is determined by the spacing provided in fault.
Args:
* fault : instance of Fault (Rectangular or Triangular).
* verbose : Talk to me
* returnSplittedPactches : Returns a triangularPatches object with the splitted
patches.
Return:
* Ids : Id of the subpatches
* Xs : UTM x-coordinate of the subpatches (km)
* Ys : UTM y-coordinate of the subpatches (km)
* Zs : UTM z-coordinate of the subpatches (km)
* Strikes : Strike angles of the subpatches (rad)
* Dips : Dip angles of the subpatches (rad)
* Areas : Area of the subpatches (km^2)
if returnSplittedPatches:
* splitFault : Fault object with the subpatches
'''
# Create lists
Id, X, Y, Z, Strike, Dip, Area = [], [], [], [], [], [], []
Splitted = []
# Check
if (not hasattr(fault, 'sourceSpacing')) and (not hasattr(fault, 'sourceNumber')) and (not hasattr(fault, 'sourceArea')):
print('EDKS: Need to provide area, spacing or number of sources...')
sys.exit(1)
if hasattr(fault, 'sourceSpacing') and hasattr(fault, 'sourceNumber') and hasattr(fault, 'sourceArea'):
print('EDKS: Please delete sourceSpacing, sourceNumber or sourceArea...')
print('EDKS: I do not judge... You decide...')
sys.exit(1)
# show me
if verbose:
print('Dropping point sources')
# Spacing
if hasattr(fault, 'sourceArea'):
area = fault.sourceArea
charArea = np.ones((len(fault.patch),))*area
if hasattr(fault, 'sourceSpacing'):
spacing = fault.sourceSpacing
if fault.patchType == 'rectangle':
charArea = np.ones((len(fault.patch),))*spacing**2
elif fault.patchType in ('triangle', 'triangletent'):
charArea = np.ones((len(fault.patch),))*spacing**2/2.
if hasattr(fault, 'sourceNumber'):
number = fault.sourceNumber
fault.computeArea()
charArea = np.array(fault.area)/float(number)
# Create a queue
output = mp.Queue()
# how many workers
try:
nworkers = int(os.environ['OMP_NUM_THREADS'])
except:
nworkers = mp.cpu_count()
# how many patches
npatches = len(fault.patch)
# Create them
workers = [pointdropper(fault, output, charArea,
int(np.floor(i*npatches/nworkers)),
int(np.floor((i+1)*npatches/nworkers))) for i in range(nworkers)]
workers[-1].iend = npatches
# Start them
for w in range(nworkers): workers[w].start()
# I don't understand why this guy does not work...
#for w in range(nworkers): workers[w].join()
# Get things from the queue
for i in range(nworkers):
ids, xs, ys, zs, strike, dip, area, splitted = output.get()
Id.extend(ids)
X.extend(xs)
Y.extend(ys)
Z.extend(zs)
Strike.extend(strike)
Dip.extend(dip)
Area.extend(area)
Splitted.extend(splitted)
# Make arrays
isort = np.argsort(Id)
Ids = np.array([Id[i] for i in isort])
Xs = np.array([X[i] for i in isort])
Ys = np.array([Y[i] for i in isort])
Zs = np.array([Z[i] for i in isort])
Strikes = np.array([Strike[i] for i in isort])
Dips = np.array([Dip[i] for i in isort])
Areas = np.array([Area[i] for i in isort])
allSplitted = [Splitted[i] for i in isort]
# All done
if returnSplittedPatches:
return Ids, Xs, Ys, Zs, Strikes, Dips, Areas, allSplitted
else:
return Ids, Xs, Ys, Zs, Strikes, Dips, Areas
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
# Compute the Green's functions for the patches
def sum_layered(xs, ys, zs, strike, dip, rake, slip, width, length,\
npw, npy,\
xr, yr, edks,\
prefix, \
BIN_EDKS = 'EDKS_BIN',
cleanUp=True, verbose=True):
'''
Compute the Green's functions for the given patches
Args:
<-- Sources --> 1-D numpy arrays
* xs : m, east coord to center of fault patch
* ys : m, north coord to center of fault patch
* zs : m,depth coord to center of fault patch (+ down)
* strike : deg, clockwise from north
* dip : deg, 90 is vertical
* rake : deg, 0 left lateral strike slip, 90 up-dip slip
* slip : m, slip in the rake direction
* width : m, width of the patch
* length : m, length of the patch
* npw : integers, number of sources along strike
* npy : integers, number of sources along dip
<-- Receivers --> 1-D numpy arrays
* xr : m, east coordinate of receivers
* yr : m, north coordinate of receivers
<-- Elastic structure -->
* edks : string, full name of edks file, e.g., halfspace.edks
<-- File Naming -->
* prefix : string, prefix for the files generated by sum_layered
Kwargs:
* BIN_EDKS : Environement variable where EDKS executables are.
* cleanUp : Remove the intermediate files
* verbose : Talk to me
Return:
<-- 2D arrays (#receivers, #fault patches) -->
* ux : m, east displacement
* uy : m, west displacement
* uz : m, up displacement (+ up)
'''
# Get executables
BIN_EDKS = os.environ[BIN_EDKS]
# Some initializations
Np = len(xs) # number of sources
nrec = len(xr) # number of receivers
A = length*width # Area of the patches
# Some formats
BIN_FILE_FMT = 'f' # python float = C/C++ float = Fortran 'real*4'
NBYTES_FILE_FMT = 4 # a Fortran (real*4) uses 4 bytes.
# convert sources from center to top edge of fault patch ("sum_layered" needs that)
sind = np.sin( dip * np.pi / 180.0 )
cosd = np.cos( dip * np.pi / 180.0 )
sins = np.sin( strike * np.pi / 180.0 )
coss = np.cos( strike * np.pi / 180.0 )
# displacement in local coordinates (phi, delta)
dZ = (width/2.0) * sind
dD = (width/2.0) * cosd
# rotation to global coordinates
xs = xs - dD * coss
ys = ys + dD * sins
zs = zs - dZ
# Define filenames:
file_rec = prefix + '.rec'
file_pat = prefix + '.pat'
file_dux = prefix + '_ux.dis'
file_duy = prefix + '_uy.dis'
file_duz = prefix + '_uz.dis'
# Clean the file if they exist
cmd = 'rm -f {} {} {} {} {}'.format(file_rec, file_pat, file_dux, file_duy, file_duz)
os.system(cmd)
# write receiver location file (observation points)
temp = [xr, yr]
file = open(file_rec, 'wb')
for k in range(0, nrec):
for i in range(0, len(temp)):
file.write( struct.pack( BIN_FILE_FMT, temp[i][k] ) )
file.close()
# write point sources information
temp = [xs, ys, zs, strike, dip, rake, width, length, slip]
file = open(file_pat, 'wb');
for k in range(0, Np):
for i in range(0, len(temp)):
file.write( struct.pack( BIN_FILE_FMT, temp[i][k] ) )
file.close()
# call sum_layered
if not os.path.exists(os.path.basename(edks)):
os.symlink(edks, os.path.basename(edks))
os.symlink(os.path.join(os.path.dirname(edks), 'hdr.'+os.path.basename(edks)), 'hdr.'+os.path.basename(edks))
removeSymLink = True
else:
removeSymLink = False
cmd = '{}/sum_layered {} {} {} {} {} {}'.format(BIN_EDKS, os.path.basename(edks), prefix, nrec, Np, npw, npy)
if verbose:
print(cmd)
os.system(cmd)
if removeSymLink:
os.unlink(os.path.basename(edks))
os.unlink('hdr.'+os.path.basename(edks))
# read sum_layered output Greens function
# ux
ux = np.fromfile(file_dux, 'f').reshape((nrec, Np), order='F')
# uy
uy = np.fromfile(file_duy, 'f').reshape((nrec, Np), order='F')
# uz
uz = np.fromfile(file_duz, 'f').reshape((nrec, Np), order='F')
# remove IO files.
if cleanUp:
cmd = 'rm -f {} {} {} {} {}'.format(file_rec, file_pat, file_dux, file_duy, file_duz)
os.system(cmd)
# return the GF matrices
return [ux, uy, uz]
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
# A class that interpolates edks Kernels (same as fortran's sum_layered,
# but with more flexibility for the interpolation part)
class interpolateEDKS(object):
'''
A class that will interpolate the EDKS Kernels and produce Green's
functions in a stratified medium. This class will only use point
sources as the summation is done in the fault object.
What goes in this class is a translation of the point source case of
EDKS. We use the case where slip perpendicular to the rake angle is
equal to zero.
Args:
* kernel : EDKS Kernel file (mykernel.edks). One needs to
provide the header file as well (hdr.mykernel.edks)
'''
def __init__(self, kernel, verbose=True):
# Set verbose
self.verbose = verbose
# Set kernel
self.kernel = kernel
# Make sure we start on the right foot
self.interpolationDone = False
# All done
return
def readHeader(self):
'''
Read the EDKS Kernel header file and stores it in {self}
Returns:
* None
'''
# Show me
if self.verbose:
print('Read Kernel Header file hdr.{}'.format(self.kernel))
# Open the header file
fhd = open('hdr.{}'.format(self.kernel), 'r')
# Read things
self.prefix = fhd.readline().split()[0]
self.nlayer = int(fhd.readline().split()[0])
# Layer characteristics
rho = []
alpha = []
beta = []
thickness = []
# Iterate over layers
for i in range(self.nlayer):
line = fhd.readline().split()
rho.append(float(line[0]))
alpha.append(float(line[1]))
beta.append(float(line[2]))
thickness.append(float(line[3]))
# Software date
self.softwareDate = fhd.readline().split()
self.softwareVersion = fhd.readline().split()
self.softwareComments = fhd.readline().split()
# Depths, Distances
depths = fhd.readline().split()
self.depthmin = float(depths[0])
self.depthmax = float(depths[1])
self.ndepth = int(depths[2])
distances = fhd.readline().split()
self.distamin = float(distances[0])
self.distamax = float(distances[1])
self.ndista = int(distances[2])
# Close file
fhd.close()
# All done
return
def readKernel(self):
'''
Read the EDKS Kernel and stores it in {self}
Returns:
* None
'''
# Show me
if self.verbose:
print('Read Kernel file {}'.format(self.kernel))
# Open the file
fedks = FortranFile(self.kernel, 'r')
# Read
kernel = fedks.read_reals(np.float32).reshape((self.ndepth*self.ndista,12))
# Save
self.depths = kernel[:,0]
self.distas = kernel[:,1]
self.zrtdsx = kernel[:,2:]
# CLose file
fedks.close()
def interpolate(self, xs, ys, zs, strike, dip, rake, area, slip, xr, yr, method='linear'):
'''
Interpolate the Green's functions for a given source in (xs, ys, zs) with
a strike, dip and rake and slip parameters and a given receiver (xr, yr)
Args:
* xs, ys, zs : Source location (floats or np.array)
* strike : strike angle (rad)
* dip : dip angle (rad)
* rake : rake angle (rad, 0 left-lateral strike slip, 2pi pure thrust)
* slip : Slip value. The unit of slip will condition the unit of the output displacement
* area : Area of the point source
* xr, yr : Receiver location (floats or np.array)
Kwargs:
* method : Interpolation scheme. Can be linear, nearest or CloughTocher.
Returns:
* G : np.array
'''
# Arrange things
if type(xs) in (float, np.float64, np.float32):
xs = np.array([xs])
ys = np.array([ys])
zs = np.array([zs])
strike = np.array([strike])
dip = np.array([dip])
rake = np.array([rake])
slip = np.array([slip])
area = np.array([area])
if type(xr) in (float, np.float64, np.float32):
xr = np.array([xr])
yr = np.array([yr])
# convert sources from center to top edge of fault patch
sind = np.sin( dip )
cosd = np.cos( dip )
sins = np.sin( strike )
coss = np.cos( strike )
# displacement in local coordinates (phi, delta)
dZ = (np.sqrt(area)/2.0) * sind
dD = (np.sqrt(area)/2.0) * cosd
# rotation to global coordinates
xs = xs - dD * coss
ys = ys + dD * sins
zs = zs - dZ
#Show me
if self.verbose:
print('Interpolate GFs for {} sources and {} receivers'.format(len(xs), len(xr)))
# Get moment (here potency)
M = self.src2mom(slip, area, strike, dip, rake)
# Create an interpolator
self.createInterpolator(method=method)
# Calculate geometry -- dim(r) is (sources, receivers)
if self.verbose:
print('Calculate geometry')
distance, depth, caz, saz, c2az, s2az = self._getGeometry(xs, ys, zs, xr, yr)
if not self.interpolationDone:
# Interpolate
if self.verbose:
print('Interpolate')
# Create holder
self.interpKernels = np.zeros((len(xs)*len(xr), 10))
# Multiprocessing
try:
nworkers = int(os.environ['OMP_NUM_THREADS'])
except:
nworkers = mp.cpu_count()
# Create a queue
output = mp.Queue()
# Create the workers
todo = len(distance.flatten())
workers = [interpolator(self.interpolators, output,
depth.flatten(), distance.flatten(),
int(np.floor(i*todo/nworkers)),
int(np.floor((i+1)*todo/nworkers))) for i in range(nworkers)]
workers[-1].iend = todo
# Start
for w in range(nworkers): workers[w].start()
# Get from the queue
for worker in workers:
values = output.get()
istart,iend = values.pop()
for iv,value in enumerate(values):
self.interpKernels[istart:iend,iv] = value
# Reshape
self.interpKernels = self.interpKernels.reshape((len(xs),len(xr),10))
# reconstruction of the actual displacement Xiaobi vs Herrman
# The coefficients created by tab5 are in Xiaobi's notation
# The equations just below follow Herrman's notation; hence
self.interpKernels[:,:,1] *= -1.
self.interpKernels[:,:,4] *= -1.
self.interpKernels[:,:,7] *= -1.
# Set it to True
self.interpolationDone = True
else:
if self.verbose:
print('Use interpolated kernels')
# Get what's done
kernels = self.interpKernels
# Vertical component (positive down)
ws = M[:,np.newaxis,1]*( kernels[:,:,2]*c2az/2. - kernels[:,:,0]/6. + kernels[:,:,8]/3.) \
+ M[:,np.newaxis,2]*(-kernels[:,:,2]*c2az/2. - kernels[:,:,0]/6. + kernels[:,:,8]/3.) \
+ M[:,np.newaxis,0]*( kernels[:,:,0] + kernels[:,:,8])/3. \
+ M[:,np.newaxis,5]* kernels[:,:,2]*s2az \
+ M[:,np.newaxis,3]* kernels[:,:,1]*caz \
+ M[:,np.newaxis,4]* kernels[:,:,1]*saz
# Radial component (positive away from the source)
qr = M[:,np.newaxis,1]*( kernels[:,:,5]*c2az/2. - kernels[:,:,3]/6. + kernels[:,:,9]/3.) \
+ M[:,np.newaxis,2]*(-kernels[:,:,5]*c2az/2. - kernels[:,:,3]/6. + kernels[:,:,9]/3.) \
+ M[:,np.newaxis,0]*( kernels[:,:,3] + kernels[:,:,9])/3. \
+ M[:,np.newaxis,5]* kernels[:,:,5]*s2az \
+ M[:,np.newaxis,3]* kernels[:,:,4]*caz \
+ M[:,np.newaxis,4]* kernels[:,:,4]*saz
# Tangential component (positive if clockwise from zenithal view)
vt = M[:,np.newaxis,1]*kernels[:,:,7]*s2az/2. \
- M[:,np.newaxis,2]*kernels[:,:,7]*s2az/2. \
- M[:,np.newaxis,5]*kernels[:,:,7]*c2az \
+ M[:,np.newaxis,3]*kernels[:,:,6]*saz \
- M[:,np.newaxis,4]*kernels[:,:,6]*caz
# Cartesian components
Ux = qr*saz + vt*caz
Uy = qr*caz - vt*saz
Uz = -ws
# All done
return Ux.T, Uy.T, Uz.T
def createInterpolator(self, method='linear'):
'''
Create the interpolation method. This is based on scipy.interpolate.LinearNDInterpolator.
Returns:
* None
'''
# Create Arrays
depths = np.unique(self.depths)
distas = np.unique(self.distas)
values = self.zrtdsx.reshape((self.ndepth, self.ndista,10))
# Create the interpolators (if points fall outside the interpolating box,
# the value will be extrapolated)
self.interpolators = [sciint.RegularGridInterpolator((depths,
distas),
values[:,:,i],
method=method,
bounds_error=False,
fill_value=None) \
for i in range(10)]
# All done
return
def src2mom(self, slip, area, strike, dip, rake):
'''
Convert slip and point source geometry to moment.
Args:
* slip : Slip value (m).
* area : Area of the point source (m^2)
* strike : Strike angle (rad)
* dip : Dip angle (rad)
* rake : Rake angle (rad, 0 left-lateral strike slip, 2pi pure thrust)
Returns:
* M : Moment tensor
'''
# Nor
nor = []
nor.append(-1.*np.sin(dip)*np.sin(strike))
nor.append(np.sin(dip)*np.cos(strike))
nor.append(-1.*np.cos(dip))
# Sli
s = []
s.append((np.cos(rake)*np.cos(strike)+np.cos(dip)*np.sin(rake)*np.sin(strike))*slip)
s.append((np.cos(rake)*np.sin(strike)-np.cos(dip)*np.sin(rake)*np.cos(strike))*slip)
s.append((-1.*np.sin(rake)*np.sin(dip))*slip)
# Iterate
Maki = np.zeros((3,3,len(dip)))
for ix in range(3):
for iy in range(3):
Maki[iy,ix,:] += nor[ix]*s[iy] + nor[iy]*s[ix]
# Order
M = np.zeros((len(dip),6))
M[:,0] = Maki[2,2,:]
M[:,1] = Maki[0,0,:]
M[:,2] = Maki[1,1,:]
M[:,3] = Maki[2,0,:]
M[:,4] = Maki[2,1,:]
M[:,5] = Maki[1,0,:]
# All done
return area[:,np.newaxis]*M
# PRIVATE METHODS
def _getGeometry(self, xs, ys, zs, xr, yr):
'''
Returns some geometrical features
Args:
* xs, ys, zs : Source location (floats or np.array)
* xr, yr : Receiver location (floats or np.array)
Returns:
* distance, depth, caz, saz, c2az, s2az
'''
# Machine precision
eps = np.finfo(float).eps
# Compute geometry
distance = np.sqrt( (xs[:,np.newaxis] - xr[np.newaxis,:])**2 +\
(ys[:,np.newaxis] - yr[np.newaxis,:])**2)
depth = zs[:,np.newaxis]*np.ones(distance.shape)
caz = (yr[np.newaxis,:] - ys[:,np.newaxis])/distance
saz = (xr[np.newaxis,:] - xs[:,np.newaxis])/distance
caz[distance<=eps] = 1.
saz[distance<=eps] = 0.
c2az = 2.*caz*caz - 1.
s2az = 2.*saz*caz
return distance, depth, caz, saz, c2az, s2az
#EOF