-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathseismic_phase.py
1547 lines (1409 loc) · 72.2 KB
/
seismic_phase.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
#!/usr/bin/env python
# -*- coding: utf-8 -*-
"""
Objects and functions dealing with seismic phases.
"""
from __future__ import (absolute_import, division, print_function,
unicode_literals)
from future.builtins import * # NOQA
from future.utils import raise_from
from itertools import count
import math
import re
import numpy as np
from obspy.core.util.obspy_types import Enum
from .helper_classes import (Arrival, SlownessModelError, TauModelError,
TimeDist)
from .c_wrappers import clibtau
REFINE_DIST_RADIAN_TOL = 0.0049 * math.pi / 180
_ACTIONS = Enum([
# Used by add_to_branch when the path turns within a segment. We assume
# that no ray will turn downward so turning implies turning from downward
# to upward, ie U.
"turn",
# Used by add_to_branch when the path reflects off the top of the end of
# a segment, ie ^.
"reflect_underside",
# Used by add_to_branch when the path reflects off the bottom of the end
# of a segment, ie v.
"reflect_topside",
# Used by add_to_branch when the path transmits up through the end of a
# segment.
"transup",
# Used by add_to_branch when the path transmits down through the end of a
# segment.
"transdown"
])
class SeismicPhase(object):
"""
Stores and transforms seismic phase names to and from their
corresponding sequence of branches. Will maybe contain "expert" mode
wherein paths may start in the core. Principal use is to calculate leg
contributions for scattered phases. Nomenclature: "K" - downgoing wave
from source in core; "k" - upgoing wave from source in core.
"""
def __init__(self, name, tau_model, receiver_depth=0.0):
# The phase name, e.g. PKiKP.
self.name = name
# The receiver depth within the TauModel that was used to generate this
# phase. Normally this is 0.0 for a surface station, but can be
# different for borehole or scattering calculations.
self.receiver_depth = receiver_depth
# TauModel to generate phase for.
self.tau_model = tau_model
# The source depth within the TauModel that was used to generate
# this phase.
self.source_depth = self.tau_model.source_depth
# List containing strings for each leg.
self.legs = leg_puller(name)
# Name with depths corrected to be actual discontinuities in the model.
self.purist_name = self.create_purist_name(tau_model)
# Settings for this instance. Should eventually be configurable.
self._settings = {
# The maximum degrees that a Pn or Sn can refract along the moho.
# Note this is not the total distance, only the segment along the
# moho.
"max_refraction_in_radians": np.radians(20.0),
# The maximum degrees that a Pdiff or Sdiff can diffract along the
# CMB. Note this is not the total distance, only the segment along
# the CMB.
"max_diffraction_in_radians": np.radians(60.0),
# The maximum number of refinements to make to an Arrival.
"max_recursion": 5
}
# Enables phases originating in core.
self.expert = False
# Minimum/maximum ray parameters that exist for this phase.
self.min_ray_param = None
self.max_ray_param = None
# Index within TauModel.ray_param that corresponds to max_ray_param.
# Note that max_ray_param_index < min_ray_param_index as ray parameter
# decreases with increasing index.
self.max_ray_param_index = -1
# Index within TauModel.ray_param that corresponds to min_ray_param.
# Note that max_ray_param_index < min_ray_param_index as ray parameter
# decreases with increasing index.
self.min_ray_param_index = -1
# Temporary branch numbers determining where to start adding to the
# branch sequence.
self.current_branch = None
# Array of distances corresponding to the ray parameters stored in
# ray_param.
self.dist = None
# Array of times corresponding to the ray parameters stored in
# ray_param.
self.time = None
# Array of possible ray parameters for this phase.
self.ray_param = None
# The minimum distance that this phase can be theoretically observed.
self.min_distance = 0.0
# The maximum distance that this phase can be theoretically observed.
self.max_distance = 1e300
# List (could make array!) of branch numbers for the given phase.
# Note that this depends upon both the planet model and the source
# depth.
self.branch_seq = []
# True if the current leg of the phase is down going. This allows a
# check to make sure the path is correct.
# Used in addToBranch() and parseName().
self.down_going = []
# ArrayList of wave types corresponding to each leg of the phase.
self.wave_type = []
self.parse_name(tau_model)
self.sum_branches(tau_model)
def create_purist_name(self, tau_model):
current_leg = self.legs[0]
# Deal with surface wave veocities first, since they are a special
# case.
if len(self.legs) == 2 and current_leg.endswith("kmps"):
purist_name = self.name
return purist_name
purist_name = ""
# Only loop to penultimate element as last leg is always "END".
for current_leg in self.legs[:-1]:
# Find out if the next leg represents a phase conversion or
# reflection depth.
if current_leg[0] in "v^":
discon_branch = closest_branch_to_depth(tau_model,
current_leg[1:])
leg_depth = tau_model.tau_branches[0, discon_branch].top_depth
purist_name += current_leg[0]
purist_name += str(int(round(leg_depth)))
else:
try:
float(current_leg)
except ValueError:
# If current_leg is just a string:
purist_name += current_leg
else:
# If it is indeed a number:
discon_branch = closest_branch_to_depth(tau_model,
current_leg)
leg_depth = \
tau_model.tau_branches[0, discon_branch].top_depth
purist_name += str(leg_depth)
return purist_name
def parse_name(self, tau_model):
"""
Construct a branch sequence from the given phase name and tau model.
"""
current_leg = self.legs[0]
next_leg = current_leg
# Deal with surface wave velocities first, since they are a special
# case.
if len(self.legs) == 2 and current_leg.endswith("kmps"):
return
# Make a check for J legs if the model doesn't allow J:
if "J" in self.name and not tau_model.s_mod.allow_inner_core_s:
raise TauModelError("J phases are not created for this model: {}"
.format(self.name))
# Set currWave to be the wave type for this leg, P or S
if current_leg in ("p", "K", "k", "I") or current_leg[0] == "P":
is_p_wave = True
is_p_wave_previous = is_p_wave
elif current_leg in ("s", "J") or current_leg[0] == "S":
is_p_wave = False
is_p_wave_previous = is_p_wave
else:
raise TauModelError('Unknown starting phase: ' + current_leg)
# First, decide whether the ray is upgoing or downgoing from the
# source. If it is up going then the first branch number would be
# model.source_branch-1 and downgoing would be model.source_branch.
upgoing_rec_branch = tau_model.find_branch(self.receiver_depth)
downgoing_rec_branch = upgoing_rec_branch - 1 # One branch shallower.
if current_leg[0] in "sS":
# Exclude S sources in fluids.
sdep = tau_model.source_depth
if tau_model.cmb_depth < sdep < tau_model.iocb_depth:
self.max_ray_param, self.min_ray_param = -1, -1
return
# Set self.max_ray_param to be a horizontal ray leaving the source and
# self.min_ray_param to be a vertical (p=0) ray.
if current_leg[0] in "PS" or (self.expert and current_leg[0] in "KIJ"):
# Downgoing from source.
self.current_branch = tau_model.source_branch
# Treat initial downgoing as if it were an underside reflection.
end_action = _ACTIONS["reflect_underside"]
try:
s_layer_num = tau_model.s_mod.layer_number_below(
tau_model.source_depth, is_p_wave_previous)
layer = tau_model.s_mod.get_slowness_layer(
s_layer_num, is_p_wave_previous)
self.max_ray_param = layer['top_p']
except SlownessModelError as e:
raise_from(RuntimeError('Please contact the developers. This '
'error should not occur.'), e)
self.max_ray_param = tau_model.get_tau_branch(
tau_model.source_branch, is_p_wave).max_ray_param
elif current_leg in ("p", "s") or (
self.expert and current_leg[0] == "k"):
# Upgoing from source: treat initial downgoing as if it were a
# topside reflection.
end_action = _ACTIONS["reflect_topside"]
try:
s_layer_num = tau_model.s_mod.layer_number_above(
tau_model.source_depth, is_p_wave_previous)
layer = tau_model.s_mod.get_slowness_layer(
s_layer_num, is_p_wave_previous)
self.max_ray_param = layer['bot_p']
except SlownessModelError as e:
raise_from(RuntimeError('Please contact the developers. This '
'error should not occur.'), e)
if tau_model.source_branch != 0:
self.current_branch = tau_model.source_branch - 1
else:
# p and s for zero source depth are only at zero distance
# and then can be called P or S.
self.max_ray_param = -1
self.min_ray_param = -1
return
else:
raise TauModelError(
'First phase not recognised {}: Must be one of P, Pg, Pn, '
'Pdiff, p, Ped or the S equivalents.'.format(current_leg))
if self.receiver_depth != 0:
if self.legs[-2] in ('Ped', 'Sed'):
# Downgoing at receiver
self.max_ray_param = min(
tau_model.get_tau_branch(
downgoing_rec_branch, is_p_wave).min_turn_ray_param,
self.max_ray_param)
else:
# upgoing at receiver
self.max_ray_param = min(
tau_model.get_tau_branch(upgoing_rec_branch,
is_p_wave).min_turn_ray_param,
self.max_ray_param)
self.min_ray_param = 0
is_leg_depth, is_next_leg_depth = False, False
# Now loop over all the phase legs and construct the proper branch
# sequence.
current_leg = "START" # So the prev_leg isn't wrong on the first pass.
for leg_num in range(len(self.legs) - 1):
prev_leg = current_leg
current_leg = next_leg
next_leg = self.legs[leg_num + 1]
is_leg_depth = is_next_leg_depth
# Find out if the next leg represents a phase conversion depth.
try:
next_leg_depth = float(next_leg)
is_next_leg_depth = True
except ValueError:
next_leg_depth = -1
is_next_leg_depth = False
# Set currWave to be the wave type for this leg, "P" or "S".
is_p_wave_previous = is_p_wave
if current_leg in ("p", "k", "I") or current_leg[0] == "P":
is_p_wave = True
elif current_leg in ("s", "J") or current_leg[0] == "S":
is_p_wave = False
elif current_leg == "K":
# Here we want to use whatever is_p_wave was on the last leg
# so do nothing. This makes sure we use the correct
# max_ray_param from the correct TauBranch within the outer
# core. In other words K has a high slowness zone if it
# entered the outer core as a mantle P wave, but doesn't if
# it entered as a mantle S wave. It shouldn't matter for
# inner core to outer core type legs.
pass
# Check to see if there has been a phase conversion.
if len(self.branch_seq) > 0 and is_p_wave_previous != is_p_wave:
self.phase_conversion(tau_model, self.branch_seq[-1],
end_action, is_p_wave_previous)
if current_leg in ('Ped', 'Sed'):
if next_leg == "END":
if self.receiver_depth > 0:
end_action = _ACTIONS["reflect_underside"]
self.add_to_branch(tau_model, self.current_branch,
downgoing_rec_branch, is_p_wave,
end_action)
else:
# This should be impossible except for 0 dist 0 source
# depth which can be called p or P.
self.max_ray_param = -1
self.min_ray_param = -1
return
else:
raise TauModelError(
"Phase not recognized: {} followed by {}".format(
current_leg, next_leg))
# Deal with p and s case.
elif current_leg in ("p", "s", "k"):
if next_leg[0] == "v":
raise TauModelError(
"p and s must always be upgoing and cannot come "
"immediately before a top-sided reflection.")
elif next_leg.startswith("^"):
discon_branch = closest_branch_to_depth(
tau_model, next_leg[1:])
if self.current_branch >= discon_branch:
end_action = _ACTIONS["reflect_underside"]
self.add_to_branch(
tau_model, self.current_branch, discon_branch,
is_p_wave, end_action)
else:
raise TauModelError(
"Phase not recognised: {} followed by {} when "
"current_branch > discon_branch".format(
current_leg, next_leg))
elif next_leg == "m" and \
self.current_branch >= tau_model.moho_branch:
end_action = _ACTIONS["transup"]
self.add_to_branch(
tau_model, self.current_branch, tau_model.moho_branch,
is_p_wave, end_action)
elif next_leg[0] in ("P", "S") or next_leg in ("K", "END"):
if next_leg == 'END':
discon_branch = upgoing_rec_branch
elif next_leg == 'K':
discon_branch = tau_model.cmb_branch
else:
discon_branch = 0
if current_leg == 'k' and next_leg != 'K':
end_action = _ACTIONS["transup"]
else:
end_action = _ACTIONS["reflect_underside"]
self.add_to_branch(
tau_model, self.current_branch, discon_branch,
is_p_wave, end_action)
elif is_next_leg_depth:
discon_branch = closest_branch_to_depth(tau_model,
next_leg)
end_action = _ACTIONS["transup"]
self.add_to_branch(
tau_model, self.current_branch, discon_branch,
is_p_wave, end_action)
else:
raise TauModelError(
"Phase not recognized: {} followed by {}".format(
current_leg, next_leg))
# Now deal with P and S case.
elif current_leg in ("P", "S"):
if next_leg in ("P", "S", "Pn", "Sn", "END"):
if end_action == _ACTIONS["transdown"] or \
end_action == _ACTIONS["reflect_underside"]:
# Was downgoing, so must first turn in mantle.
end_action = _ACTIONS["turn"]
self.add_to_branch(tau_model, self.current_branch,
tau_model.cmb_branch - 1, is_p_wave,
end_action)
if next_leg == 'END':
end_action = _ACTIONS["reflect_underside"]
self.add_to_branch(tau_model, self.current_branch,
upgoing_rec_branch, is_p_wave,
end_action)
else:
end_action = _ACTIONS["reflect_underside"]
self.add_to_branch(
tau_model, self.current_branch, 0, is_p_wave,
end_action)
elif next_leg[0] == "v":
discon_branch = closest_branch_to_depth(tau_model,
next_leg[1:])
if self.current_branch <= discon_branch - 1:
end_action = _ACTIONS["reflect_topside"]
self.add_to_branch(tau_model, self.current_branch,
discon_branch - 1, is_p_wave,
end_action)
else:
raise TauModelError(
"Phase not recognised: {} followed by {} when "
"current_branch > discon_branch".format(
current_leg, next_leg))
elif next_leg[0] == "^":
discon_branch = closest_branch_to_depth(tau_model,
next_leg[1:])
if prev_leg == "K":
end_action = _ACTIONS["reflect_underside"]
self.add_to_branch(tau_model, self.current_branch,
discon_branch, is_p_wave,
end_action)
elif prev_leg[0] == "^" or prev_leg in ("P", "S", "p", "s",
"START"):
end_action = _ACTIONS["turn"]
self.add_to_branch(tau_model, self.current_branch,
tau_model.cmb_branch - 1, is_p_wave,
end_action)
end_action = _ACTIONS["reflect_underside"]
self.add_to_branch(
tau_model, self.current_branch, discon_branch,
is_p_wave, end_action)
elif ((prev_leg[0] == "v" and
discon_branch < closest_branch_to_depth(
tau_model, prev_leg[1:]) or
(prev_leg == "m" and
discon_branch < tau_model.moho_branch) or
(prev_leg == "c" and
discon_branch < tau_model.cmb_branch))):
end_action = _ACTIONS["reflect_underside"]
self.add_to_branch(
tau_model, self.current_branch, discon_branch,
is_p_wave, end_action)
else:
raise TauModelError(
"Phase not recognised: {} followed by {} when "
"current_branch > discon_branch".format(
current_leg, next_leg))
elif next_leg == "c":
end_action = _ACTIONS["reflect_topside"]
self.add_to_branch(
tau_model, self.current_branch,
tau_model.cmb_branch - 1, is_p_wave, end_action)
elif next_leg == "K":
end_action = _ACTIONS["transdown"]
self.add_to_branch(
tau_model, self.current_branch,
tau_model.cmb_branch - 1, is_p_wave, end_action)
elif next_leg == "m" or (is_next_leg_depth and
next_leg_depth < tau_model.cmb_depth):
# Treat the Moho in the same way as 410 type
# discontinuities.
discon_branch = closest_branch_to_depth(tau_model,
next_leg)
if end_action == _ACTIONS["turn"] \
or end_action == _ACTIONS["reflect_topside"] \
or end_action == _ACTIONS["transup"]:
# Upgoing section
if discon_branch > self.current_branch:
# Check the discontinuity below the current
# branch when the ray should be upgoing
raise TauModelError(
"Phase not recognised: {} followed by {} when "
"current_branch > discon_branch".format(
current_leg, next_leg))
end_action = _ACTIONS["transup"]
self.add_to_branch(
tau_model, self.current_branch, discon_branch,
is_p_wave, end_action)
else:
# Downgoing section, must look at leg after next to
# determine whether to convert on the downgoing or
# upgoing part of the path.
next_next_leg = self.legs[leg_num + 2]
if next_next_leg == "p" or next_next_leg == "s":
# Convert on upgoing section
end_action = _ACTIONS["turn"]
self.add_to_branch(
tau_model, self.current_branch,
tau_model.cmb_branch - 1, is_p_wave,
end_action)
end_action = _ACTIONS["transup"]
self.add_to_branch(tau_model, self.current_branch,
discon_branch, is_p_wave,
end_action)
elif next_next_leg == "P" or next_next_leg == "S":
if discon_branch > self.current_branch:
# discon is below current loc
end_action = _ACTIONS["transdown"]
self.add_to_branch(
tau_model, self.current_branch,
discon_branch - 1, is_p_wave, end_action)
else:
# Discontinuity is above current location,
# but we have a downgoing ray, so this is an
# illegal ray for this source depth.
self.max_ray_param = -1
return
else:
raise TauModelError(
"Phase not recognized: {} followed by {} "
"followed by {}".format(current_leg, next_leg,
next_next_leg))
else:
raise TauModelError(
"Phase not recognized: {} followed by {}".format(
current_leg, next_leg))
elif current_leg[0] in "PS":
if current_leg == "Pdiff" or current_leg == "Sdiff":
# In the diffracted case we trick addtoBranch into
# thinking we are turning, but then make max_ray_param
# equal to min_ray_param, which is the deepest turning ray.
if (self.max_ray_param >= tau_model.get_tau_branch(
tau_model.cmb_branch - 1,
is_p_wave).min_turn_ray_param >=
self.min_ray_param):
end_action = _ACTIONS["turn"]
self.add_to_branch(tau_model, self.current_branch,
tau_model.cmb_branch - 1, is_p_wave,
end_action)
self.max_ray_param = self.min_ray_param
if next_leg == "END":
end_action = _ACTIONS["reflect_underside"]
self.add_to_branch(tau_model, self.current_branch,
upgoing_rec_branch, is_p_wave,
end_action)
elif next_leg[0] in "PS":
end_action = _ACTIONS["reflect_underside"]
self.add_to_branch(
tau_model, self.current_branch, 0, is_p_wave,
end_action)
else:
# Can't have head wave as ray param is not within
# range.
self.max_ray_param = -1
return
elif current_leg in ("Pg", "Sg", "Pn", "Sn"):
if self.current_branch >= tau_model.moho_branch:
# Pg, Pn, Sg and Sn must be above the moho and so is
# not valid for rays coming upwards from below,
# possibly due to the source depth. Setting
# max_ray_param = -1 effectively disallows this phase.
self.max_ray_param = -1
return
if current_leg in ("Pg", "Sg"):
end_action = _ACTIONS["turn"]
self.add_to_branch(tau_model, self.current_branch,
tau_model.moho_branch - 1,
is_p_wave, end_action)
end_action = _ACTIONS["reflect_underside"]
self.add_to_branch(tau_model, self.current_branch,
upgoing_rec_branch, is_p_wave,
end_action)
elif current_leg in ("Pn", "Sn"):
# In the diffracted case we trick addtoBranch into
# thinking we are turning below the Moho, but then
# make the min_ray_param equal to max_ray_param,
# which is the head wave ray.
if (self.max_ray_param >= tau_model.get_tau_branch(
tau_model.moho_branch,
is_p_wave).max_ray_param >=
self.min_ray_param):
end_action = _ACTIONS["turn"]
self.add_to_branch(
tau_model, self.current_branch,
tau_model.moho_branch, is_p_wave, end_action)
end_action = _ACTIONS["transup"]
self.add_to_branch(
tau_model, self.current_branch,
tau_model.moho_branch, is_p_wave, end_action)
self.min_ray_param = self.max_ray_param
if next_leg == "END":
end_action = _ACTIONS["reflect_underside"]
self.add_to_branch(
tau_model, self.current_branch,
upgoing_rec_branch, is_p_wave, end_action)
elif next_leg[0] in "PS":
end_action = _ACTIONS["reflect_underside"]
self.add_to_branch(
tau_model, self.current_branch, 0,
is_p_wave, end_action)
else:
# Can't have head wave as ray param is not
# within range.
self.max_ray_param = -1
return
else:
raise TauModelError(
"Phase not recognized: {} followed by {}".format(
current_leg, next_leg))
elif current_leg == "K":
if next_leg in ("P", "S"):
if prev_leg in ("P", "S", "K", "k", "START"):
end_action = _ACTIONS["turn"]
self.add_to_branch(
tau_model, self.current_branch,
tau_model.iocb_branch - 1, is_p_wave, end_action)
end_action = _ACTIONS["transup"]
self.add_to_branch(
tau_model, self.current_branch, tau_model.cmb_branch,
is_p_wave, end_action)
elif next_leg == "K":
if prev_leg in ("P", "S", "K"):
end_action = _ACTIONS["turn"]
self.add_to_branch(
tau_model, self.current_branch,
tau_model.iocb_branch - 1, is_p_wave, end_action)
end_action = _ACTIONS["reflect_underside"]
self.add_to_branch(
tau_model, self.current_branch, tau_model.cmb_branch,
is_p_wave, end_action)
elif next_leg in ("I", "J"):
end_action = _ACTIONS["transdown"]
self.add_to_branch(
tau_model, self.current_branch,
tau_model.iocb_branch - 1, is_p_wave, end_action)
elif next_leg == "i":
end_action = _ACTIONS["reflect_topside"]
self.add_to_branch(
tau_model, self.current_branch,
tau_model.iocb_branch - 1, is_p_wave, end_action)
else:
raise TauModelError(
"Phase not recognized: {} followed by {}".format(
current_leg, next_leg))
elif current_leg in ("I", "J"):
end_action = _ACTIONS["turn"]
self.add_to_branch(
tau_model, self.current_branch,
tau_model.tau_branches.shape[1] - 1, is_p_wave, end_action)
if next_leg in ("I", "J"):
end_action = _ACTIONS["reflect_underside"]
self.add_to_branch(
tau_model, self.current_branch, tau_model.iocb_branch,
is_p_wave, end_action)
elif next_leg == "K":
end_action = _ACTIONS["transup"]
self.add_to_branch(
tau_model, self.current_branch, tau_model.iocb_branch,
is_p_wave, end_action)
elif current_leg in ("m", "c", "i") or current_leg[0] == "^":
pass
elif current_leg[0] == "v":
b = closest_branch_to_depth(tau_model, current_leg[1:])
if b == 0:
raise TauModelError(
"Phase not recognized: {} looks like a top side "
"reflection at the free surface.".format(current_leg))
elif is_leg_depth:
# Check for phase like P0s, but could also be P2s if first
# discontinuity is deeper.
b = closest_branch_to_depth(tau_model, current_leg)
if b == 0 and next_leg in ("p", "s"):
raise TauModelError(
"Phase not recognized: {} followed by {} looks like "
"an upgoing wave from the free surface as closest "
"discontinuity to {} is zero depth.".format(
current_leg, next_leg, current_leg))
else:
raise TauModelError(
"Phase not recognized: {} followed by {}".format(
current_leg, next_leg))
if self.max_ray_param != -1:
if (end_action == _ACTIONS["reflect_underside"] and
downgoing_rec_branch == self.branch_seq[-1]):
# Last action was upgoing, so last branch should be
# upgoing_rec_branch
self.min_ray_param = -1
self.max_ray_param = -1
elif (end_action == _ACTIONS["reflect_topside"] and
upgoing_rec_branch == self.branch_seq[-1]):
# Last action was downgoing, so last branch should be
# downgoing_rec_branch
self.min_ray_param = -1
self.max_ray_param = -1
def phase_conversion(self, tau_model, from_branch, end_action, is_p_to_s):
"""
Change max_ray_param and min_ray_param where there is a phase
conversion.
For instance, SKP needs to change the max_ray_param because there are
SKS ray parameters that cannot propagate from the CMB into the mantle
as a P wave.
"""
if end_action == _ACTIONS["turn"]:
# Can't phase convert for just a turn point
raise TauModelError("Bad end_action: phase conversion is not "
"allowed at turn points.")
elif end_action == _ACTIONS["reflect_underside"]:
self.max_ray_param = \
min(self.max_ray_param,
tau_model.get_tau_branch(from_branch,
is_p_to_s).max_ray_param,
tau_model.get_tau_branch(from_branch,
not is_p_to_s).max_ray_param)
elif end_action == _ACTIONS["reflect_topside"]:
self.max_ray_param = min(
self.max_ray_param,
tau_model.get_tau_branch(from_branch,
is_p_to_s).min_turn_ray_param,
tau_model.get_tau_branch(from_branch,
not is_p_to_s).min_turn_ray_param)
elif end_action == _ACTIONS["transup"]:
self.max_ray_param = min(
self.max_ray_param,
tau_model.get_tau_branch(from_branch, is_p_to_s).max_ray_param,
tau_model.get_tau_branch(from_branch - 1,
not is_p_to_s).min_turn_ray_param)
elif end_action == _ACTIONS["transdown"]:
self.max_ray_param = min(
self.max_ray_param,
tau_model.get_tau_branch(from_branch, is_p_to_s).min_ray_param,
tau_model.get_tau_branch(from_branch + 1,
not is_p_to_s).max_ray_param)
else:
raise TauModelError("Illegal end_action = {}".format(end_action))
def add_to_branch(self, tau_model, start_branch, end_branch, is_p_wave,
end_action):
"""
Add branch numbers to branch_seq.
Branches from start_branch to end_branch, inclusive, are added in
order. Also, current_branch is set correctly based on the value of
end_action. end_action can be one of transup, transdown,
reflect_underside, reflect_topside, or turn.
"""
if end_branch < 0 or end_branch > tau_model.tau_branches.shape[1]:
raise ValueError('End branch outside range: %d' % (end_branch, ))
if end_action == _ACTIONS["turn"]:
end_offset = 0
is_down_going = True
self.min_ray_param = max(
self.min_ray_param,
tau_model.get_tau_branch(end_branch,
is_p_wave).min_turn_ray_param)
elif end_action == _ACTIONS["reflect_underside"]:
end_offset = 0
is_down_going = False
self.max_ray_param = min(
self.max_ray_param,
tau_model.get_tau_branch(end_branch, is_p_wave).max_ray_param)
elif end_action == _ACTIONS["reflect_topside"]:
end_offset = 0
is_down_going = True
self.max_ray_param = min(
self.max_ray_param,
tau_model.get_tau_branch(end_branch,
is_p_wave).min_turn_ray_param)
elif end_action == _ACTIONS["transup"]:
end_offset = -1
is_down_going = False
self.max_ray_param = min(
self.max_ray_param,
tau_model.get_tau_branch(end_branch, is_p_wave).max_ray_param)
elif end_action == _ACTIONS["transdown"]:
end_offset = 1
is_down_going = True
self.max_ray_param = min(
self.max_ray_param,
tau_model.get_tau_branch(end_branch, is_p_wave).min_ray_param)
else:
raise TauModelError("Illegal end_action: {}".format(end_action))
if is_down_going:
if start_branch > end_branch:
# Can't be downgoing as we are already below.
self.min_ray_param = -1
self.max_ray_param = -1
else:
# Must be downgoing, so increment i.
for i in range(start_branch, end_branch + 1):
self.branch_seq.append(i)
self.down_going.append(is_down_going)
self.wave_type.append(is_p_wave)
else:
if start_branch < end_branch:
# Can't be upgoing as we are already above.
self.min_ray_param = -1
self.max_ray_param = -1
else:
# Upgoing, so decrement i.
for i in range(start_branch, end_branch - 1, -1):
self.branch_seq.append(i)
self.down_going.append(is_down_going)
self.wave_type.append(is_p_wave)
self.current_branch = end_branch + end_offset
def sum_branches(self, tau_model):
"""Sum the appropriate branches for this phase."""
# Special case for surface waves.
if self.name.endswith("kmps"):
self.dist = np.zeros(2)
self.time = np.zeros(2)
self.ray_param = np.empty(2)
self.ray_param[0] = \
tau_model.radius_of_planet / float(self.name[:-4])
self.dist[1] = 2 * math.pi
self.time[1] = \
2 * math.pi * tau_model.radius_of_planet / \
float(self.name[:-4])
self.ray_param[1] = self.ray_param[0]
self.min_distance = 0
self.max_distance = 2 * math.pi
self.down_going.append(True)
return
if self.max_ray_param < 0 or self.min_ray_param > self.max_ray_param:
# Phase has no arrivals, possibly due to source depth.
self.ray_param = np.empty(0)
self.min_ray_param = -1
self.max_ray_param = -1
self.dist = np.empty(0)
self.time = np.empty(0)
self.max_distance = -1
return
# Find the ray parameter index that corresponds to the min_ray_param
# and max_ray_param.
index = np.where(tau_model.ray_params >= self.min_ray_param)[0]
if len(index):
self.min_ray_param_index = index[-1]
index = np.where(tau_model.ray_params >= self.max_ray_param)[0]
if len(index):
self.max_ray_param_index = index[-1]
if self.max_ray_param_index == 0 \
and self.min_ray_param_index == len(tau_model.ray_params) - 1:
# All ray parameters are valid so just copy:
self.ray_param = tau_model.ray_param.copy()
elif self.max_ray_param_index == self.min_ray_param_index:
# if "Sdiff" in self.name or "Pdiff" in self.name:
# self.ray_param = [self.min_ray_param, self.min_ray_param]
# elif "Pn" in self.name or "Sn" in self.name:
# self.ray_param = [self.min_ray_param, self.min_ray_param]
if self.name.endswith("kmps"):
self.ray_param = np.array([0, self.max_ray_param])
else:
self.ray_param = np.array([self.min_ray_param,
self.min_ray_param])
else:
# Only a subset of the ray parameters is valid so use these.
self.ray_param = \
tau_model.ray_params[self.max_ray_param_index:
self.min_ray_param_index + 1].copy()
self.dist = np.zeros(shape=self.ray_param.shape)
self.time = np.zeros(shape=self.ray_param.shape)
# Counter for passes through each branch. 0 is P and 1 is S.
times_branches = self.calc_branch_mult(tau_model)
# Sum the branches with the appropriate multiplier.
size = self.min_ray_param_index - self.max_ray_param_index + 1
index = slice(self.max_ray_param_index, self.min_ray_param_index + 1)
for i in range(tau_model.tau_branches.shape[1]):
tb = times_branches[0, i]
tbs = times_branches[1, i]
taub = tau_model.tau_branches[0, i]
taubs = tau_model.tau_branches[1, i]
if tb != 0:
self.dist[:size] += tb * taub.dist[index]
self.time[:size] += tb * taub.time[index]
if tbs != 0:
self.dist[:size] += tbs * taubs.dist[index]
self.time[:size] += tbs * taubs.time[index]
if "Sdiff" in self.name or "Pdiff" in self.name:
if tau_model.s_mod.depth_in_high_slowness(
tau_model.cmb_depth - 1e-10, self.min_ray_param,
self.name[0] == "P"):
# No diffraction if there is a high slowness zone at the CMB.
self.min_ray_param = -1
self.max_ray_param = -1
self.max_distance = -1
self.time = np.empty(0)
self.dist = np.empty(0)
self.ray_param = np.empty(0)
return
else:
self.dist[1] = self.dist[0] + \
self._settings["max_diffraction_in_radians"]
self.time[1] = (
self.time[0] +
self._settings["max_diffraction_in_radians"] *
self.min_ray_param)
elif "Pn" in self.name or "Sn" in self.name:
self.dist[1] = self.dist[0] + \
self._settings["max_refraction_in_radians"]
self.time[1] = (
self.time[0] + self._settings["max_refraction_in_radians"] *
self.min_ray_param)
elif self.max_ray_param_index == self.min_ray_param_index:
self.dist[1] = self.dist[0]
self.time[1] = self.time[0]
self.min_distance = np.min(self.dist)
self.max_distance = np.max(self.dist)
# Now check to see if our ray parameter range includes any ray
# parameters that are associated with high slowness zones. If so,
# then we will need to insert a "shadow zone" into our time and
# distance arrays. It is represented by a repeated ray parameter.
for is_pwave in [True, False]:
hsz = tau_model.s_mod.high_slowness_layer_depths_p \
if is_pwave \
else tau_model.s_mod.high_slowness_layer_depths_s
index_offset = 0
for hszi in hsz:
if self.max_ray_param > hszi.ray_param > self.min_ray_param:
# There is a high slowness zone within our ray parameter
# range so might need to add a shadow zone. Need to
# check if the current wave type is part of the phase at
# this depth/ray parameter.
branch_num = tau_model.find_branch(hszi.top_depth)
found_overlap = False
for leg_num in range(len(self.branch_seq)):
# Check for downgoing legs that cross the high
# slowness zone with the same wave type.
if (self.branch_seq[leg_num] == branch_num and
self.wave_type[leg_num] == is_pwave and
self.down_going[leg_num] is True and
self.branch_seq[leg_num - 1] ==
branch_num - 1 and
self.wave_type[leg_num - 1] == is_pwave and
self.down_going[leg_num - 1] is True):
found_overlap = True
break
if found_overlap:
hsz_index = np.where(self.ray_param == hszi.ray_param)
hsz_index = hsz_index[0][0]
newlen = self.ray_param.shape[0] + 1
new_ray_params = np.empty(newlen)
newdist = np.empty(newlen)
newtime = np.empty(newlen)
new_ray_params[:hsz_index] = self.ray_param[:hsz_index]
newdist[:hsz_index] = self.dist[:hsz_index]
newtime[:hsz_index] = self.time[:hsz_index]
# Sum the branches with an appropriate multiplier.
new_ray_params[hsz_index] = hszi.ray_param
newdist[hsz_index] = 0
newtime[hsz_index] = 0
for tb, tbs, taub, taubs in zip(
times_branches[0], times_branches[1],
tau_model.tau_branches[0],
tau_model.tau_branches[1]):
if tb != 0 and taub.top_depth < hszi.top_depth:
newdist[hsz_index] += tb * taub.dist[
self.max_ray_param_index + hsz_index -
index_offset]
newtime[hsz_index] += tb * taub.time[
self.max_ray_param_index + hsz_index -
index_offset]
if tbs != 0 and taubs.top_depth < hszi.top_depth:
newdist[hsz_index] += tbs * taubs.dist[
self.max_ray_param_index + hsz_index -
index_offset]
newtime[hsz_index] += tbs * taubs.time[
self.max_ray_param_index + hsz_index -
index_offset]
newdist[hsz_index + 1:] = self.dist[hsz_index:]
newtime[hsz_index + 1:] = self.time[hsz_index:]
new_ray_params[hsz_index + 1:] = \
self.ray_param[hsz_index:]
index_offset += 1
self.dist = newdist
self.time = newtime
self.ray_param = new_ray_params
def calc_branch_mult(self, tau_model):
"""
Calculate how many times the phase passes through a branch, up or down.
With this result, we can just multiply instead of doing the ray calc
for each time.
"""
# Initialise the counter for each branch to 0. 0 is P and 1 is S.
times_branches = np.zeros((2, tau_model.tau_branches.shape[1]))
# Count how many times each branch appears in the path.