-
Notifications
You must be signed in to change notification settings - Fork 32
/
TriangularPatches.py
3230 lines (2629 loc) · 110 KB
/
TriangularPatches.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
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
'''
A parent class that deals with triangular patches fault
Written by Bryan Riel, Z. Duputel and R. Jolivet November 2013
'''
# Externals
import numpy as np
import pyproj as pp
import matplotlib.pyplot as plt
import scipy.interpolate as sciint
from . import triangularDisp as tdisp
from scipy.linalg import block_diag
import copy
import sys
import os
# Personals
from .Fault import Fault
from .geodeticplot import geodeticplot as geoplot
from .gps import gps as gpsclass
class TriangularPatches(Fault):
'''
Classes implementing a fault made of triangular patches. Inherits from Fault
Args:
* name : Name of the fault.
Kwargs:
* utmzone : UTM zone (optional, default=None)
* lon0 : Longitude of the center of the UTM zone
* lat0 : Latitude of the center of the UTM zone
* ellps : ellipsoid (optional, default='WGS84')
* verbose : Speak to me (default=True)
'''
# ----------------------------------------------------------------------
def __init__(self, name, utmzone=None, ellps='WGS84', verbose=True, lon0=None, lat0=None):
# Base class init
super(TriangularPatches,self).__init__(name, utmzone=utmzone, ellps=ellps, lon0=lon0, lat0=lat0, verbose=verbose)
# Specify the type of patch
self.patchType = 'triangle'
# The case of vertical faults with triangular patches is tricky, so we leave this option here
self.vertical = False
# All done
return
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def setdepth(self):
'''
Set depth patch attributes
Returns:
* None
'''
# Set depth
self.depth = np.max([v[2] for v in self.Vertices])
self.z_patches = np.linspace(self.depth, 0.0, 5)
# All done
return
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def computeArea(self):
'''
Computes the area of all triangles.
Returns:
* None
'''
# Area
self.area = []
# Loop over patches
for patch in self.patch:
self.area.append(self.patchArea(patch))
# all done
return
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def patchArea(self, patch):
'''
Returns the area of one patch.
Args:
* patch : one item of the patch list.
Returns:
* Area : float
'''
# Get vertices of patch
p1, p2, p3 = list(patch)
# Compute side lengths
a = np.sqrt((p2[0] - p1[0])**2 + (p2[1] - p1[1])**2 + (p2[2] - p1[2])**2)
b = np.sqrt((p3[0] - p2[0])**2 + (p3[1] - p2[1])**2 + (p3[2] - p2[2])**2)
c = np.sqrt((p1[0] - p3[0])**2 + (p1[1] - p3[1])**2 + (p1[2] - p3[2])**2)
# Compute area using numerically stable Heron's formula
c,b,a = np.sort([a, b, c])
area = 0.25 * np.sqrt((a + (b + c)) * (c - (a - b))
* (c + (a - b)) * (a + (b - c)))
# All Done
return area
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def splitPatch(self, patch):
'''
Splits a patch into 4 patches, based on the mid-point of each side.
Args:
* patch : item of the patch list.
Returns:
* t1, t2, t3, t4 : 4 patches
'''
# Get corners
p1, p2, p3 = list(patch)
if type(p1) is not list:
p1 = p1.tolist()
if type(p2) is not list:
p2 = p2.tolist()
if type(p3) is not list:
p3 = p3.tolist()
# Compute mid-points
p12 = [p1[0] + (p2[0]-p1[0])/2.,
p1[1] + (p2[1]-p1[1])/2.,
p1[2] + (p2[2]-p1[2])/2.]
p23 = [p2[0] + (p3[0]-p2[0])/2.,
p2[1] + (p3[1]-p2[1])/2.,
p2[2] + (p3[2]-p2[2])/2.]
p31 = [p3[0] + (p1[0]-p3[0])/2.,
p3[1] + (p1[1]-p3[1])/2.,
p3[2] + (p1[2]-p3[2])/2.]
# make 4 triangles
t1 = [p1, p12, p31]
t2 = [p12, p2, p23]
t3 = [p31, p23, p3]
t4 = [p31, p12, p23]
# All done
return t1, t2, t3, t4
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def selectPatches(self,minlon,maxlon,minlat,maxlat,mindep,maxdep):
'''
Removes patches that are outside of a 3D box.
Args:
* minlon : west longitude
* maxlon : east longitude
* minlat : south latitude
* maxlat : north latitude
* mindep : Minimum depth
* maxdep : Maximum depth
Returns:
* None
'''
xmin,ymin = self.ll2xy(minlon,minlat)
xmax,ymax = self.ll2xy(maxlon,maxlat)
for p in range(len(self.patch)-1,-1,-1):
x1, x2, x3, width, length, strike, dip = self.getpatchgeometry(p)
if x1<xmin or x1>xmax or x2<ymin or x2>ymax or x3<mindep or x3>maxdep:
self.deletepatch(p)
for i in range(len(self.xf)-1,-1,-1):
x1 = self.xf[i]
x2 = self.yf[i]
if x1<xmin or x1>xmax or x2<ymin or x2>ymax:
self.xf = np.delete(self.xf,i)
self.yf = np.delete(self.yf,i)
return
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def vertices2ll(self):
'''
Converts all the vertices into lonlat coordinates.
Returns:
* None
'''
# Create a list
vertices_ll = []
# iterate
for vertex in self.Vertices:
lon, lat = self.xy2ll(vertex[0], vertex[1])
vertices_ll.append([lon, lat, vertex[2]])
# Save
self.Vertices_ll = np.array(vertices_ll)
# All done
return
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def setVerticesFromPatches(self):
'''
Takes the patches and constructs a list of Vertices and Faces
Returns:
* None
'''
# Get patches
patches = self.patch
# Create lists
vertices = []
faces = []
# Iterate to build vertices
for patch in patches:
for vertex in patch.tolist():
if vertex not in vertices:
vertices.append(vertex)
# Iterate to build Faces
for patch in patches:
face = []
for vertex in patch.tolist():
uu = np.flatnonzero(np.array([vertex==v for v in vertices]))[0]
face.append(uu)
faces.append(face)
# Set them
self.Vertices = np.array(vertices)
self.Faces = np.array(faces)
# 2 lon lat
self.vertices2ll()
# All done
return
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def patches2triangles(self, fault, numberOfTriangles=4):
'''
Takes a fault with rectangular patches and splits them into triangles to
initialize self.
Args:
* fault : instance of rectangular patches.
Kwargs:
* numberOfTriangles : Split each patch in 2 or 4 (default) triangle
Returns:
* None
'''
# Initialize the lists of patches
self.patch = []
self.patchll = []
# Initialize vertices and faces
vertices = []
faces = []
# Each patch is being splitted in 2 or 4 triangles
for patch in fault.patch:
# Add vertices
for vertex in patch.tolist():
if vertex not in vertices:
vertices.append(vertex)
# Find the vertices in the list
i0 = np.flatnonzero(np.array([patch[0].tolist()==v for v in vertices]))[0]
i1 = np.flatnonzero(np.array([patch[1].tolist()==v for v in vertices]))[0]
i2 = np.flatnonzero(np.array([patch[2].tolist()==v for v in vertices]))[0]
i3 = np.flatnonzero(np.array([patch[3].tolist()==v for v in vertices]))[0]
if numberOfTriangles==4:
# Get the center
center = np.array(fault.getpatchgeometry(patch, center=True)[:3])
vertices.append(list(center))
ic = np.flatnonzero(np.array([center.tolist()==v for v in vertices]))[0]
# Split in 4
t1 = np.array([patch[0], patch[1], center])
t2 = np.array([patch[1], patch[2], center])
t3 = np.array([patch[2], patch[3], center])
t4 = np.array([patch[3], patch[0], center])
# faces
fs = ([i0, i1, ic],
[i1, i2, ic],
[i2, i3, ic],
[i3, i0, ic])
# patches
ps = [t1, t2, t3, t4]
elif numberOfTriangles==2:
# Split in 2
t1 = np.array([patch[0], patch[1], patch[2]])
t2 = np.array([patch[2], patch[3], patch[0]])
# faces
fs = ([i0, i1, i2], [i2, i3, i0])
# patches
ps = (t1, t2)
else:
assert False, 'numberOfTriangles should be 2 or 4'
for f,p in zip(fs, ps):
faces.append(f)
self.patch.append(p)
# Save
self.Vertices = np.array(vertices)
self.Faces = np.array(faces)
# Convert
self.vertices2ll()
self.patch2ll()
# Initialize slip
self.initializeslip()
# Set fault trace
self.xf = fault.xf
self.yf = fault.yf
self.lon = fault.lon
self.lat = fault.lat
if hasattr(fault, 'xi'):
self.xi = fault.xi
self.yi = fault.yi
self.loni = fault.loni
self.lati = fault.lati
# Set depth
self.setdepth()
# All done
return
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def setMu(self,model_file):
'''
Set shear modulus values for seismic moment calculation
from an edks model_file
Args:
* model_file: EDKS .model file
Returns:
* None
'''
# Read model file
mu = []
depth = 0.
depths = []
with open(model_file) as f:
l1 = f.readline()
items = l1.strip().split()
assert len(items)==2, 'Incorrect first line format in %s'%(model_file)
nd = int(items[0])
fc = float(items[1])
for l in f:
items = l.strip().split()
assert len(items)==4, 'Incorrect line format in %s'%(model_file)
RHO = float(items[0])*fc
VP = float(items[1])*fc
VS = float(items[2])*fc
H = float(items[3])
mu.append(VS*VS*RHO)
if H==0.:
H = np.inf
depths.append([depth,depth+H])
depth += H
Nd = len(depths)
assert Nd==nd, 'Incorrect number of layes in %s (%d vs %d)'%(model_file,Nd,nd)
Np = len(self.patch)
# Set Mu for each patch
self.mu = np.zeros((Np,))
for p in range(Np):
p_x, p_y, p_z, width, length, strike_rad, dip_rad = self.getpatchgeometry(p,center=True)
for d in range(Nd):
if p_z>=depths[d][0] and p_z<depths[d][1]:
self.mu[p] = mu[d]
# All done
return
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def readPatchesFromFile(self, filename, readpatchindex=True,
donotreadslip=False, gmtslip=True,
inputCoordinates='lonlat'):
'''
Reads patches from a GMT formatted file.
Args:
* filename : Name of the file
Kwargs:
* inputCoordinates : Default is 'lonlat'. Can be 'utm'
* readpatchindex : Default True.
* donotreadslip : Default is False. If True, does not read the slip
* gmtslip : A -Zxxx is in the header of each patch
* inputCoordinates : Default is 'lonlat', can be 'xyz'
Returns:
* None
'''
# create the lists
self.patch = []
self.patchll = []
if readpatchindex:
self.index_parameter = []
if not donotreadslip:
Slip = []
# open the files
fin = open(filename, 'r')
# read all the lines
A = fin.readlines()
# depth
D = 0.0
d = 10000.
# Index
if gmtslip:
ipatch = 3
else:
ipatch = 2
# Loop over the file
i = 0
while i<len(A):
# Assert it works
assert A[i].split()[0]=='>', 'Not a patch, reformat your file...'
# Get the Patch Id
if readpatchindex:
self.index_parameter.append([int(A[i].split()[ipatch]),int(A[i].split()[ipatch+1]),int(A[i].split()[ipatch+2])])
# Get the slip value
if not donotreadslip:
if len(A[i].split())>7:
slip = np.array([float(A[i].split()[ipatch+4]), float(A[i].split()[ipatch+5]), float(A[i].split()[ipatch+6])])
else:
slip = np.array([0.0, 0.0, 0.0])
Slip.append(slip)
# get the values
if inputCoordinates in ('lonlat'):
lon1, lat1, z1 = A[i+1].split()
lon2, lat2, z2 = A[i+2].split()
lon3, lat3, z3 = A[i+3].split()
# Pass as floating point
lon1 = float(lon1); lat1 = float(lat1); z1 = float(z1)
lon2 = float(lon2); lat2 = float(lat2); z2 = float(z2)
lon3 = float(lon3); lat3 = float(lat3); z3 = float(z3)
# translate to utm
x1, y1 = self.ll2xy(lon1, lat1)
x2, y2 = self.ll2xy(lon2, lat2)
x3, y3 = self.ll2xy(lon3, lat3)
elif inputCoordinates in ('xyz'):
x1, y1, z1 = A[i+1].split()
x2, y2, z2 = A[i+2].split()
x3, y3, z3 = A[i+3].split()
# Pass as floating point
x1 = float(x1); y1 = float(y1); z1 = float(z1)
x2 = float(x2); y2 = float(y2); z2 = float(z2)
x3 = float(x3); y3 = float(y3); z3 = float(z3)
# translate to utm
lon1, lat1 = self.xy2ll(x1, y1)
lon2, lat2 = self.xy2ll(x2, y2)
lon3, lat3 = self.xy2ll(x3, y3)
# Depth
mm = min([float(z1), float(z2), float(z3)])
if D<mm:
D=mm
if d>mm:
d=mm
# Set points
p1 = [x1, y1, z1]; p1ll = [lon1, lat1, z1]
p2 = [x2, y2, z2]; p2ll = [lon2, lat2, z2]
p3 = [x3, y3, z3]; p3ll = [lon3, lat3, z3]
# Store these
p = [p1, p2, p3]
pll = [p1ll, p2ll, p3ll]
p = np.array(p)
pll = np.array(pll)
# Store these in the lists
self.patch.append(p)
self.patchll.append(pll)
# increase i
i += 4
# Close the file
fin.close()
# depth
self.depth = D
self.top = d
self.z_patches = np.linspace(0,D,5)
self.factor_depth = 1.
# Patches 2 vertices
self.setVerticesFromPatches()
self.numpatch = self.Faces.shape[0]
# Translate slip to np.array
if not donotreadslip:
self.initializeslip(values=np.array(Slip))
else:
self.initializeslip()
if readpatchindex:
self.index_parameter = np.array(self.index_parameter)
# All done
return
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def readPatchesFromTriFile(self, filename, readpatchindex=False,
donotreadslip=True,
inputCoordinates='lonlat'):
'''
Reads patches from a file with triangles having all three vertices on each line
modified from the GMT file reader
Args:
* filename : Name of the file
Kwargs:
* inputCoordinates : Default is 'lonlat'. Can be 'utm'
* readpatchindex : Default False.
* donotreadslip : Default is True. If True, does not read the slip
Returns:
* None
'''
negDepth = False
# create the lists
self.patch = []
self.patchll = []
if readpatchindex:
self.index_parameter = []
if not donotreadslip:
Slip = []
# open the files
fin = open(filename, 'r')
# read all the lines
A = fin.readlines()
# depth
D = 0.0
d = 10000.
# Loop over the file
i = 0
while i<len(A):
# Get the Patch Id
if readpatchindex:
self.index_parameter.append([int(A[i].split()[ipatch]),int(A[i].split()[ipatch+1]),int(A[i].split()[ipatch+2])])
# Get the slip value
if not donotreadslip:
if len(A[i].split())>7:
slip = np.array([float(A[i].split()[ipatch+4]), float(A[i].split()[ipatch+5]), float(A[i].split()[ipatch+6])])
else:
slip = np.array([0.0, 0.0, 0.0])
Slip.append(slip)
# get the values
if inputCoordinates in ('lonlat'):
# all three vertices on one line
lon1, lat1, z1, lon2, lat2, z2, lon3, lat3, z3 = A[i].split()
# Pass as floating point
lon1 = float(lon1); lat1 = float(lat1); z1 = float(z1)
lon2 = float(lon2); lat2 = float(lat2); z2 = float(z2)
lon3 = float(lon3); lat3 = float(lat3); z3 = float(z3)
if negDepth and z1 > 0.0: # check depth sign, assuming all three are the same sign
z1 = -1*z1
z2 = -1*z2
z3 = -1*z3
# translate to utm
x1, y1 = self.ll2xy(lon1, lat1)
x2, y2 = self.ll2xy(lon2, lat2)
x3, y3 = self.ll2xy(lon3, lat3)
elif inputCoordinates in ('xyz'):
# all three vertices on one line
lon1, lat1, z1, lon2, lat2, z2, lon3, lat3, z3 = A[i].split()
# Pass as floating point
x1 = float(x1); y1 = float(y1); z1 = float(z1)
x2 = float(x2); y2 = float(y2); z2 = float(z2)
x3 = float(x3); y3 = float(y3); z3 = float(z3)
# translate to utm
lon1, lat1 = self.xy2ll(x1, y1)
lon2, lat2 = self.xy2ll(x2, y2)
lon3, lat3 = self.xy2ll(x3, y3)
# Depth
mm = min([float(z1), float(z2), float(z3)])
if D<mm:
D=mm
if d>mm:
d=mm
# Set points
p1 = [x1, y1, z1]; p1ll = [lon1, lat1, z1]
p2 = [x2, y2, z2]; p2ll = [lon2, lat2, z2]
p3 = [x3, y3, z3]; p3ll = [lon3, lat3, z3]
# Store these
p = [p1, p2, p3]
pll = [p1ll, p2ll, p3ll]
p = np.array(p)
pll = np.array(pll)
# Store these in the lists
self.patch.append(p)
self.patchll.append(pll)
# increase i
i += 1
# Close the file
fin.close()
# depth
self.depth = D
self.top = d
self.z_patches = np.linspace(0,D,5)
self.factor_depth = 1.
# Patches 2 vertices
self.setVerticesFromPatches()
self.numpatch = self.Faces.shape[0]
# Translate slip to np.array
if not donotreadslip:
self.initializeslip(values=np.array(Slip))
else:
self.initializeslip()
if readpatchindex:
self.index_parameter = np.array(self.index_parameter)
# All done
return
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def readGocadPatches(self, filename, neg_depth=False, utm=False, factor_xy=1.0,
factor_depth=1.0, verbose=False):
'''
Load a triangulated surface from a Gocad formatted file. Vertices
must be in geographical coordinates.
Args:
* filename: tsurf file to read
Kwargs:
* neg_depth: if true, use negative depth
* utm: if true, input file is given as utm coordinates (if false -> lon/lat)
* factor_xy: if utm==True, multiplication factor for x and y
* factor_depth: multiplication factor for z
* verbose: Speak to me
Returns:
* None
'''
# Initialize the lists of patches
self.patch = []
self.patchll = []
# Factor to correct input negative depths (we want depths to be positive)
if neg_depth:
negFactor = -1.0
else:
negFactor = 1.0
# Get the geographic vertices and connectivities from the Gocad file
with open(filename, 'r') as fid:
vertices = []
vids = []
faces = []
for line in fid:
if line.startswith('VRTX'):
items = line.split()
name, vid, x, y, z = items[:5]
vids.append(vid)
vertices.append([float(x), float(y), negFactor*float(z)])
elif line.startswith('TRGL'):
name, p1, p2, p3 = line.split()
faces.append([int(p1), int(p2), int(p3)])
fid.close()
vids = np.array(vids,dtype=int)
i0 = np.min(vids)
vids = vids - i0
i = np.argsort(vids)
vertices = np.array(vertices, dtype=float)[i,:]
faces = np.array(faces, dtype=int) - i0
# Resample vertices to UTM
if utm:
vx = vertices[:,0].copy()*factor_xy
vy = vertices[:,1].copy()*factor_xy
vertices[:,0],vertices[:,1] = self.xy2ll(vx,vy)
else:
vx, vy = self.ll2xy(vertices[:,0], vertices[:,1])
vz = vertices[:,2]*factor_depth
self.factor_depth = factor_depth
self.Vertices = np.column_stack((vx, vy, vz))
self.Vertices_ll = vertices
self.Faces = faces
if verbose:
print('min/max depth: {} km/ {} km'.format(vz.min(),vz.max()))
print('min/max lat: {} deg/ {} deg'.format(vertices[:,1].min(),vertices[:,1].max()))
print('min/max lon: {} deg/ {} deg'.format(vertices[:,0].min(),vertices[:,0].max()))
print('min/max x: {} km/ {} km'.format(vx.min(),vx.max()))
print('min/max y: {} km/ {} km'.format(vy.min(),vy.max()))
# Loop over faces and create a triangular patch consisting of coordinate tuples
self.numpatch = faces.shape[0]
for i in range(self.numpatch):
# Get the indices of the vertices
v1, v2, v3 = faces[i,:]
# Get the coordinates
x1, y1, lon1, lat1, z1 = vx[v1], vy[v1], vertices[v1,0], vertices[v1,1], vz[v1]
x2, y2, lon2, lat2, z2 = vx[v2], vy[v2], vertices[v2,0], vertices[v2,1], vz[v2]
x3, y3, lon3, lat3, z3 = vx[v3], vy[v3], vertices[v3,0], vertices[v3,1], vz[v3]
# Make the coordinate tuples
p1 = [x1, y1, z1]; pll1 = [lon1, lat1, z1]
p2 = [x2, y2, z2]; pll2 = [lon2, lat2, z2]
p3 = [x3, y3, z3]; pll3 = [lon3, lat3, z3]
# Store the patch
self.patch.append(np.array([p1, p2, p3]))
self.patchll.append(np.array([pll1, pll2, pll3]))
# Update the depth of the bottom of the fault
if neg_depth:
self.top = np.max(vz)
self.depth = np.min(vz)
else:
self.top = np.min(vz)
self.depth = np.max(vz)
self.z_patches = np.linspace(self.depth, 0.0, 5)
self.initializeslip()
# All done
return
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def writeGocadPatches(self, filename, utm=False):
'''
Write a triangulated Gocad surface file.
Args:
* filename : output file name
Kwargs:
* utm : Write in utm coordinates if True
Returns:
* None
'''
# Get the geographic vertices and connectivities from the Gocad file
fid = open(filename, 'w')
if utm:
vertices = self.Vertices*1.0e3
else:
vertices = self.Vertices_ll
for i in range(vertices.shape[0]):
v = vertices[i]
fid.write('VRTX {} {} {} {}\n'.format(i+1,v[0],v[1],v[2]))
for i in range(self.Faces.shape[0]):
vid = self.Faces[i,:]+1
fid.write('TRGL {} {} {}\n'.format(vid[0],vid[1],vid[2]))
fid.close()
# All done
return
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def getStrikes(self):
'''
Returns an array of strikes.
'''
# all done in one line
return np.array([self.getpatchgeometry(p)[5] for p in self.patch])
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def getDips(self):
'''
Returns an array of dips.
'''
# all done in one line
return np.array([self.getpatchgeometry(p)[6] for p in self.patch])
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def getDepths(self, center=True):
'''
Returns an array of depths.
Kwargs:
* center : If True, returns the center of the patches
'''
# All done in one line
return np.array([self.getpatchgeometry(p, center=True)[2] for p in self.patch])
# ----------------------------------------------------------------------
def writePatches2File(self, filename, add_slip=None, scale=1.0, stdh5=None, decim=1):
'''
Writes the patch corners in a file that can be used in psxyz.
Args:
* filename : Name of the file.
Kwargs:
* add_slip : Put the slip as a value for the color.
Can be None, strikeslip, dipslip, total, coupling
* scale : Multiply the slip value by a factor.
* patch : Can be 'normal' or 'equiv'
* stdh5 : Get the standard deviation from a h5 file
* decim : Decimate the h5 file
Returns:
* None
'''
# Check size
if self.N_slip!=None and self.N_slip!=len(self.patch) and add_slip is not None:
raise NotImplementedError('Only works for len(slip)==len(patch)')
# Write something
print('Writing geometry to file {}'.format(filename))
# Open the file
fout = open(filename, 'w')
# If an h5 file is specified, open it
if stdh5 is not None:
import h5py
h5fid = h5py.File(stdh5, 'r')
samples = h5fid['samples'].value[::decim,:]
# Loop over the patches
nPatches = len(self.patch)
for pIndex in range(nPatches):
# Select the string for the color
string = ' '
if add_slip is not None:
if add_slip=='coupling':
slp = self.coupling[pIndex]
string = '-Z{}'.format(slp)
if add_slip=='strikeslip':
if stdh5 is not None:
slp = np.std(samples[:,pIndex])
else:
slp = self.slip[pIndex,0]*scale
string = '-Z{}'.format(slp)
elif add_slip=='dipslip':
if stdh5 is not None:
slp = np.std(samples[:,pIndex+nPatches])
else:
slp = self.slip[pIndex,1]*scale
string = '-Z{}'.format(slp)
elif add_slip=='tensile':
if stdh5 is not None:
slp = np.std(samples[:,pIndex+2*nPatches])
else:
slp = self.slip[pIndex,2]*scale
string = '-Z{}'.format(slp)
elif add_slip=='total':
if stdh5 is not None:
slp = np.std(samples[:,pIndex]**2 + samples[:,pIndex+nPatches]**2)
else:
slp = np.sqrt(self.slip[pIndex,0]**2 + self.slip[pIndex,1]**2)*scale
string = '-Z{}'.format(slp)
# Put the parameter number in the file as well if it exists
parameter = ' '
if hasattr(self,'index_parameter') and add_slip is not None:
i = int(self.index_parameter[pIndex,0])
j = int(self.index_parameter[pIndex,1])
k = int(self.index_parameter[pIndex,2])
parameter = '# {} {} {} '.format(i,j,k)
else:
parameter = '# 9999 9999 9999 '
# Put the slip value
if add_slip is not None:
if add_slip=='coupling':
slipstring = ' # {}'.format(self.coupling[pIndex])
else:
slipstring = ' # {} {} {} '.format(self.slip[pIndex,0],
self.slip[pIndex,1], self.slip[pIndex,2])
# Write the string to file
if add_slip is None:
fout.write('> {} {} \n'.format(string, parameter))
else:
fout.write('> {} {} {} \n'.format(string,parameter,slipstring))
# Write the 3 patch corners (the order is to be GMT friendly)
p = self.patchll[pIndex]
pp = p[0]; fout.write('{} {} {} \n'.format(np.round(pp[0], decimals=4),
np.round(pp[1], decimals=4),
np.round(pp[2], decimals=4)))
pp = p[1]; fout.write('{} {} {} \n'.format(np.round(pp[0], decimals=4),
np.round(pp[1], decimals=4),
np.round(pp[2], decimals=4)))
pp = p[2]; fout.write('{} {} {} \n'.format(np.round(pp[0], decimals=4),
np.round(pp[1], decimals=4),
np.round(pp[2], decimals=4)))
# Close the file
fout.close()
# Close h5 file if it is open
if stdh5 is not None:
h5fid.close()
# All done
return
# ----------------------------------------------------------------------
# ----------------------------------------------------------------------
def writeSlipDirection2File(self, filename, scale=1.0, factor=1.0,
neg_depth=False, ellipse=False,nsigma=1.):
'''
Write a psxyz compatible file to draw lines starting from the center
of each patch, indicating the direction of slip. Scale can be a real
number or a string in 'total', 'strikeslip', 'dipslip' or 'tensile'
Args:
* filename : Name of the output file
Kwargs:
* scale : Scale of the line
* factor : Multiply slip by a factor
* neg_depth : if True, depth is a negative nmber
* ellipse : Write the ellipse
* nsigma : Nxsigma for the ellipse
Returns:
* None
'''
# Copmute the slip direction
self.computeSlipDirection(scale=scale, factor=factor, ellipse=ellipse,nsigma=nsigma, neg_depth=neg_depth)
# Write something
print('Writing slip direction to file {}'.format(filename))
# Open the file
fout = open(filename, 'w')
# Loop over the patches
for p in self.slipdirection:
# Write the > sign to the file
fout.write('> \n')