1# http://pyrocko.org - GPLv3 

2# 

3# The Pyrocko Developers, 21st Century 

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

5 

6''' 

7QuakeML 1.2 input, output, and data model. 

8 

9This modules provides support to read and write `QuakeML version 1.2 

10<https://quake.ethz.ch/quakeml>`_. It defines a hierarchy of Python objects, 

11closely following the QuakeML data model. 

12 

13QuakeML is a flexible, extensible and modular XML representation of 

14seismological data which is intended to cover a broad range of fields of 

15application in modern seismology. It covers basic seismic event description, 

16including moment tensors. 

17 

18For convenience and ease of use, this documentation contains excerpts from the 

19`QuakeML Manual 

20<https://quake.ethz.ch/quakeml/docs/REC?action=AttachFile&do=get&target=QuakeML-BED-20130214b.pdf>`_. 

21However, this information is not intended to be complete. Please refer to the 

22QuakeML Manual for details. 

23''' 

24 

25from __future__ import absolute_import 

26import logging 

27from pyrocko.guts import StringPattern, StringChoice, String, Float, Int,\ 

28 Timestamp, Object, List, Union, Bool, Unicode 

29from pyrocko.model import event 

30from pyrocko.gui import marker 

31from pyrocko import moment_tensor 

32import numpy as num 

33 

34logger = logging.getLogger('pyrocko.io.quakeml') 

35 

36guts_prefix = 'quakeml' 

37guts_xmlns = 'http://quakeml.org/xmlns/bed/1.2' 

38polarity_choices = {'positive': 1, 'negative': -1, 'undecidable': None} 

39 

40 

41class QuakeMLError(Exception): 

42 pass 

43 

44 

45class NoPreferredOriginSet(QuakeMLError): 

46 pass 

47 

48 

49def one_element_or_none(li): 

50 if len(li) == 1: 

51 return li[0] 

52 elif len(li) == 0: 

53 return None 

54 else: 

55 logger.warning('More than one element in list: {}'.format(li)) 

56 return None 

57 

58 

59class ResourceIdentifier(StringPattern): 

60 ''' 

61 Identifies resource origin. 

62 

63 They consist of an authority identifier, a unique resource identifier, and 

64 an optional local identifier. The URI schema name smi stands for 

65 seismological meta-information, thus indicating a connection to a set of 

66 metadata associated with the resource. 

67 ''' 

68 pattern = "^(smi|quakeml):[\\w\\d][\\w\\d\\-\\.\\*\\(\\)_~']{2,}/[\\w" +\ 

69 "\\d\\-\\.\\*\\(\\)_~'][\\w\\d\\-\\.\\*\\(\\)\\+\\?_~'=,;#/&]*$" 

70 

71 

72class WhitespaceOrEmptyStringType(StringPattern): 

73 pattern = '^\\s*$' 

74 

75 

76class OriginUncertaintyDescription(StringChoice): 

77 ''' 

78 Preferred uncertainty description. 

79 ''' 

80 choices = [ 

81 'horizontal uncertainty', 

82 'uncertainty ellipse', 

83 'confidence ellipsoid'] 

84 

85 

86class AmplitudeCategory(StringChoice): 

87 ''' 

88 Description of the way the waveform trace is evaluated to get an amplitude 

89 value. 

90 

91 This can be just reading a single value for a given point in time (point), 

92 taking a mean value over a time interval (mean), integrating the trace 

93 over a time interval (integral), specifying just a time interval 

94 (duration), or evaluating a period (period). 

95 ''' 

96 choices = ['point', 'mean', 'duration', 'period', 'integral', 'other'] 

97 

98 

99class OriginDepthType(StringChoice): 

100 ''' 

101 Type of depth determination. 

102 ''' 

103 choices = [ 

104 'from location', 

105 'from moment tensor inversion', 

106 'from modeling of broad-band P waveforms', 

107 'constrained by depth phases', 

108 'constrained by direct phases', 

109 'constrained by depth and direct phases', 

110 'operator assigned', 

111 'other'] 

112 

113 

114class OriginType(StringChoice): 

115 ''' 

116 Describes the origin type. 

117 ''' 

118 choices = [ 

119 'hypocenter', 

120 'centroid', 

121 'amplitude', 

122 'macroseismic', 

123 'rupture start', 

124 'rupture end'] 

125 

126 

127class MTInversionType(StringChoice): 

128 ''' 

129 Type of moment tensor inversion. 

130 

131 Users should avoid to give contradictory information in 

132 :py:class:`MTInversionType` and :py:gattr:`MomentTensor.method_id`. 

133 ''' 

134 choices = ['general', 'zero trace', 'double couple'] 

135 

136 

137class EvaluationMode(StringChoice): 

138 ''' 

139 Mode of an evaluation. 

140 

141 Used in :py:class:`Pick`, :py:class:`Amplitude`, :py:class:`Magnitude`, 

142 :py:class:`Origin`, :py:class:`FocalMechanism`. 

143 ''' 

144 choices = ['manual', 'automatic'] 

145 

146 

147class EvaluationStatus(StringChoice): 

148 ''' 

149 Status of an evaluation. 

150 

151 Used in :py:class:`Pick`, :py:class:`Amplitude`, :py:class:`Magnitude`, 

152 :py:class:`Origin`, :py:class:`FocalMechanism`. 

153 ''' 

154 choices = ['preliminary', 'confirmed', 'reviewed', 'final', 'rejected'] 

155 

156 

157class PickOnset(StringChoice): 

158 ''' 

159 Flag that roughly categorizes the sharpness of the onset. 

160 ''' 

161 choices = ['emergent', 'impulsive', 'questionable'] 

162 

163 

164class EventType(StringChoice): 

165 ''' 

166 Describes the type of an event. 

167 ''' 

168 choices = [ 

169 'not existing', 

170 'not reported', 

171 'earthquake', 

172 'anthropogenic event', 

173 'collapse', 

174 'cavity collapse', 

175 'mine collapse', 

176 'building collapse', 

177 'explosion', 

178 'accidental explosion', 

179 'chemical explosion', 

180 'controlled explosion', 

181 'experimental explosion', 

182 'industrial explosion', 

183 'mining explosion', 

184 'quarry blast', 

185 'road cut', 

186 'blasting levee', 

187 'nuclear explosion', 

188 'induced or triggered event', 

189 'rock burst', 

190 'reservoir loading', 

191 'fluid injection', 

192 'fluid extraction', 

193 'crash', 

194 'plane crash', 

195 'train crash', 

196 'boat crash', 

197 'other event', 

198 'atmospheric event', 

199 'sonic boom', 

200 'sonic blast', 

201 'acoustic noise', 

202 'thunder', 

203 'avalanche', 

204 'snow avalanche', 

205 'debris avalanche', 

206 'hydroacoustic event', 

207 'ice quake', 

208 'slide', 

209 'landslide', 

210 'rockslide', 

211 'meteorite', 

212 'volcanic eruption', 

213 'duplicate earthquake', 

214 'rockburst'] 

215 

216 

217class DataUsedWaveType(StringChoice): 

218 ''' 

219 Type of waveform data. 

220 ''' 

221 choices = [ 

222 'P waves', 

223 'body waves', 

224 'surface waves', 

225 'mantle waves', 

226 'combined', 

227 'unknown'] 

228 

229 

230class AmplitudeUnit(StringChoice): 

231 ''' 

232 Provides the most likely measurement units. 

233 

234 The measurement units for physical quantity are described in the 

235 :py:gattr:`Amplitude.generic_amplitude` attribute. Possible values are 

236 specified as combination of SI base units. 

237 ''' 

238 choices = ['m', 's', 'm/s', 'm/(s*s)', 'm*s', 'dimensionless', 'other'] 

239 

240 

241class EventDescriptionType(StringChoice): 

242 ''' 

243 Category of earthquake description. 

244 ''' 

245 choices = [ 

246 'felt report', 

247 'Flinn-Engdahl region', 

248 'local time', 

249 'tectonic summary', 

250 'nearest cities', 

251 'earthquake name', 

252 'region name'] 

253 

254 

255class MomentTensorCategory(StringChoice): 

256 ''' 

257 Category of moment tensor inversion. 

258 ''' 

259 choices = ['teleseismic', 'regional'] 

260 

261 

262class EventTypeCertainty(StringChoice): 

263 ''' 

264 Denotes how certain the information on event type is. 

265 ''' 

266 choices = ['known', 'suspected'] 

267 

268 

269class SourceTimeFunctionType(StringChoice): 

270 ''' 

271 Type of source time function. 

272 ''' 

273 choices = ['box car', 'triangle', 'trapezoid', 'unknown'] 

274 

275 

276class PickPolarity(StringChoice): 

277 choices = list(polarity_choices.keys()) 

278 

279 

280class AgencyID(String): 

281 pass 

282 

283 

284class Author(Unicode): 

285 pass 

286 

287 

288class Version(String): 

289 pass 

290 

291 

292class Phase(Object): 

293 value = String.T(xmlstyle='content') 

294 

295 

296class GroundTruthLevel(String): 

297 pass 

298 

299 

300class AnonymousNetworkCode(String): 

301 pass 

302 

303 

304class AnonymousStationCode(String): 

305 pass 

306 

307 

308class AnonymousChannelCode(String): 

309 pass 

310 

311 

312class AnonymousLocationCode(String): 

313 pass 

314 

315 

316class Type(String): 

317 pass 

318 

319 

320class MagnitudeHint(String): 

321 pass 

322 

323 

324class Region(Unicode): 

325 pass 

326 

327 

328class RealQuantity(Object): 

329 value = Float.T() 

330 uncertainty = Float.T(optional=True) 

331 lower_uncertainty = Float.T(optional=True) 

332 upper_uncertainty = Float.T(optional=True) 

333 confidence_level = Float.T(optional=True) 

334 

335 

336class IntegerQuantity(Object): 

337 value = Int.T() 

338 uncertainty = Int.T(optional=True) 

339 lower_uncertainty = Int.T(optional=True) 

340 upper_uncertainty = Int.T(optional=True) 

341 confidence_level = Float.T(optional=True) 

342 

343 

344class ConfidenceEllipsoid(Object): 

345 ''' 

346 Representation of the location uncertainty as a confidence ellipsoid with 

347 arbitrary orientation in space. 

348 ''' 

349 semi_major_axis_length = Float.T() 

350 semi_minor_axis_length = Float.T() 

351 semi_intermediate_axis_length = Float.T() 

352 major_axis_plunge = Float.T() 

353 major_axis_azimuth = Float.T() 

354 major_axis_rotation = Float.T() 

355 

356 

357class TimeQuantity(Object): 

358 ''' 

359 Representation of a point in time. 

360 

361 It's given in ISO 8601 format, with optional symmetric or asymmetric 

362 uncertainties given in seconds. The time has to be specified in UTC. 

363 ''' 

364 value = Timestamp.T() 

365 uncertainty = Float.T(optional=True) 

366 lower_uncertainty = Float.T(optional=True) 

367 upper_uncertainty = Float.T(optional=True) 

368 confidence_level = Float.T(optional=True) 

369 

370 

371class TimeWindow(Object): 

372 ''' 

373 Representation of a time window for amplitude measurements. 

374 

375 Which is given by a central point in time, and points in time before and 

376 after this central point. Both points before and after may coincide with 

377 the central point. 

378 ''' 

379 begin = Float.T() 

380 end = Float.T() 

381 reference = Timestamp.T() 

382 

383 

384class ResourceReference(ResourceIdentifier): 

385 ''' 

386 This type is used to refer to QuakeML resources as described in Sect. 3.1 

387 in the `QuakeML manual <https://quake.ethz.ch/quakeml/docs/REC?action=Att\ 

388 achFile&do=get&target=QuakeML-BED-20130214b.pdf>`_. 

389 ''' 

390 pass 

391 

392 

393class DataUsed(Object): 

394 ''' 

395 Description of the type of data used in a moment-tensor inversion. 

396 ''' 

397 wave_type = DataUsedWaveType.T() 

398 station_count = Int.T(optional=True) 

399 component_count = Int.T(optional=True) 

400 shortest_period = Float.T(optional=True) 

401 longest_period = Float.T(optional=True) 

402 

403 

404class EventDescription(Object): 

405 ''' 

406 Free-form string with additional event description. 

407 

408 This can be a well-known name, like 1906 San Francisco Earthquake. 

409 A number of categories can be given in :py:gattr:`type`. 

410 ''' 

411 text = Unicode.T() 

412 type = EventDescriptionType.T(optional=True) 

413 

414 

415class SourceTimeFunction(Object): 

416 ''' 

417 Source time function used in moment-tensor inversion. 

418 ''' 

419 type = SourceTimeFunctionType.T() 

420 duration = Float.T() 

421 rise_time = Float.T(optional=True) 

422 decay_time = Float.T(optional=True) 

423 

424 

425class OriginQuality(Object): 

426 ''' 

427 Description of an origin's quality. 

428 

429 It contains various attributes commonly used to describe the 

430 quality of an origin, e. g., errors, azimuthal coverage, etc. 

431 :py:class:`Origin` objects have an optional attribute of the type 

432 :py:gattr:`OriginQuality`. 

433 ''' 

434 associated_phase_count = Int.T(optional=True) 

435 used_phase_count = Int.T(optional=True) 

436 associated_station_count = Int.T(optional=True) 

437 used_station_count = Int.T(optional=True) 

438 depth_phase_count = Int.T(optional=True) 

439 standard_error = Float.T(optional=True) 

440 azimuthal_gap = Float.T(optional=True) 

441 secondary_azimuthal_gap = Float.T(optional=True) 

442 ground_truth_level = GroundTruthLevel.T(optional=True) 

443 maximum_distance = Float.T(optional=True) 

444 minimum_distance = Float.T(optional=True) 

445 median_distance = Float.T(optional=True) 

446 

447 

448class Axis(Object): 

449 ''' 

450 Representation of an eigenvector of a moment tensor. 

451 

452 Which is expressed in its principal-axes system and uses the angles 

453 :py:gattr:`azimuth`, :py:gattr:`plunge`, and the eigenvalue 

454 :py:gattr:`length`. 

455 ''' 

456 azimuth = RealQuantity.T() 

457 plunge = RealQuantity.T() 

458 length = RealQuantity.T() 

459 

460 

461class Tensor(Object): 

462 ''' 

463 Representation of the six independent moment-tensor elements in spherical 

464 coordinates. 

465 

466 These are the moment-tensor elements Mrr, Mtt, Mpp, Mrt, Mrp, Mtp in the 

467 spherical coordinate system defined by local upward vertical (r), 

468 North-South (t), and West-East (p) directions. 

469 ''' 

470 mrr = RealQuantity.T(xmltagname='Mrr') 

471 mtt = RealQuantity.T(xmltagname='Mtt') 

472 mpp = RealQuantity.T(xmltagname='Mpp') 

473 mrt = RealQuantity.T(xmltagname='Mrt') 

474 mrp = RealQuantity.T(xmltagname='Mrp') 

475 mtp = RealQuantity.T(xmltagname='Mtp') 

476 

477 

478class NodalPlane(Object): 

479 ''' 

480 Description of a nodal plane of a focal mechanism. 

481 ''' 

482 strike = RealQuantity.T() 

483 dip = RealQuantity.T() 

484 rake = RealQuantity.T() 

485 

486 

487class CompositeTime(Object): 

488 ''' 

489 Representation of a time instant. 

490 

491 If the specification is given with no greater accuracy than days (i.e., no 

492 time components are given), the date refers to local time. However, if 

493 time components are given, they have to refer to UTC. 

494 ''' 

495 year = IntegerQuantity.T(optional=True) 

496 month = IntegerQuantity.T(optional=True) 

497 day = IntegerQuantity.T(optional=True) 

498 hour = IntegerQuantity.T(optional=True) 

499 minute = IntegerQuantity.T(optional=True) 

500 second = RealQuantity.T(optional=True) 

501 

502 

503class OriginUncertainty(Object): 

504 ''' 

505 Description of the location uncertainties of an origin. 

506 

507 The uncertainty can be described either as a simple circular horizontal 

508 uncertainty, an uncertainty ellipse according to IMS1.0, or a confidence 

509 ellipsoid. If multiple uncertainty models are given, the preferred variant 

510 can be specified in the attribute :py:gattr:`preferred_description`. 

511 ''' 

512 horizontal_uncertainty = Float.T(optional=True) 

513 min_horizontal_uncertainty = Float.T(optional=True) 

514 max_horizontal_uncertainty = Float.T(optional=True) 

515 azimuth_max_horizontal_uncertainty = Float.T(optional=True) 

516 confidence_ellipsoid = ConfidenceEllipsoid.T(optional=True) 

517 preferred_description = OriginUncertaintyDescription.T(optional=True) 

518 confidence_level = Float.T(optional=True) 

519 

520 

521class ResourceReferenceOptional(Union): 

522 members = [ResourceReference.T(), WhitespaceOrEmptyStringType.T()] 

523 

524 

525class CreationInfo(Object): 

526 ''' 

527 Description of creation metadata (author, version, and creation time) of a 

528 resource. 

529 ''' 

530 agency_id = AgencyID.T(optional=True, xmltagname='agencyID') 

531 agency_uri = ResourceReference.T(optional=True, xmltagname='agencyURI') 

532 author = Author.T(optional=True) 

533 author_uri = ResourceReference.T(optional=True, xmltagname='authorURI') 

534 creation_time = Timestamp.T(optional=True) 

535 version = Version.T(optional=True) 

536 

537 

538class StationMagnitudeContribution(Object): 

539 ''' 

540 Description of the weighting of magnitude values from several 

541 :py:class:`StationMagnitude` objects for computing a network magnitude 

542 estimation. 

543 ''' 

544 station_magnitude_id = ResourceReference.T(xmltagname='stationMagnitudeID') 

545 residual = Float.T(optional=True) 

546 weight = Float.T(optional=True) 

547 

548 

549class PrincipalAxes(Object): 

550 ''' 

551 Representation of the principal axes of a moment tensor. 

552 

553 :py:gattr:`t_axis` and :py:gattr:`p_axis` are required, while 

554 :py:gattr:`n_axis` is optional. 

555 ''' 

556 t_axis = Axis.T() 

557 p_axis = Axis.T() 

558 n_axis = Axis.T(optional=True) 

559 

560 

561class NodalPlanes(Object): 

562 ''' 

563 Representation of the nodal planes of a moment tensor. 

564 

565 The attribute :py:gattr:`preferred_plane` can be used to define which plane 

566 is the preferred one. 

567 ''' 

568 preferred_plane = Int.T(optional=True, xmlstyle='attribute') 

569 nodal_plane1 = NodalPlane.T(optional=True) 

570 nodal_plane2 = NodalPlane.T(optional=True) 

571 

572 

573class WaveformStreamID(Object): 

574 ''' 

575 Reference to a stream description in an inventory. 

576 

577 This is mostly equivalent to the combination of networkCode, stationCode, 

578 locationCode, and channelCode. However, additional information, e.g., 

579 sampling rate, can be referenced by the resourceURI. It is recommended to 

580 use resourceURI as a flexible, abstract, and unique stream ID that allows 

581 to describe different processing levels, or resampled/filtered products of 

582 the same initial stream, without violating the intrinsic meaning of the 

583 legacy identifiers (:py:gattr:`network_code`, :py:gattr:`station_code`, 

584 :py:gattr:`channel_code`, and :py:gattr:`location_code`). However, for 

585 operation in the context of legacy systems, the classical identifier 

586 components are supported. 

587 ''' 

588 value = ResourceReferenceOptional.T(xmlstyle='content') 

589 network_code = AnonymousNetworkCode.T(xmlstyle='attribute') 

590 station_code = AnonymousStationCode.T(xmlstyle='attribute') 

591 channel_code = AnonymousChannelCode.T(optional=True, xmlstyle='attribute') 

592 location_code = AnonymousLocationCode.T( 

593 optional=True, xmlstyle='attribute') 

594 

595 @property 

596 def nslc_id(self): 

597 return (self.network_code, self.station_code, self.location_code, 

598 self.channel_code) 

599 

600 

601class Comment(Object): 

602 ''' 

603 Comment to a resource together with author and creation time information. 

604 ''' 

605 id = ResourceReference.T(optional=True, xmlstyle='attribute') 

606 text = Unicode.T() 

607 creation_info = CreationInfo.T(optional=True) 

608 

609 

610class MomentTensor(Object): 

611 ''' 

612 Representation of a moment tensor solution for an event. 

613 

614 It is an optional part of a :py:class:`FocalMechanism` description. 

615 ''' 

616 public_id = ResourceReference.T( 

617 xmlstyle='attribute', xmltagname='publicID') 

618 data_used_list = List.T(DataUsed.T()) 

619 comment_list = List.T(Comment.T()) 

620 derived_origin_id = ResourceReference.T( 

621 optional=True, xmltagname='derivedOriginID') 

622 moment_magnitude_id = ResourceReference.T( 

623 optional=True, xmltagname='momentMagnitudeID') 

624 scalar_moment = RealQuantity.T(optional=True) 

625 tensor = Tensor.T(optional=True) 

626 variance = Float.T(optional=True) 

627 variance_reduction = Float.T(optional=True) 

628 double_couple = Float.T(optional=True) 

629 clvd = Float.T(optional=True) 

630 iso = Float.T(optional=True) 

631 greens_function_id = ResourceReference.T( 

632 optional=True, xmltagname='greensFunctionID') 

633 filter_id = ResourceReference.T(optional=True, xmltagname='filterID') 

634 source_time_function = SourceTimeFunction.T(optional=True) 

635 method_id = ResourceReference.T(optional=True, xmltagname='methodID') 

636 category = MomentTensorCategory.T(optional=True) 

637 inversion_type = MTInversionType.T(optional=True) 

638 creation_info = CreationInfo.T(optional=True) 

639 

640 def pyrocko_moment_tensor(self): 

641 mrr = self.tensor.mrr.value 

642 mtt = self.tensor.mtt.value 

643 mpp = self.tensor.mpp.value 

644 mrt = self.tensor.mrt.value 

645 mrp = self.tensor.mrp.value 

646 mtp = self.tensor.mtp.value 

647 mt = moment_tensor.MomentTensor(m_up_south_east=num.matrix([ 

648 [mrr, mrt, mrp], [mrt, mtt, mtp], [mrp, mtp, mpp]])) 

649 

650 return mt 

651 

652 

653class Amplitude(Object): 

654 ''' 

655 Quantification of the waveform anomaly. 

656 

657 Usually it consists of a single amplitude measurement or a measurement of 

658 the visible signal duration for duration magnitudes. 

659 ''' 

660 public_id = ResourceReference.T( 

661 xmlstyle='attribute', xmltagname='publicID') 

662 comment_list = List.T(Comment.T()) 

663 generic_amplitude = RealQuantity.T() 

664 type = Type.T(optional=True) 

665 category = AmplitudeCategory.T(optional=True) 

666 unit = AmplitudeUnit.T(optional=True) 

667 method_id = ResourceReference.T(optional=True, xmltagname='methodID') 

668 period = RealQuantity.T(optional=True) 

669 snr = Float.T(optional=True) 

670 time_window = TimeWindow.T(optional=True) 

671 pick_id = ResourceReference.T(optional=True, xmltagname='pickID') 

672 waveform_id = WaveformStreamID.T(optional=True, xmltagname='waveformID') 

673 filter_id = ResourceReference.T(optional=True, xmltagname='filterID') 

674 scaling_time = TimeQuantity.T(optional=True) 

675 magnitude_hint = MagnitudeHint.T(optional=True) 

676 evaluation_mode = EvaluationMode.T(optional=True) 

677 evaluation_status = EvaluationStatus.T(optional=True) 

678 creation_info = CreationInfo.T(optional=True) 

679 

680 

681class Magnitude(Object): 

682 ''' 

683 Description of a magnitude value. 

684 

685 It can, but does not need to be associated with an origin. Association 

686 with an origin is expressed with the optional attribute 

687 :py:gattr:`origin_id`. It is either a combination of different magnitude 

688 estimations, or it represents the reported magnitude for the given event. 

689 ''' 

690 public_id = ResourceReference.T( 

691 xmlstyle='attribute', xmltagname='publicID') 

692 comment_list = List.T(Comment.T()) 

693 station_magnitude_contribution_list = List.T( 

694 StationMagnitudeContribution.T()) 

695 mag = RealQuantity.T() 

696 type = Type.T(optional=True) 

697 origin_id = ResourceReference.T(optional=True, xmltagname='originID') 

698 method_id = ResourceReference.T(optional=True, xmltagname='methodID') 

699 station_count = Int.T(optional=True) 

700 azimuthal_gap = Float.T(optional=True) 

701 evaluation_mode = EvaluationMode.T(optional=True) 

702 evaluation_status = EvaluationStatus.T(optional=True) 

703 creation_info = CreationInfo.T(optional=True) 

704 

705 

706class StationMagnitude(Object): 

707 ''' 

708 Description of a magnitude derived from a single waveform stream. 

709 ''' 

710 public_id = ResourceReference.T( 

711 xmlstyle='attribute', xmltagname='publicID') 

712 comment_list = List.T(Comment.T()) 

713 origin_id = ResourceReference.T(optional=True, xmltagname='originID') 

714 mag = RealQuantity.T() 

715 type = Type.T(optional=True) 

716 amplitude_id = ResourceReference.T(optional=True, xmltagname='amplitudeID') 

717 method_id = ResourceReference.T(optional=True, xmltagname='methodID') 

718 waveform_id = WaveformStreamID.T(optional=True, xmltagname='waveformID') 

719 creation_info = CreationInfo.T(optional=True) 

720 

721 

722class Arrival(Object): 

723 ''' 

724 Successful association of a pick with an origin qualifies this pick as an 

725 arrival. 

726 

727 An arrival thus connects a pick with an origin and provides 

728 additional attributes that describe this relationship. Usually 

729 qualification of a pick as an arrival for a given origin is a hypothesis, 

730 which is based on assumptions about the type of arrival (phase) as well as 

731 observed and (on the basis of an earth model) computed arrival times, or 

732 the residual, respectively. Additional pick attributes like the horizontal 

733 slowness and backazimuth of the observed wave—especially if derived from 

734 array data—may further constrain the nature of the arrival. 

735 ''' 

736 public_id = ResourceReference.T( 

737 xmlstyle='attribute', xmltagname='publicID') 

738 comment_list = List.T(Comment.T()) 

739 pick_id = ResourceReference.T(xmltagname='pickID') 

740 phase = Phase.T() 

741 time_correction = Float.T(optional=True) 

742 azimuth = Float.T(optional=True) 

743 distance = Float.T(optional=True) 

744 takeoff_angle = RealQuantity.T(optional=True) 

745 time_residual = Float.T(optional=True) 

746 horizontal_slowness_residual = Float.T(optional=True) 

747 backazimuth_residual = Float.T(optional=True) 

748 time_weight = Float.T(optional=True) 

749 time_used = Int.T(optional=True) 

750 horizontal_slowness_weight = Float.T(optional=True) 

751 backazimuth_weight = Float.T(optional=True) 

752 earth_model_id = ResourceReference.T( 

753 optional=True, xmltagname='earthModelID') 

754 creation_info = CreationInfo.T(optional=True) 

755 

756 

757class Pick(Object): 

758 ''' 

759 A pick is the observation of an amplitude anomaly in a seismogram at a 

760 specific point in time. 

761 

762 It is not necessarily related to a seismic event. 

763 ''' 

764 public_id = ResourceReference.T( 

765 xmlstyle='attribute', xmltagname='publicID') 

766 comment_list = List.T(Comment.T()) 

767 time = TimeQuantity.T() 

768 waveform_id = WaveformStreamID.T(xmltagname='waveformID') 

769 filter_id = ResourceReference.T(optional=True, xmltagname='filterID') 

770 method_id = ResourceReference.T(optional=True, xmltagname='methodID') 

771 horizontal_slowness = RealQuantity.T(optional=True) 

772 backazimuth = RealQuantity.T(optional=True) 

773 slowness_method_id = ResourceReference.T( 

774 optional=True, xmltagname='slownessMethodID') 

775 onset = PickOnset.T(optional=True) 

776 phase_hint = Phase.T(optional=True) 

777 polarity = PickPolarity.T(optional=True) 

778 evaluation_mode = EvaluationMode.T(optional=True) 

779 evaluation_status = EvaluationStatus.T(optional=True) 

780 creation_info = CreationInfo.T(optional=True) 

781 

782 @property 

783 def pyrocko_polarity(self): 

784 return polarity_choices.get(self.polarity, None) 

785 

786 def get_pyrocko_phase_marker(self, event=None): 

787 if not self.phase_hint: 

788 logger.warn('Pick %s: phase_hint undefined' % self.public_id) 

789 phasename = 'undefined' 

790 else: 

791 phasename = self.phase_hint.value 

792 

793 return marker.PhaseMarker( 

794 event=event, nslc_ids=[self.waveform_id.nslc_id], 

795 tmin=self.time.value, tmax=self.time.value, 

796 phasename=phasename, 

797 polarity=self.pyrocko_polarity, 

798 automatic=self.evaluation_mode) 

799 

800 

801class FocalMechanism(Object): 

802 ''' 

803 Description of the focal mechanism of an event. 

804 

805 It includes different descriptions like nodal planes, principal axes, and 

806 a moment tensor. The moment tensor description is provided by objects of 

807 the class :py:class:`MomentTensor` which can be specified as 

808 child elements of :py:class:`FocalMechanism`. 

809 ''' 

810 public_id = ResourceReference.T( 

811 xmlstyle='attribute', xmltagname='publicID') 

812 waveform_id_list = List.T(WaveformStreamID.T(xmltagname='waveformID')) 

813 comment_list = List.T(Comment.T()) 

814 moment_tensor_list = List.T(MomentTensor.T()) 

815 triggering_origin_id = ResourceReference.T( 

816 optional=True, xmltagname='triggeringOriginID') 

817 nodal_planes = NodalPlanes.T(optional=True) 

818 principal_axes = PrincipalAxes.T(optional=True) 

819 azimuthal_gap = Float.T(optional=True) 

820 station_polarity_count = Int.T(optional=True) 

821 misfit = Float.T(optional=True) 

822 station_distribution_ratio = Float.T(optional=True) 

823 method_id = ResourceReference.T(optional=True, xmltagname='methodID') 

824 evaluation_mode = EvaluationMode.T(optional=True) 

825 evaluation_status = EvaluationStatus.T(optional=True) 

826 creation_info = CreationInfo.T(optional=True) 

827 

828 

829class Origin(Object): 

830 ''' 

831 Representation of the focal time and geographical location of an earthquake 

832 hypocenter, as well as additional meta-information. 

833 

834 :py:class:`Origin` can have objects of type :py:class:`OriginUncertainty` 

835 and :py:class:`Arrival` as child elements. 

836 ''' 

837 public_id = ResourceReference.T( 

838 xmlstyle='attribute', xmltagname='publicID') 

839 composite_time_list = List.T(CompositeTime.T()) 

840 comment_list = List.T(Comment.T()) 

841 origin_uncertainty_list = List.T(OriginUncertainty.T()) 

842 arrival_list = List.T(Arrival.T()) 

843 time = TimeQuantity.T() 

844 longitude = RealQuantity.T() 

845 latitude = RealQuantity.T() 

846 depth = RealQuantity.T(optional=True) 

847 depth_type = OriginDepthType.T(optional=True) 

848 time_fixed = Bool.T(optional=True) 

849 epicenter_fixed = Bool.T(optional=True) 

850 reference_system_id = ResourceReference.T( 

851 optional=True, xmltagname='referenceSystemID') 

852 method_id = ResourceReference.T(optional=True, xmltagname='methodID') 

853 earth_model_id = ResourceReference.T( 

854 optional=True, xmltagname='earthModelID') 

855 quality = OriginQuality.T(optional=True) 

856 type = OriginType.T(optional=True) 

857 region = Region.T(optional=True) 

858 evaluation_mode = EvaluationMode.T(optional=True) 

859 evaluation_status = EvaluationStatus.T(optional=True) 

860 creation_info = CreationInfo.T(optional=True) 

861 

862 def position_values(self): 

863 lat = self.latitude.value 

864 lon = self.longitude.value 

865 if not self.depth: 

866 logger.warn( 

867 'Origin %s: Depth is undefined. Set to depth=0.' % 

868 self.public_id) 

869 depth = 0. 

870 else: 

871 depth = self.depth.value 

872 

873 return lat, lon, depth 

874 

875 def get_pyrocko_event(self): 

876 lat, lon, depth = self.position_values() 

877 otime = self.time.value 

878 if self.creation_info: 

879 cat = self.creation_info.agency_id 

880 else: 

881 cat = None 

882 

883 return event.Event( 

884 name=self.public_id, 

885 lat=lat, 

886 lon=lon, 

887 time=otime, 

888 depth=depth, 

889 catalog=cat, 

890 region=self.region) 

891 

892 

893class Event(Object): 

894 ''' 

895 Representation of a seismic event. 

896 

897 The Event does not necessarily need to be a tectonic earthquake. An event 

898 is usually associated with one or more origins, which contain information 

899 about focal time and geographical location of the event. Multiple origins 

900 can cover automatic and manual locations, a set of location from different 

901 agencies, locations generated with different location programs and earth 

902 models, etc. Furthermore, an event is usually associated with one or more 

903 magnitudes, and with one or more focal mechanism determinations. In 

904 standard QuakeML-BED, :py:class:`Origin`, :py:class:`Magnitude`, 

905 :py:class:`StationMagnitude`, and :py:class:`FocalMechanism` are child 

906 elements of Event. In BED-RT all these elements are on the same hierarchy 

907 level as child elements of :py:class:`EventParameters`. The association of 

908 origins, magnitudes, and focal mechanisms to a particular event is 

909 expressed using references inside :py:class:`Event`. 

910 ''' 

911 public_id = ResourceReference.T( 

912 xmlstyle='attribute', xmltagname='publicID') 

913 description_list = List.T(EventDescription.T()) 

914 comment_list = List.T(Comment.T()) 

915 focal_mechanism_list = List.T(FocalMechanism.T()) 

916 amplitude_list = List.T(Amplitude.T()) 

917 magnitude_list = List.T(Magnitude.T()) 

918 station_magnitude_list = List.T(StationMagnitude.T()) 

919 origin_list = List.T(Origin.T()) 

920 pick_list = List.T(Pick.T()) 

921 preferred_origin_id = ResourceReference.T( 

922 optional=True, xmltagname='preferredOriginID') 

923 preferred_magnitude_id = ResourceReference.T( 

924 optional=True, xmltagname='preferredMagnitudeID') 

925 preferred_focal_mechanism_id = ResourceReference.T( 

926 optional=True, xmltagname='preferredFocalMechanismID') 

927 type = EventType.T( 

928 optional=True) 

929 type_certainty = EventTypeCertainty.T( 

930 optional=True) 

931 creation_info = CreationInfo.T( 

932 optional=True) 

933 region = Region.T( 

934 optional=True) 

935 

936 def describe(self): 

937 return '''%s: 

938 origins: %i %s 

939 magnitudes: %i %s 

940 focal_machanisms: %i %s 

941 picks: %i 

942 amplitudes: %i 

943 station_magnitudes: %i''' % ( 

944 self.public_id, 

945 len(self.origin_list), 

946 '@' if self.preferred_origin_id else '-', 

947 len(self.magnitude_list), 

948 '@' if self.preferred_magnitude_id else '-', 

949 len(self.focal_mechanism_list), 

950 '@' if self.preferred_focal_mechanism_id else '-', 

951 len(self.pick_list), 

952 len(self.amplitude_list), 

953 len(self.station_magnitude_list)) 

954 

955 def get_pyrocko_phase_markers(self): 

956 event = self.get_pyrocko_event() 

957 return [ 

958 p.get_pyrocko_phase_marker(event=event) for p in self.pick_list] 

959 

960 def get_pyrocko_event(self): 

961 ''' 

962 Convert into Pyrocko event object. 

963 

964 Uses *preferred* origin, magnitude, and moment tensor. If no preferred 

965 item is specified, it picks the first from the list and emits a 

966 warning. 

967 ''' 

968 

969 origin = self.preferred_origin 

970 if not origin and self.origin_list: 

971 origin = self.origin_list[0] 

972 if len(self.origin_list) > 1: 

973 logger.warn( 

974 'Event %s: No preferred origin set, ' 

975 'more than one available, using first' % self.public_id) 

976 

977 if not origin: 

978 raise QuakeMLError( 

979 'No origin available for event: %s' % self.public_id) 

980 

981 ev = origin.get_pyrocko_event() 

982 

983 foc_mech = self.preferred_focal_mechanism 

984 if not foc_mech and self.focal_mechanism_list: 

985 foc_mech = self.focal_mechanism_list[0] 

986 if len(self.focal_mechanism_list) > 1: 

987 logger.warn( 

988 'Event %s: No preferred focal mechanism set, ' 

989 'more than one available, using first' % ev.name) 

990 

991 if foc_mech and foc_mech.moment_tensor_list: 

992 ev.moment_tensor = \ 

993 foc_mech.moment_tensor_list[0].pyrocko_moment_tensor() 

994 

995 if len(foc_mech.moment_tensor_list) > 1: 

996 logger.warn( 

997 'more than one moment tensor available, using first') 

998 

999 mag = None 

1000 pref_mag = self.preferred_magnitude 

1001 if pref_mag: 

1002 mag = pref_mag 

1003 elif self.magnitude_list: 

1004 mag = self.magnitude_list[0] 

1005 if len(self.magnitude_list) > 1: 

1006 logger.warn( 

1007 'Event %s: No preferred magnitude set, ' 

1008 'more than one available, using first' % ev.name) 

1009 

1010 if mag: 

1011 ev.magnitude = mag.mag.value 

1012 ev.magnitude_type = mag.type 

1013 

1014 ev.region = self.get_effective_region() 

1015 

1016 return ev 

1017 

1018 def get_effective_region(self): 

1019 if self.region: 

1020 return self.region 

1021 

1022 for desc in self.description_list: 

1023 if desc.type in ('Flinn-Engdahl region', 'region name'): 

1024 return desc.text 

1025 

1026 return None 

1027 

1028 @property 

1029 def preferred_origin(self): 

1030 return one_element_or_none( 

1031 [x for x in self.origin_list 

1032 if x.public_id == self.preferred_origin_id]) 

1033 

1034 @property 

1035 def preferred_magnitude(self): 

1036 return one_element_or_none( 

1037 [x for x in self.magnitude_list 

1038 if x.public_id == self.preferred_magnitude_id]) 

1039 

1040 @property 

1041 def preferred_focal_mechanism(self): 

1042 return one_element_or_none( 

1043 [x for x in self.focal_mechanism_list 

1044 if x.public_id == self.preferred_focal_mechanism_id]) 

1045 

1046 

1047class EventParameters(Object): 

1048 ''' 

1049 In the bulletin-type (non real-time) model, this class serves as a 

1050 container for Event objects. 

1051 

1052 In the real-time version, it can hold objects of type :py:class:`Event`, 

1053 :py:class:`Origin`, :py:class:`Magnitude`, :py:class:`StationMagnitude`, 

1054 :py:class:`FocalMechanism`, Reading, :py:class:`Amplitude` and 

1055 :py:class:`Pick` (real-time mode is not covered by this module at the 

1056 moment). 

1057 ''' 

1058 public_id = ResourceReference.T( 

1059 xmlstyle='attribute', xmltagname='publicID') 

1060 comment_list = List.T(Comment.T()) 

1061 event_list = List.T(Event.T(xmltagname='event')) 

1062 description = Unicode.T(optional=True) 

1063 creation_info = CreationInfo.T(optional=True) 

1064 

1065 

1066class QuakeML(Object): 

1067 ''' 

1068 QuakeML data container. 

1069 ''' 

1070 xmltagname = 'quakeml' 

1071 xmlns = 'http://quakeml.org/xmlns/quakeml/1.2' 

1072 guessable_xmlns = [xmlns, guts_xmlns] 

1073 

1074 event_parameters = EventParameters.T(optional=True) 

1075 

1076 def get_events(self): 

1077 return self.event_parameters.event_list 

1078 

1079 def get_pyrocko_events(self): 

1080 ''' 

1081 Get event information in Pyrocko's basic event format. 

1082 

1083 :rtype: 

1084 List of :py:class:`pyrocko.model.event.Event` objects. 

1085 ''' 

1086 events = [] 

1087 for e in self.event_parameters.event_list: 

1088 events.append(e.get_pyrocko_event()) 

1089 

1090 return events 

1091 

1092 def get_pyrocko_phase_markers(self): 

1093 ''' 

1094 Get pick information in Pyrocko's basic marker format. 

1095 

1096 :rtype: 

1097 List of :py:class:`pyrocko.gui.marker.PhaseMarker` objects. 

1098 ''' 

1099 markers = [] 

1100 for e in self.event_parameters.event_list: 

1101 markers.extend(e.get_pyrocko_phase_markers()) 

1102 

1103 return markers 

1104 

1105 @classmethod 

1106 def load_xml(cls, stream=None, filename=None, string=None): 

1107 ''' 

1108 Load QuakeML data from stream, file or string. 

1109 

1110 :param stream: 

1111 Stream open for reading in binary mode. 

1112 :type stream: 

1113 file-like object, optional 

1114 

1115 :param filename: 

1116 Path to file to be opened for reading. 

1117 :type filename: 

1118 str, optional 

1119 

1120 :param string: 

1121 String with QuakeML data to be deserialized. 

1122 :type string: 

1123 str, optional 

1124 

1125 The arguments ``stream``, ``filename``, and ``string`` are mutually 

1126 exclusive. 

1127 

1128 :returns: 

1129 Parsed QuakeML data structure. 

1130 :rtype: 

1131 :py:class:`QuakeML` object 

1132 

1133 ''' 

1134 

1135 return super(QuakeML, cls).load_xml( 

1136 stream=stream, 

1137 filename=filename, 

1138 string=string, 

1139 ns_hints=[ 

1140 'http://quakeml.org/xmlns/quakeml/1.2', 

1141 'http://quakeml.org/xmlns/bed/1.2'], 

1142 ns_ignore=True)