Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

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

1001

1002

1003

1004

1005

1006

1007

1008

1009

1010

1011

1012

1013

1014

1015

1016

1017

1018

1019

1020

1021

1022

1023

1024

1025

1026

1027

1028

1029

1030

1031

1032

1033

1034

1035

1036

1037

1038

1039

1040

1041

1042

1043

1044

1045

1046

1047

1048

1049

1050

1051

1052

1053

1054

1055

1056

1057

1058

1059

1060

1061

1062

1063

1064

1065

1066

1067

1068

1069

1070

1071

1072

1073

1074

1075

1076

1077

1078

1079

1080

1081

1082

1083

1084

1085

1086

1087

1088

1089

1090

1091

1092

1093

1094

1095

1096

1097

1098

1099

1100

1101

1102

1103

1104

1105

1106

1107

1108

1109

1110

1111

1112

1113

1114

1115

1116

1117

1118

1119

1120

1121

1122

1123

1124

1125

1126

1127

1128

1129

1130

1131

1132

1133

1134

1135

1136

1137

1138

1139

1140

1141

1142

1143

1144

1145

1146

1147

1148

1149

1150

1151

1152

1153

1154

1155

1156

1157

1158

1159

1160

1161

1162

1163

1164

1165

1166

1167

1168

1169

1170

1171

1172

1173

1174

1175

1176

1177

1178

1179

1180

1181

1182

1183

1184

1185

1186

1187

1188

1189

1190

1191

1192

1193

1194

1195

1196

1197

1198

1199

1200

1201

1202

1203

1204

1205

1206

1207

1208

1209

1210

1211

1212

1213

1214

1215

1216

1217

1218

1219

1220

1221

1222

1223

1224

1225

1226

1227

1228

1229

1230

1231

1232

1233

1234

1235

1236

1237

1238

1239

1240

1241

1242

1243

1244

1245

1246

1247

1248

1249

1250

1251

1252

1253

1254

1255

1256

1257

1258

1259

1260

1261

1262

1263

1264

1265

1266

1267

1268

1269

1270

1271

1272

1273

1274

1275

1276

1277

1278

1279

1280

1281

1282

1283

1284

1285

1286

1287

1288

1289

1290

1291

1292

1293

1294

1295

1296

1297

1298

1299

1300

1301

1302

1303

1304

1305

1306

1307

1308

1309

1310

1311

1312

1313

1314

1315

1316

1317

1318

1319

1320

1321

1322

1323

1324

1325

1326

1327

1328

1329

1330

1331

1332

1333

1334

1335

1336

1337

1338

1339

1340

1341

1342

1343

1344

1345

1346

1347

1348

1349

1350

1351

1352

1353

1354

1355

1356

1357

1358

1359

1360

1361

1362

1363

1364

1365

1366

1367

1368

1369

1370

1371

1372

1373

1374

1375

1376

1377

1378

1379

1380

1381

1382

1383

1384

1385

1386

1387

1388

1389

1390

1391

1392

1393

1394

1395

1396

1397

1398

1399

1400

1401

1402

1403

1404

1405

1406

1407

1408

1409

1410

1411

1412

1413

1414

1415

1416

1417

1418

1419

1420

1421

1422

1423

1424

1425

1426

1427

1428

1429

1430

1431

1432

1433

1434

1435

1436

1437

1438

1439

1440

1441

1442

1443

1444

1445

1446

1447

1448

1449

1450

1451

1452

1453

1454

1455

1456

1457

1458

1459

1460

1461

1462

1463

1464

1465

1466

1467

1468

1469

1470

1471

1472

1473

1474

1475

1476

1477

1478

1479

1480

1481

1482

1483

1484

1485

1486

1487

1488

1489

1490

1491

1492

1493

1494

1495

1496

1497

1498

1499

1500

1501

1502

1503

1504

1505

1506

1507

1508

1509

1510

1511

1512

1513

1514

1515

1516

1517

1518

1519

1520

1521

1522

1523

1524

1525

1526

1527

1528

1529

1530

1531

1532

1533

1534

1535

1536

1537

1538

1539

1540

1541

1542

1543

1544

1545

1546

1547

1548

1549

1550

1551

1552

1553

1554

1555

1556

1557

1558

1559

1560

1561

1562

1563

1564

1565

1566

#!/usr/bin/env python 

#------------------------------------------------------------------- 

# Filename: sacio.py 

#  Purpose: Read & Write Seismograms, Format SAC. 

#   Author: Yannik Behr, C. J. Ammon's 

#    Email: yannik.behr@vuw.ac.nz 

# 

# Copyright (C) 2008-2012 Yannik Behr, C. J. Ammon's 

#------------------------------------------------------------------- 

from obspy import UTCDateTime, Trace 

from obspy.core.util import gps2DistAzimuth, loadtxt, AttribDict 

import numpy as np 

import os 

import string 

import time 

import warnings 

""" 

Low-level module internally used for handling SAC files 

 

An object-oriented version of C. J. Ammon's SAC I/O module. 

 

:copyright: 

    The ObsPy Development Team (devs@obspy.org) & C. J. Ammon 

:license: 

    GNU Lesser General Public License, Version 3 

    (http://www.gnu.org/copyleft/lesser.html) 

""" 

 

 

# we put here everything but the time, they are going to stats.starttime 

# left SAC attributes, right trace attributes, see also 

# http://www.iris.edu/KB/questions/13/SAC+file+format 

convert_dict = {'npts': 'npts', 

                'delta': 'delta', 

                'kcmpnm': 'channel', 

                'kstnm': 'station', 

                'scale': 'calib', 

                'knetwk': 'network', 

                'khole': 'location'} 

 

# all the sac specific extras, the SAC reference time specific headers are 

# handled separately and are directly controlled by trace.stats.starttime. 

SAC_EXTRA = ('depmin', 'depmax', 'odelta', 'o', 'a', 't0', 't1', 't2', 't3', 

             't4', 't5', 't6', 't7', 't8', 't9', 'f', 'stla', 'stlo', 'stel', 

             'stdp', 'evla', 'evlo', 'evdp', 'mag', 'user0', 'user1', 'user2', 

             'user3', 'user4', 'user5', 'user6', 'user7', 'user8', 'user9', 

             'dist', 'az', 'baz', 'gcarc', 'depmen', 'cmpaz', 'cmpinc', 

             'nvhdr', 'norid', 'nevid', 'nwfid', 'iftype', 'idep', 'iztype', 

             'iinst', 'istreg', 'ievreg', 'ievtype', 'iqual', 'isynth', 

             'imagtyp', 'imagsrc', 'leven', 'lpspol', 'lovrok', 'lcalda', 

             'kevnm', 'ko', 'ka', 'kt0', 'kt1', 'kt2', 'kt3', 'kt4', 'kt5', 

             'kt6', 'kt7', 'kt8', 'kt9', 'kf', 'kuser0', 'kuser1', 'kuser2', 

             'kdatrd', 'kinst', 'cmpinc', 'xminimum', 'xmaximum', 'yminimum', 

             'ymaximum', 'unused6', 'unused7', 'unused8', 'unused9', 

             'unused10', 'unused11', 'unused12') 

 

FDICT = {'delta': 0, 'depmin': 1, 'depmax': 2, 'scale': 3, 

         'odelta': 4, 'b': 5, 'e': 6, 'o': 7, 'a': 8, 'int1': 9, 

         't0': 10, 't1': 11, 't2': 12, 't3': 13, 't4': 14, 

         't5': 15, 't6': 16, 't7': 17, 't8': 18, 't9': 19, 

         'f': 20, 'stla': 31, 'stlo': 32, 'stel': 33, 'stdp': 34, 

         'evla': 35, 'evlo': 36, 'evdp': 38, 'mag': 39, 

         'user0': 40, 'user1': 41, 'user2': 42, 'user3': 43, 

         'user4': 44, 'user5': 45, 'user6': 46, 'user7': 47, 

         'user8': 48, 'user9': 49, 'dist': 50, 'az': 51, 

         'baz': 52, 'gcarc': 53, 'depmen': 56, 'cmpaz': 57, 

         'cmpinc': 58, 'xminimum': 59, 'xmaximum': 60, 

         'yminimum': 61, 'ymaximum': 62, 'unused6': 63, 

         'unused7': 64, 'unused8': 65, 'unused9': 66, 

         'unused10': 67, 'unused11': 68, 'unused12': 69} 

 

IDICT = {'nzyear': 0, 'nzjday': 1, 'nzhour': 2, 'nzmin': 3, 

         'nzsec': 4, 'nzmsec': 5, 'nvhdr': 6, 'norid': 7, 

         'nevid': 8, 'npts': 9, 'nwfid': 11, 

         'iftype': 15, 'idep': 16, 'iztype': 17, 'iinst': 19, 

         'istreg': 20, 'ievreg': 21, 'ievtype': 22, 'iqual': 23, 

         'isynth': 24, 'imagtyp': 25, 'imagsrc': 26, 

         'leven': 35, 'lpspol': 36, 'lovrok': 37, 

         'lcalda': 38} 

 

SDICT = {'kstnm': 0, 'kevnm': 1, 'khole': 2, 'ko': 3, 'ka': 4, 

         'kt0': 5, 'kt1': 6, 'kt2': 7, 'kt3': 8, 'kt4': 9, 

         'kt5': 10, 'kt6': 11, 'kt7': 12, 'kt8': 13, 

         'kt9': 14, 'kf': 15, 'kuser0': 16, 'kuser1': 17, 

         'kuser2': 18, 'kcmpnm': 19, 'knetwk': 20, 

         'kdatrd': 21, 'kinst': 22} 

 

 

class SacError(Exception): 

    """ 

    Raised if the SAC file is corrupt or if necessary information 

    in the SAC file is missing. 

    """ 

    pass 

 

 

class SacIOError(Exception): 

    """ 

    Raised if the given SAC file can't be read. 

    """ 

    pass 

 

 

def _isText(filename, blocksize=512): 

    """ 

    Check if it is a text or a binary file. 

    """ 

    # This should  always be true if a file is a text-file and only true for a 

    # binary file in rare occasions (see Recipe 173220 found on 

    # http://code.activestate.com/) 

    text_characters = "".join(map(chr, range(32, 127)) + list("\n\r\t\b")) 

    _null_trans = string.maketrans("", "") 

    s = open(filename).read(blocksize) 

    if "\0" in s: 

        return False 

 

    if not s:  # Empty files are considered text 

        return True 

 

    # Get the non-text characters (maps a character to itself then 

    # use the 'remove' option to get rid of the text characters.) 

    t = s.translate(_null_trans, text_characters) 

 

    # If more than 30% non-text characters, then 

    # this is considered a binary file 

    if len(t) / len(s) > 0.30: 

        return 0 

    return True 

 

 

class SacIO(object): 

    """ 

    Class for SAC file IO. 

 

    Functions are given below, attributes/header 

    fields (described below) can be directly accessed (via the 

    :meth:`~obspy.sac.sacio.SacIO.__getattr__` method, see the link for 

    an example). 

 

    .. rubric::Description of attributes/header fields (based on SacIris_). 

 

    .. _SacIris: http://www.iris.edu/manuals/sac/SAC_Manuals/FileFormatPt2.html 

 

    ============ ==== ========================================================= 

    Field Name   Type Description 

    ============ ==== ========================================================= 

    npts         N    Number of points per data component. [required] 

    nvhdr        N    Header version number. Current value is the integer 6. 

                      Older version data (NVHDR < 6) are automatically updated 

                      when read into sac. [required] 

    b            F    Beginning value of the independent variable. [required] 

    e            F    Ending value of the independent variable. [required] 

    iftype       I    Type of file [required]: 

                          * ITIME {Time series file} 

                          * IRLIM {Spectral file---real and imaginary} 

                          * IAMPH {Spectral file---amplitude and phase} 

                          * IXY {General x versus y data} 

                          * IXYZ {General XYZ (3-D) file} 

    leven        L    TRUE if data is evenly spaced. [required] 

    delta        F    Increment between evenly spaced samples (nominal value). 

                      [required] 

    odelta       F    Observed increment if different from nominal value. 

    idep         I    Type of dependent variable: 

                          * IUNKN (Unknown) 

                          * IDISP (Displacement in nm) 

                          * IVEL (Velocity in nm/sec) 

                          * IVOLTS (Velocity in volts) 

                          * IACC (Acceleration in nm/sec/sec) 

    scale        F    Multiplying scale factor for dependent variable 

                      [not currently used] 

    depmin       F    Minimum value of dependent variable. 

    depmax       F    Maximum value of dependent variable. 

    depmen       F    Mean value of dependent variable. 

    nzyear       N    GMT year corresponding to reference (zero) time in file. 

    nyjday       N    GMT julian day. 

    nyhour       N    GMT hour. 

    nzmin        N    GMT minute. 

    nzsec        N    GMT second. 

    nzmsec       N    GMT millisecond. 

    iztype       I    Reference time equivalence: 

                          * IUNKN (5): Unknown 

                          * IB (9): Begin time 

                          * IDAY (10): Midnight of refernece GMT day 

                          * IO (11): Event origin time 

                          * IA (12): First arrival time 

                          * ITn (13-22): User defined time pick n, n=0,9 

    o            F    Event origin time (seconds relative to reference time.) 

    a            F    First arrival time (seconds relative to reference time.) 

    ka           K    First arrival time identification. 

    f            F    Fini or end of event time (seconds relative to reference 

                      time.) 

    tn           F    User defined time {ai n}=0,9 (seconds picks or markers 

                      relative to reference time). 

    kt{ai n}     K    A User defined time {ai n}=0,9.  pick identifications 

    kinst        K    Generic name of recording instrument 

    iinst        I    Type of recording instrument. [currently not used] 

    knetwk       K    Name of seismic network. 

    kstnm        K    Station name. 

    istreg       I    Station geographic region. [not currently used] 

    stla         F    Station latitude (degrees, north positive) 

    stlo         F    Station longitude (degrees, east positive). 

    stel         F    Station elevation (meters). [not currently used] 

    stdp         F    Station depth below surface (meters). [not currently 

                      used] 

    cmpaz        F    Component azimuth (degrees, clockwise from north). 

    cmpinc       F    Component incident angle (degrees, from vertical). 

    kcmpnm       K    Component name. 

    lpspol       L    TRUE if station components have a positive polarity 

                      (left-hand rule). 

    kevnm        K    Event name. 

    ievreg       I    Event geographic region. [not currently used] 

    evla         F    Event latitude (degrees north positive). 

    evlo         F    Event longitude (degrees east positive). 

    evel         F    Event elevation (meters). [not currently used] 

    evdp         F    Event depth below surface (meters). [not currently used] 

    mag          F    Event magnitude. 

    imagtyp      I    Magnitude type: 

                          * IMB (Bodywave Magnitude) 

                          * IMS (Surfacewave Magnitude) 

                          * IML (Local Magnitude) 

                          * IMW (Moment Magnitude) 

                          * IMD (Duration Magnitude) 

                          * IMX (User Defined Magnitude) 

    imagsrc      I    Source of magnitude information: 

                          * INEIC (National Earthquake Information Center) 

                          * IPDE (Preliminary Determination of Epicenter) 

                          * IISC (Internation Seismological Centre) 

                          * IREB (Reviewed Event Bulletin) 

                          * IUSGS (US Geological Survey) 

                          * IBRK (UC Berkeley) 

                          * ICALTECH (California Institute of Technology) 

                          * ILLNL (Lawrence Livermore National Laboratory) 

                          * IEVLOC (Event Location (computer program) ) 

                          * IJSOP (Joint Seismic Observation Program) 

                          * IUSER (The individual using SAC2000) 

                          * IUNKNOWN (unknown) 

    ievtyp       I    Type of event: 

                          * IUNKN (Unknown) 

                          * INUCL (Nuclear event) 

                          * IPREN (Nuclear pre-shot event) 

                          * IPOSTN (Nuclear post-shot event) 

                          * IQUAKE (Earthquake) 

                          * IPREQ (Foreshock) 

                          * IPOSTQ (Aftershock) 

                          * ICHEM (Chemical explosion) 

                          * IQB (Quarry or mine blast confirmed by quarry) 

                          * IQB1 (Quarry/mine blast with designed shot 

                            info-ripple fired) 

                          * IQB2 (Quarry/mine blast with observed shot 

                            info-ripple fired) 

                          * IQMT (Quarry/mining-induced events: 

                            tremors and rockbursts) 

                          * IEQ (Earthquake) 

                          * IEQ1 (Earthquakes in a swarm or aftershock 

                            sequence) 

                          * IEQ2 (Felt earthquake) 

                          * IME (Marine explosion) 

                          * IEX (Other explosion) 

                          * INU (Nuclear explosion) 

                          * INC (Nuclear cavity collapse) 

                          * IO\_ (Other source of known origin) 

                          * IR (Regional event of unknown origin) 

                          * IT (Teleseismic event of unknown origin) 

                          * IU (Undetermined or conflicting information) 

                          * IOTHER (Other) 

    nevid        N    Event ID (CSS 3.0) 

    norid        N    Origin ID (CSS 3.0) 

    nwfid        N    Waveform ID (CSS 3.0) 

    khole        k    Hole identification if nuclear event. 

    dist         F    Station to event distance (km). 

    az           F    Event to station azimuth (degrees). 

    baz          F    Station to event azimuth (degrees). 

    gcarc        F    Station to event great circle arc length (degrees). 

    lcalda       L    TRUE if DIST AZ BAZ and GCARC are to be calculated from 

                      st event coordinates. 

    iqual        I    Quality of data [not currently used]: 

                          * IGOOD (Good data) 

                          * IGLCH (Glitches) 

                          * IDROP (Dropouts) 

                          * ILOWSN (Low signal to noise ratio) 

                          * IOTHER (Other) 

    isynth       I    Synthetic data flag [not currently used]: 

                          * IRLDTA (Real data) 

                          * ????? (Flags for various synthetic seismogram 

                            codes) 

    user{ai n}   F    User defined variable storage area {ai n}=0,9. 

    kuser{ai n}  K    User defined variable storage area {ai n}=0,2. 

    lovrok       L    TRUE if it is okay to overwrite this file on disk. 

    ============ ==== ========================================================= 

    """ 

 

    def __init__(self, filen=False, headonly=False, alpha=False, 

                 debug_headers=False): 

        self.byteorder = 'little' 

        self.InitArrays() 

        self.debug_headers = debug_headers 

        if filen is False: 

            return 

        # parse Trace object if we get one 

        if isinstance(filen, Trace): 

            self.readTrace(filen) 

            return 

        if alpha: 

            if headonly: 

                self.ReadSacXYHeader(filen) 

            else: 

                self.ReadSacXY(filen) 

        elif headonly: 

            self.ReadSacHeader(filen) 

        else: 

            self.ReadSacFile(filen) 

 

    def InitArrays(self): 

        """ 

        Function to initialize the floating, character and integer 

        header arrays (self.hf, self.hs, self.hi) with dummy values. This 

        function is useful for writing SAC files from artificial data, 

        thus the header arrays are not filled by a read method 

        beforehand 

 

        :return: Nothing 

        """ 

        # The SAC header has 70 floats, then 40 integers, then 192 bytes 

        # in strings. Store them in array (an convert the char to a 

        # list). That's a total of 632 bytes. 

        # 

        # allocate the array for header floats 

        self.hf = np.ndarray(70, dtype='<f4') 

        self.hf[:] = -12345.0 

        # 

        # allocate the array for header integers 

        self.hi = np.ndarray(40, dtype='<i4') 

        self.hi[:] = -12345 

        # 

        # allocate the array for header characters 

        self.hs = np.ndarray(24, dtype='|S8') 

        self.hs[:] = '-12345  '  # setting default value 

        # allocate the array for the points 

        self.seis = np.ndarray([], dtype='<f4') 

 

    def fromarray(self, trace, begin=0.0, delta=1.0, distkm=0, 

                  starttime=UTCDateTime("1970-01-01T00:00:00.000000")): 

        """ 

        Create a SAC file from an numpy.ndarray instance 

 

        >>> t = SacIO() 

        >>> b = np.arange(10) 

        >>> t.fromarray(b) 

        >>> t.GetHvalue('npts') 

        10 

        """ 

        if not isinstance(trace, np.ndarray): 

            raise SacError("input needs to be of instance numpy.ndarray") 

        else: 

            # Only copy the data if they are not of the required type 

            self.seis = np.require(trace, '<f4') 

        # convert stattime to sac reference time, if it is not default 

        if begin == -12345: 

            reftime = starttime 

        else: 

            reftime = starttime - begin 

        # if there are any micro-seconds, use begin to store them 

        # integer arithmetic 

        millisecond = reftime.microsecond // 1000 

        # integer arithmetic 

        microsecond = (reftime.microsecond - millisecond * 1000) 

        if microsecond != 0: 

            begin += microsecond * 1e-6 

        # set a few values that are required to create a valid SAC-file 

        self.SetHvalue('int1', 2) 

        self.SetHvalue('cmpaz', 0) 

        self.SetHvalue('cmpinc', 0) 

        self.SetHvalue('nvhdr', 6) 

        self.SetHvalue('leven', 1) 

        self.SetHvalue('lpspol', 1) 

        self.SetHvalue('lcalda', 0) 

        self.SetHvalue('nzyear', reftime.year) 

        self.SetHvalue('nzjday', reftime.strftime("%j")) 

        self.SetHvalue('nzhour', reftime.hour) 

        self.SetHvalue('nzmin', reftime.minute) 

        self.SetHvalue('nzsec', reftime.second) 

        self.SetHvalue('nzmsec', millisecond) 

        self.SetHvalue('kcmpnm', 'Z') 

        self.SetHvalue('evla', 0) 

        self.SetHvalue('evlo', 0) 

        self.SetHvalue('iftype', 1) 

        self.SetHvalue('npts', len(trace)) 

        self.SetHvalue('delta', delta) 

        self.SetHvalue('b', begin) 

        self.SetHvalue('e', begin + (len(trace) - 1) * delta) 

        self.SetHvalue('iztype', 9) 

        self.SetHvalue('dist', distkm) 

 

    def GetHvalue(self, item): 

        """ 

        Read SAC-header variable. 

 

        :param item: header variable name (e.g. 'npts' or 'delta') 

 

        >>> from obspy.sac import SacIO # doctest: +SKIP 

        >>> tr = SacIO('test.sac') # doctest: +SKIP 

        >>> tr.GetHvalue('npts') # doctest: +SKIP 

        100 

 

        This is equivalent to: 

 

        >>> SacIO().GetHvalueFromFile('test.sac','npts') # doctest: +SKIP 

        100 

 

        Or: 

 

        >>> tr = SacIO('test.sac') # doctest: +SKIP 

        >>> tr.npts # doctest: +SKIP 

        100 

 

        """ 

        key = item.lower()  # convert the item to lower case 

        if key in FDICT: 

            index = FDICT[key] 

            return(self.hf[index]) 

        elif key in IDICT: 

            index = IDICT[key] 

            return(self.hi[index]) 

        elif key in SDICT: 

            index = SDICT[key] 

            if index == 0: 

                myarray = self.hs[0] 

            elif index == 1: 

                myarray = self.hs[1] + self.hs[2] 

            else: 

                myarray = self.hs[index + 1]  # extra 1 is from item #2 

            return myarray 

        else: 

            raise SacError("Cannot find header entry for: " + item) 

 

    def SetHvalue(self, item, value): 

        """ 

        Assign new value to SAC-header variable. 

 

        :param item: SAC-header variable name 

        :param value: numeric or string value to be assigned to header 

                      variable. 

 

        >>> from obspy.sac import SacIO 

        >>> tr = SacIO() 

        >>> tr.GetHvalue('kstnm') 

        '-12345  ' 

        >>> tr.SetHvalue('kstnm', 'STA_NEW') 

        >>> tr.GetHvalue('kstnm') 

        'STA_NEW ' 

 

 

        """ 

        key = item.lower()  # convert the item to lower case 

        # 

        if key in FDICT: 

            index = FDICT[key] 

            self.hf[index] = float(value) 

        elif key in IDICT: 

            index = IDICT[key] 

            self.hi[index] = int(value) 

        elif key in SDICT: 

            index = SDICT[key] 

            if value: 

                value = '%-8s' % value 

            else: 

                value = '-12345  ' 

            if index == 0: 

                self.hs[0] = value 

            elif index == 1: 

                value1 = '%-8s' % value[0:8] 

                value2 = '%-8s' % value[8:16] 

                self.hs[1] = value1 

                self.hs[2] = value2 

            else: 

                self.hs[index + 1] = value 

        else: 

            raise SacError("Cannot find header entry for: " + item) 

 

    def IsSACfile(self, name, fsize=True, lenchk=False): 

        """ 

        Test for a valid SAC file using arrays. 

 

        :param f: filename (Sac binary). 

        """ 

        try: 

            npts = self.GetHvalue('npts') 

        except: 

            raise SacError("Unable to read number of points from header") 

        if lenchk and npts != len(self.seis): 

            raise SacError("Number of points in header and " + \ 

                           "length of trace inconsistent!") 

        if fsize: 

            st = os.stat(name)  # file's size = st[6] 

            sizecheck = st[6] - (632 + 4 * int(npts)) 

            # size check info 

            if sizecheck != 0: 

                msg = "File-size and theoretical size are inconsistent: %s\n" \ 

                      "Check that headers are consistent with time series." 

                raise SacError(msg % name) 

        # get the SAC file version number 

        version = self.GetHvalue('nvhdr') 

        if version < 0 or version > 20: 

            raise SacError("Unknown header version!") 

        if self.GetHvalue('delta') <= 0: 

            raise SacError("Delta < 0 is not a valid header entry!") 

 

    def ReadSacHeader(self, fname): 

        """ 

        Reads only the header portion of a binary SAC-file. 

 

        :param f: filename (SAC binary). 

 

        >>> from obspy.sac import SacIO # doctest: +SKIP 

        >>> tr = SacIO() # doctest: +SKIP 

        >>> tr.ReadSacHeader('test.sac') # doctest: +SKIP 

 

        This is equivalent to: 

 

        >>> tr = SacIO('test.sac', headonly=True)  # doctest: +SKIP 

        """ 

        #### check if file exists 

        try: 

            #### open the file 

            f = open(fname, 'rb') 

        except IOError: 

            raise SacIOError("No such file: " + fname) 

        #-------------------------------------------------------------- 

        # parse the header 

        # 

        # The sac header has 70 floats, 40 integers, then 192 bytes 

        #    in strings. Store them in array (an convert the char to a 

        #    list). That's a total of 632 bytes. 

        #-------------------------------------------------------------- 

        self.hf = np.fromfile(f, dtype='<f4', count=70) 

        self.hi = np.fromfile(f, dtype='<i4', count=40) 

        # read in the char values 

        self.hs = np.fromfile(f, dtype='|S8', count=24) 

        if len(self.hf) != 70 or len(self.hi) != 40 or len(self.hs) != 24: 

            self.hf = self.hi = self.hs = None 

            f.close() 

            raise SacIOError("Cannot read all header values") 

        try: 

            self.IsSACfile(fname) 

        except SacError, e: 

            try: 

                # if it is not a valid SAC-file try with big endian 

                # byte order 

                f.seek(0, 0) 

                self.hf = np.fromfile(f, dtype='>f4', count=70) 

                self.hi = np.fromfile(f, dtype='>i4', count=40) 

                # read in the char values 

                self.hs = np.fromfile(f, dtype='|S8', count=24) 

                self.IsSACfile(fname) 

                self.byteorder = 'big' 

            except SacError, e: 

                self.hf = self.hi = self.hs = None 

                f.close() 

                raise SacError(e) 

        try: 

            self._get_date() 

        except SacError: 

            warnings.warn('Cannot determine date') 

        if self.GetHvalue('lcalda'): 

            try: 

                self._get_dist() 

            except SacError: 

                pass 

        f.close() 

 

    def WriteSacHeader(self, fname): 

        """ 

        Writes an updated header to an 

        existing binary SAC-file. 

 

        :param f: filename (SAC binary). 

 

        >>> from obspy.sac import SacIO # doctest: +SKIP 

        >>> tr = SacIO('test.sac') # doctest: +SKIP 

        >>> tr.WriteSacBinary('test2.sac') # doctest: +SKIP 

        >>> u = SacIO('test2.sac') # doctest: +SKIP 

        >>> u.SetHvalue('kevnm','hullahulla') # doctest: +SKIP 

        >>> u.WriteSacHeader('test2.sac') # doctest: +SKIP 

        >>> u.GetHvalueFromFile('test2.sac',"kevnm") # doctest: +SKIP 

        'hullahulla      ' 

        """ 

        #-------------------------------------------------------------- 

        # open the file 

        # 

        try: 

            os.path.exists(fname) 

        except IOError: 

            warnings.warn("No such file: " + fname, category=Warning) 

        else: 

            f = open(fname, 'rb+')  # open file for modification 

            f.seek(0, 0)  # set pointer to the file beginning 

            try: 

                # write the header 

                self.hf.tofile(f) 

                self.hi.tofile(f) 

                self.hs.tofile(f) 

            except Exception, e: 

                f.close() 

                raise SacError("Cannot write header to file: " + fname, e) 

        f.close() 

 

    def ReadSacFile(self, fname): 

        """ 

        Read read in the header and data in a SAC file 

 

        The header is split into three arrays - floats, ints, and strings and 

        the data points are returned in the array seis 

 

        :param f: filename (SAC binary) 

 

        >>> from obspy.sac import SacIO # doctest: +SKIP 

        >>> tr = SacIO() # doctest: +SKIP 

        >>> tr.ReadSacFile('test.sac') # doctest: +SKIP 

 

        This is equivalent to: 

 

        >>> tr = SacIO('test.sac')  # doctest: +SKIP 

        """ 

        try: 

            #### open the file 

            f = open(fname, 'rb') 

        except IOError: 

            raise SacIOError("No such file: " + fname) 

        #-------------------------------------------------------------- 

        # parse the header 

        # 

        # The sac header has 70 floats, 40 integers, then 192 bytes 

        #    in strings. Store them in array (an convert the char to a 

        #    list). That's a total of 632 bytes. 

        #-------------------------------------------------------------- 

        self.hf = np.fromfile(f, dtype='<f4', count=70) 

        self.hi = np.fromfile(f, dtype='<i4', count=40) 

        # read in the char values 

        self.hs = np.fromfile(f, dtype='|S8', count=24) 

        if len(self.hf) != 70 or len(self.hi) != 40 or len(self.hs) != 24: 

            self.hf = self.hi = self.hs = None 

            f.close() 

            raise SacIOError("Cannot read all header values") 

        ##### only continue if it is a SAC file 

        try: 

            self.IsSACfile(fname) 

        except SacError: 

            try: 

                # if it is not a valid SAC-file try with big endian 

                # byte order 

                f.seek(0, 0) 

                self.hf = np.fromfile(f, dtype='>f4', count=70) 

                self.hi = np.fromfile(f, dtype='>i4', count=40) 

                # read in the char values 

                self.hs = np.fromfile(f, dtype='|S8', count=24) 

                self.IsSACfile(fname) 

                self.byteorder = 'big' 

            except SacError, e: 

                f.close() 

                raise SacError(e) 

        #-------------------------------------------------------------- 

        # read in the seismogram points 

        #-------------------------------------------------------------- 

        # you just have to know it's in the 10th place 

        # actually, it's in the SAC manual 

        npts = self.hi[9] 

        if self.byteorder == 'big': 

            self.seis = np.fromfile(f, dtype='>f4', count=npts) 

        else: 

            self.seis = np.fromfile(f, dtype='<f4', count=npts) 

        if len(self.seis) != npts: 

            self.hf = self.hi = self.hs = self.seis = None 

            f.close() 

            raise SacIOError("Cannot read any or only some data points") 

        try: 

            self._get_date() 

        except SacError: 

            warnings.warn('Cannot determine date') 

        if self.GetHvalue('lcalda'): 

            try: 

                self._get_dist() 

            except SacError: 

                pass 

        f.close() 

 

    def ReadSacXY(self, fname): 

        """ 

        Read SAC XY files (ascii) 

 

        :param f: filename (SAC ascii). 

 

        >>> from obspy.sac import SacIO # doctest: +SKIP 

        >>> tr = SacIO() # doctest: +SKIP 

        >>> tr.ReadSacXY('testxy.sac') # doctest: +SKIP 

        >>> tr.GetHvalue('npts') # doctest: +SKIP 

        100 

 

        This is equivalent to: 

 

        >>> tr = SacIO('testxy.sac',alpha=True)  # doctest: +SKIP 

 

        Reading only the header portion of alphanumeric SAC-files is currently 

        not supported. 

        """ 

        ###### open the file 

        try: 

            f = open(fname, 'r') 

        except IOError: 

            raise SacIOError("No such file: " + fname) 

        try: 

            #-------------------------------------------------------------- 

            # parse the header 

            # 

            # The sac header has 70 floats, 40 integers, then 192 bytes 

            #    in strings. Store them in array (an convert the char to a 

            #    list). That's a total of 632 bytes. 

            #-------------------------------------------------------------- 

            # read in the float values 

            self.hf = np.fromfile(f, dtype='<f4', count=70, sep=" ") 

            # read in the int values 

            self.hi = np.fromfile(f, dtype='<i4', count=40, sep=" ") 

            # reading in the string part is a bit more complicated 

            # because every string field has to be 8 characters long 

            # apart from the second field which is 16 characters long 

            # resulting in a total length of 192 characters 

            for i in xrange(0, 24, 3): 

                self.hs[i:i + 3] = np.fromfile(f, dtype='|S8', count=3) 

                f.readline()  # strip the newline 

            #-------------------------------------------------------------- 

            # read in the seismogram points 

            #-------------------------------------------------------------- 

            self.seis = loadtxt(f, dtype='<f4', ndlim=1).ravel() 

        except IOError, e: 

            self.hf = self.hs = self.hi = self.seis = None 

            f.close() 

            raise SacIOError("%s is not a valid SAC file:" % fname, e) 

        try: 

            self._get_date() 

        except SacError: 

            warnings.warn('Cannot determine date') 

        if self.GetHvalue('lcalda'): 

            try: 

                self._get_dist() 

            except SacError: 

                pass 

        f.close() 

 

    def ReadSacXYHeader(self, fname): 

        """ 

        Read SAC XY files (ascii) 

 

        :param f: filename (SAC ascii). 

 

        >>> from obspy.sac import SacIO # doctest: +SKIP 

        >>> tr = SacIO() # doctest: +SKIP 

        >>> tr.ReadSacXY('testxy.sac') # doctest: +SKIP 

        >>> tr.GetHvalue('npts') # doctest: +SKIP 

        100 

 

        This is equivalent to: 

 

        >>> tr = SacIO('testxy.sac',alpha=True)  # doctest: +SKIP 

 

        Reading only the header portion of alphanumeric SAC-files is currently 

        not supported. 

        """ 

        ###### open the file 

        try: 

            f = open(fname, 'r') 

        except IOError: 

            raise SacIOError("No such file: " + fname) 

        try: 

            #-------------------------------------------------------------- 

            # parse the header 

            # 

            # The sac header has 70 floats, 40 integers, then 192 bytes 

            #    in strings. Store them in array (an convert the char to a 

            #    list). 

            #-------------------------------------------------------------- 

            # read in the float values 

            self.hf = np.fromfile(f, dtype='<f4', count=70, sep=" ") 

            # read in the int values 

            self.hi = np.fromfile(f, dtype='<i4', count=40, sep=" ") 

            # reading in the string part is a bit more complicated 

            # because every string field has to be 8 characters long 

            # apart from the second field which is 16 characters long 

            # resulting in a total length of 192 characters 

            for i in xrange(0, 24, 3): 

                self.hs[i:i + 3] = np.fromfile(f, dtype='|S8', count=3) 

                f.readline()  # strip the newline 

        except IOError, e: 

            self.hf = self.hs = self.hi = self.seis = None 

            f.close() 

            raise SacIOError("%s is not a valid SAC file:" % fname, e) 

        try: 

            self.IsSACfile(fname, fsize=False) 

        except SacError, e: 

            f.close() 

            raise SacError(e) 

        try: 

            self._get_date() 

        except SacError: 

            warnings.warn('Cannot determine date') 

        if self.GetHvalue('lcalda'): 

            try: 

                self._get_dist() 

            except SacError: 

                pass 

        f.close() 

 

    def readTrace(self, trace): 

        """ 

        Fill in SacIO object with data from obspy trace. 

        Warning: Currently only the case of a previously empty SacIO object is 

        safe! 

        """ 

        # extracting relative SAC time as specified with b 

        try: 

            b = float(trace.stats['sac']['b']) 

        except KeyError: 

            b = 0.0 

        # filling in SAC/sacio specific defaults 

        self.fromarray(trace.data, begin=b, delta=trace.stats.delta, 

                       starttime=trace.stats.starttime) 

        # overwriting with ObsPy defaults 

        for _j, _k in convert_dict.iteritems(): 

            self.SetHvalue(_j, trace.stats[_k]) 

        # overwriting up SAC specific values 

        # note that the SAC reference time values (including B and E) are 

        # not used in here any more, they are already set by t.fromarray 

        # and directly deduce from tr.starttime 

        for _i in SAC_EXTRA: 

            try: 

                self.SetHvalue(_i, trace.stats.sac[_i]) 

            except KeyError: 

                pass 

        return 

 

    def WriteSacXY(self, ofname): 

        """ 

        Write SAC XY file (ascii) 

 

        :param f: filename (SAC ascii) 

 

        >>> from obspy.sac import SacIO # doctest: +SKIP 

        >>> tr = SacIO('test.sac') # doctest: +SKIP 

        >>> tr.WriteSacXY('test2.sac') # doctest: +SKIP 

        >>> tr.IsValidXYSacFile('test2.sac') # doctest: +SKIP 

        True 

        """ 

        try: 

            f = open(ofname, 'w') 

        except IOError: 

            raise SacIOError("Cannot open file: " + ofname) 

        # header 

        try: 

            np.savetxt(f, np.reshape(self.hf, (14, 5)), 

                       fmt="%#15.7g%#15.7g%#15.7g%#15.7g%#15.7g") 

            np.savetxt(f, np.reshape(self.hi, (8, 5)), 

                       fmt="%10d%10d%10d%10d%10d") 

            for i in xrange(0, 24, 3): 

                self.hs[i:i + 3].tofile(f) 

                f.write('\n') 

        except: 

            f.close() 

            raise SacIOError("Cannot write header values: " + ofname) 

        # traces 

        npts = self.GetHvalue('npts') 

        if npts == -12345: 

            return 

        try: 

            rows = npts / 5 

            np.savetxt(f, np.reshape(self.seis[0:5 * rows], (rows, 5)), 

                       fmt="%#15.7g%#15.7g%#15.7g%#15.7g%#15.7g") 

            np.savetxt(f, self.seis[5 * rows:], delimiter='\t') 

        except: 

            f.close() 

            raise SacIOError("Cannot write trace values: " + ofname) 

        f.close() 

 

    def WriteSacBinary(self, ofname): 

        """ 

        Write a SAC binary file using the head arrays and array seis. 

 

        :param f: filename (SAC binary). 

 

        >>> from obspy.sac import SacIO # doctest: +SKIP 

        >>> tr = SacIO('test.sac') # doctest: +SKIP 

        >>> tr.WriteSacBinary('test2.sac') # doctest: +SKIP 

        >>> os.stat('test2.sac')[6] == os.stat('test.sac')[6] # doctest: +SKIP 

        True 

        """ 

        try: 

            f = open(ofname, 'wb+') 

        except IOError: 

            raise SacIOError("Cannot open file: " + ofname) 

        try: 

            self._chck_header() 

            self.hf.tofile(f) 

            self.hi.tofile(f) 

            self.hs.tofile(f) 

            self.seis.tofile(f) 

        except Exception, e: 

            f.close() 

            msg = "Cannot write SAC-buffer to file: " 

            raise SacIOError(msg, ofname, e) 

 

        f.close() 

 

    def PrintIValue(self, label='=', value=-12345): 

        """ 

        Convenience function for printing undefined integer header values. 

        """ 

        if value != -12345: 

            print(label, value) 

 

    def PrintFValue(self, label='=', value=-12345.0): 

        """ 

        Convenience function for printing undefined float header values. 

        """ 

        if value != -12345.0: 

            print('%s %.8g' % (label, value)) 

 

    def PrintSValue(self, label='=', value='-12345'): 

        """ 

        Convenience function for printing undefined string header values. 

        """ 

        if value.find('-12345') == -1: 

            print(label, value) 

 

    def ListStdValues(self):  # h is a header list, s is a float list 

        """ 

        Convenience function for printing common header values. 

 

        :param: None 

        :return: None 

 

        >>> from obspy.sac import SacIO  # doctest: +SKIP 

        >>> t = SacIO('test.sac')  # doctest: +SKIP 

        >>> t.ListStdValues()  # doctest: +SKIP +NORMALIZE_WHITESPACE 

        <BLANKLINE> 

        Reference Time = 07/18/1978 (199) 8:0:0.0 

        Npts  =  100 

        Delta =  1 

        Begin =  10 

        End   =  109 

        Min   =  -1 

        Mean  =  8.7539462e-08 

        Max   =  1 

        Header Version =  6 

        Station =  STA 

        Channel =  Q 

        Event       =  FUNCGEN: SINE 

 

        If no header values are defined (i.e. all are equal 12345) than this 

        function won't do anything. 

        """ 

        # 

        # Seismogram Info: 

        # 

        try: 

            nzyear = self.GetHvalue('nzyear') 

            nzjday = self.GetHvalue('nzjday') 

            month = time.strptime(repr(nzyear) + " " + repr(nzjday), 

                                  "%Y %j").tm_mon 

            date = time.strptime(repr(nzyear) + " " + repr(nzjday), 

                                 "%Y %j").tm_mday 

            pattern = '\nReference Time = %2.2d/%2.2d/%d (%d) %d:%d:%d.%d' 

            print(pattern % (month, date, 

                             self.GetHvalue('nzyear'), 

                             self.GetHvalue('nzjday'), 

                             self.GetHvalue('nzhour'), 

                             self.GetHvalue('nzmin'), 

                             self.GetHvalue('nzsec'), 

                             self.GetHvalue('nzmsec'))) 

        except ValueError: 

            pass 

        self.PrintIValue('Npts  = ', self.GetHvalue('npts')) 

        self.PrintFValue('Delta = ', self.GetHvalue('delta')) 

        self.PrintFValue('Begin = ', self.GetHvalue('b')) 

        self.PrintFValue('End   = ', self.GetHvalue('e')) 

        self.PrintFValue('Min   = ', self.GetHvalue('depmin')) 

        self.PrintFValue('Mean  = ', self.GetHvalue('depmen')) 

        self.PrintFValue('Max   = ', self.GetHvalue('depmax')) 

        # 

        self.PrintIValue('Header Version = ', self.GetHvalue('nvhdr')) 

        # 

        # station Info: 

        # 

        self.PrintSValue('Station = ', self.GetHvalue('kstnm')) 

        self.PrintSValue('Channel = ', self.GetHvalue('kcmpnm')) 

        self.PrintFValue('Station Lat  = ', self.GetHvalue('stla')) 

        self.PrintFValue('Station Lon  = ', self.GetHvalue('stlo')) 

        self.PrintFValue('Station Elev = ', self.GetHvalue('stel')) 

        # 

        # Event Info: 

        # 

        self.PrintSValue('Event       = ', self.GetHvalue('kevnm')) 

        self.PrintFValue('Event Lat   = ', self.GetHvalue('evla')) 

        self.PrintFValue('Event Lon   = ', self.GetHvalue('evlo')) 

        self.PrintFValue('Event Depth = ', self.GetHvalue('evdp')) 

        self.PrintFValue('Origin Time = ', self.GetHvalue('o')) 

        # 

        self.PrintFValue('Azimuth        = ', self.GetHvalue('az')) 

        self.PrintFValue('Back Azimuth   = ', self.GetHvalue('baz')) 

        self.PrintFValue('Distance (km)  = ', self.GetHvalue('dist')) 

        self.PrintFValue('Distance (deg) = ', self.GetHvalue('gcarc')) 

 

    def GetHvalueFromFile(self, thePath, theItem): 

        """ 

        Quick access to a specific header item in specified file. 

 

        :param f: filename (SAC binary) 

        :type hn: string 

        :param hn: header variable name 

 

        >>> from obspy.sac import SacIO # doctest: +SKIP 

        >>> t = SacIO() # doctest: +SKIP 

        >>> t.GetHvalueFromFile('test.sac','kcmpnm').rstrip() # doctest: +SKIP 

        'Q' 

 

        String header values have a fixed length of 8 or 16 characters. This 

        can lead to errors for example if you concatenate strings and forget to 

        strip off the trailing whitespace. 

        """ 

        # 

        #  Read in the Header 

        # 

        self.ReadSacHeader(thePath) 

        # 

        return(self.GetHvalue(theItem)) 

 

    def SetHvalueInFile(self, thePath, theItem, theValue): 

        """ 

        Quick access to change a specific header item in a specified file. 

 

        :param f: filename (SAC binary) 

        :type hn: string 

        :param hn: header variable name 

        :type hv: string, float or integer 

        :param hv: header variable value (numeric or string value to be 

            assigned to hn) 

        :return: None 

 

        >>> from obspy.sac import SacIO # doctest: +SKIP 

        >>> t = SacIO() # doctest: +SKIP 

        >>> t.GetHvalueFromFile('test.sac','kstnm').rstrip() # doctest: +SKIP 

        'STA' 

        >>> t.SetHvalueInFile('test.sac','kstnm','blub') # doctest: +SKIP 

        >>> t.GetHvalueFromFile('test.sac','kstnm').rstrip() # doctest: +SKIP 

        'blub' 

        """ 

        # 

        #  Read in the Header 

        # 

        self.ReadSacHeader(thePath) 

        # 

        self.SetHvalue(theItem, theValue) 

        self.WriteSacHeader(thePath) 

 

    def IsValidSacFile(self, thePath): 

        """ 

        Quick test for a valid SAC binary file (wraps 'IsSACfile'). 

 

        :param f: filename (SAC binary) 

        :rtype: boolean (True or False) 

 

        >>> from obspy.sac import SacIO # doctest: +SKIP 

        >>> SacIO().IsValidSacFile('test.sac') # doctest: +SKIP 

        True 

        >>> SacIO().IsValidSacFile('testxy.sac') # doctest: +SKIP 

        False 

        """ 

        # 

        #  Read in the Header 

        # 

        try: 

            self.ReadSacHeader(thePath) 

        except SacError: 

            return False 

        except SacIOError: 

            return False 

        else: 

            return True 

 

    def IsValidXYSacFile(self, filename): 

        """ 

        Quick test for a valid SAC ascii file. 

 

        :param file: filename (SAC ascii) 

        :rtype: boolean (True or False) 

 

        >>> from obspy.sac import SacIO # doctest: +SKIP 

        >>> SacIO().IsValidXYSacFile('testxy.sac') # doctest: +SKIP 

        True 

        >>> SacIO().IsValidXYSacFile('test.sac') # doctest: +SKIP 

        False 

        """ 

        # 

        #  Read in the Header 

        # 

        if _isText(filename): 

            try: 

                self.ReadSacXY(filename) 

            except: 

                return False 

            return True 

        else: 

            return False 

 

    def _get_date(self): 

        """ 

        If date header values are set calculate date in julian seconds 

 

        >>> t = SacIO() 

        >>> t.fromarray(np.random.randn(100), delta=1.0, \ 

                        starttime=UTCDateTime(1970,01,01)) 

        >>> t._get_date() 

        >>> t.reftime.timestamp 

        0.0 

        >>> t.endtime.timestamp - t.reftime.timestamp 

        100.0 

        """ 

        ### if any of the time-header values are still set to 

        ### -12345 then UTCDateTime raises an exception and 

        ### reftime is set to 0.0 

        try: 

            ms = self.GetHvalue('nzmsec') * 1000 

            self.reftime = UTCDateTime(year=self.GetHvalue('nzyear'), 

                                         julday=self.GetHvalue('nzjday'), 

                                         hour=self.GetHvalue('nzhour'), 

                                         minute=self.GetHvalue('nzmin'), 

                                         second=self.GetHvalue('nzsec'), 

                                         microsecond=ms) 

            b = float(self.GetHvalue('b')) 

            if b != -12345.0: 

                self.starttime = self.reftime + b 

            else: 

                self.starttime = self.reftime 

            self.endtime = self.starttime + \ 

                self.GetHvalue('npts') * float(self.GetHvalue('delta')) 

        except: 

            try: 

                self.reftime = UTCDateTime(0.0) 

                b = float(self.GetHvalue('b')) 

                if b != -12345.0: 

                    self.starttime = self.reftime + b 

                else: 

                    self.starttime = self.reftime 

                self.endtime = self.reftime + \ 

                    self.GetHvalue('npts') * float(self.GetHvalue('delta')) 

            except: 

                raise SacError("Cannot calculate date") 

 

    def _chck_header(self): 

        """ 

        If trace changed since read, adapt header values 

        """ 

        self.seis = np.require(self.seis, '<f4') 

        self.SetHvalue('npts', self.seis.size) 

        if self.seis.size == 0: 

            return 

        self.SetHvalue('depmin', self.seis.min()) 

        self.SetHvalue('depmax', self.seis.max()) 

        self.SetHvalue('depmen', self.seis.mean()) 

 

    def _get_dist(self): 

        """ 

        calculate distance from station and event coordinates 

 

        >>> t = SacIO() 

        >>> t.SetHvalue('evla',48.15) 

        >>> t.SetHvalue('evlo',11.58333) 

        >>> t.SetHvalue('stla',-41.2869) 

        >>> t.SetHvalue('stlo',174.7746) 

        >>> t._get_dist() 

        >>> print(round(t.GetHvalue('dist'), 2)) 

        18486.53 

        >>> print(round(t.GetHvalue('az'), 5)) 

        65.65415 

        >>> print(round(t.GetHvalue('baz'), 4)) 

        305.9755 

 

        The original SAC-program calculates the distance assuming a 

        average radius of 6371 km. Therefore, our routine should be more 

        accurate. 

        """ 

        eqlat = self.GetHvalue('evla') 

        eqlon = self.GetHvalue('evlo') 

        stlat = self.GetHvalue('stla') 

        stlon = self.GetHvalue('stlo') 

        d = self.GetHvalue('dist') 

        if eqlat == -12345.0 or eqlon == -12345.0 or \ 

           stlat == -12345.0 or stlon == -12345.0: 

            raise SacError('Insufficient information to calculate distance.') 

        if d != -12345.0: 

            raise SacError('Distance is already set.') 

        dist, az, baz = gps2DistAzimuth(eqlat, eqlon, stlat, stlon) 

        self.SetHvalue('dist', dist / 1000.) 

        self.SetHvalue('az', az) 

        self.SetHvalue('baz', baz) 

 

    def swap_byte_order(self): 

        """ 

        Swap byte order of SAC-file in memory. 

 

        Currently seems to work only for conversion from big-endian to 

        little-endian. 

 

        :param: None 

        :return: None 

 

        >>> from obspy.sac import SacIO # doctest: +SKIP 

        >>> t = SacIO('test.sac') # doctest: +SKIP 

        >>> t.swap_byte_order() # doctest: +SKIP 

        """ 

        if self.byteorder == 'big': 

            bs = 'L' 

        elif self.byteorder == 'little': 

            bs = 'B' 

        self.seis.byteswap(True) 

        self.hf.byteswap(True) 

        self.hi.byteswap(True) 

        self.seis = self.seis.newbyteorder(bs) 

        self.hf = self.hf.newbyteorder(bs) 

        self.hi = self.hi.newbyteorder(bs) 

 

    def __getattr__(self, hname): 

        """ 

        convenience function to access header values 

 

        :param hname: header variable name 

 

        >>> tr = SacIO() 

        >>> tr.fromarray(np.random.randn(100)) 

        >>> tr.npts == tr.GetHvalue('npts') # doctest: +SKIP 

        True 

        """ 

        return self.GetHvalue(hname) 

 

    def get_obspy_header(self): 

        """ 

        Return a dictionary that can be used as a header in creating a new 

        :class:`~obspy.core.trace.Trace` object. 

        Currently most likely an Exception will be raised if no SAC file was 

        read beforehand! 

        """ 

        header = {} 

        # convert common header types of the ObsPy trace object 

        for i, j in convert_dict.iteritems(): 

            value = self.GetHvalue(i) 

            if isinstance(value, str): 

                null_term = value.find('\x00') 

                if null_term >= 0: 

                    value = value[:null_term] 

                value = value.strip() 

                if value == '-12345': 

                    value = '' 

            # fix for issue #156 

            if i == 'delta': 

                header['sampling_rate'] = \ 

                        np.float32(1.0) / np.float32(self.hf[0]) 

            else: 

                header[j] = value 

        if header['calib'] == -12345.0: 

            header['calib'] = 1.0 

        # assign extra header types of SAC 

        header['sac'] = {} 

        for i in SAC_EXTRA: 

            header['sac'][i] = self.GetHvalue(i) 

        # convert time to UTCDateTime 

        header['starttime'] = self.starttime 

        # always add the begin time (if it's defined) to get the given 

        # SAC reference time, no matter which iztype is given 

        # note that the B and E times should not be in the SAC_EXTRA 

        # dictionary, as they would overwrite the self.fromarray which sets 

        # them according to the starttime, npts and delta. 

        header['sac']['b'] = float(self.GetHvalue('b')) 

        header['sac']['e'] = float(self.GetHvalue('e')) 

        # ticket #390 

        if self.debug_headers: 

            for i in ['nzyear', 'nzjday', 'nzhour', 'nzmin', 'nzsec', 'nzmsec', 

                      'delta', 'scale', 'npts', 'knetwk', 'kstnm', 'kcmpnm']: 

                header['sac'][i] = self.GetHvalue(i) 

        return header 

 

 

############# UTILITIES ################################################### 

def attach_paz(tr, paz_file, todisp=False, tovel=False, torad=False, 

               tohz=False): 

    ''' 

    Attach tr.stats.paz AttribDict to trace from SAC paz_file 

 

    This is experimental code, taken from 

    obspy.gse2.libgse2.attach_paz and adapted to the SAC-pole-zero 

    conventions. Especially the conversion from velocity to 

    displacement and vice versa is still under construction. It works 

    but I cannot guarantee that the values are correct. For more 

    information on the SAC-pole-zero format see: 

    http://www.iris.edu/software/sac/commands/transfer.html. For a 

    useful discussion on polezero files and transfer functions in 

    general see: 

    http://www.le.ac.uk/seis-uk/downloads/seisuk_instrument_resp_removal.pdf. 

    Also bear in mind that according to the SAC convention for 

    pole-zero files CONSTANT is defined as: 

    digitizer_gain*seismometer_gain*A0. This means that it does not 

    have explicit information on the digitizer gain and seismometer 

    gain which we therefore set to 1.0. 

 

    Attaches to a trace a paz AttribDict containing poles zeros and gain. 

 

    :param tr: An ObsPy :class:`~obspy.core.trace.Trace` object 

    :param paz_file: path to pazfile or file pointer 

    :param todisp: change a velocity transfer function to a displacement 

                   transfer function by adding another zero 

    :param tovel: change a displacement transfer function to a velocity 

                  transfer function by removing one 0,0j zero 

    :param torad: change to radians 

    :param tohz: change to Hertz 

 

    >>> from obspy import Trace 

    >>> tr = Trace() 

    >>> import StringIO 

    >>> f = StringIO.StringIO("""ZEROS 3 

    ... -5.032 0.0 

    ... POLES 6 

    ... -0.02365 0.02365 

    ... -0.02365 -0.02365 

    ... -39.3011 0. 

    ... -7.74904 0. 

    ... -53.5979 21.7494 

    ... -53.5979 -21.7494 

    ... CONSTANT 2.16e18""") 

    >>> attach_paz(tr, f,torad=True) 

    >>> print(tr.stats.paz['zeros'][0]) 

    (-31.6169884657+0j) 

    ''' 

 

    poles = [] 

    zeros = [] 

    found_zero = False 

 

    if isinstance(paz_file, str): 

        paz_file = open(paz_file, 'r') 

 

    while True: 

        line = paz_file.readline() 

        if not line: 

            break 

        ### lines starting with * are comments 

        if line.startswith('*'): 

            continue 

        if line.find('ZEROS') != -1: 

            a = line.split() 

            noz = int(a[1]) 

            for _k in xrange(noz): 

                line = paz_file.readline() 

                a = line.split() 

                if line.find('POLES') != -1 or line.find('CONSTANT') != -1 or \ 

                   line.startswith('*') or not line: 

                    while len(zeros) < noz: 

                        zeros.append(complex(0, 0j)) 

                    break 

                else: 

                    zeros.append(complex(float(a[0]), float(a[1]))) 

 

        if line.find('POLES') != -1: 

            a = line.split() 

            nop = int(a[1]) 

            for _k in xrange(nop): 

                line = paz_file.readline() 

                a = line.split() 

                if line.find('CONSTANT') != -1 or line.find('ZEROS') != -1 or \ 

                   line.startswith('*') or not line: 

                    while len(poles) < nop: 

                        poles.append(complex(0, 0j)) 

                    break 

                else: 

                    poles.append(complex(float(a[0]), float(a[1]))) 

        if line.find('CONSTANT') != -1: 

            a = line.split() 

            # in the observatory this is the seismometer gain [muVolt/nm/s] 

            # the A0_normalization_factor is hardcoded to 1.0 

            constant = float(a[1]) 

    paz_file.close() 

 

    ### To convert the velocity response to the displacement response, 

    ### multiplication with jw is used. This is equivalent to one more 

    ### zero in the pole-zero representation 

    if todisp: 

        zeros.append(complex(0, 0j)) 

 

    ### To convert the displacement response to the velocity response, 

    ### division with jw is used. This is equivalent to one less zero 

    ### in the pole-zero representation 

    if tovel: 

        for i, zero in enumerate(list(zeros)): 

            if zero == complex(0, 0j): 

                zeros.pop(i) 

                found_zero = True 

                break 

        if not found_zero: 

            raise Exception("Could not remove (0,0j) zero to change \ 

            displacement response to velocity response") 

 

    ### convert poles, zeros and gain in Hertz to radians 

    if torad: 

        tmp = [z * 2. * np.pi for z in zeros] 

        zeros = tmp 

        tmp = [p * 2. * np.pi for p in poles] 

        poles = tmp 

        # When extracting RESP files and SAC_PZ files 

        # from a dataless SEED using the rdseed program 

        # where the former is in Hz and the latter in radians, 

        # there gains seem to be unaffected by this. 

        # According to this document: 

        # http://www.le.ac.uk/ 

        #         seis-uk/downloads/seisuk_instrument_resp_removal.pdf 

        # the gain should also be converted when changing from 

        # hertz to radians or vice versa. However, the rdseed programs 

        # does not do this. I'm not entirely sure at this stage which one is 

        # correct or if I have missed something. I've therefore decided 

        # to leave it out for now, in order to stay compatible with the 

        # rdseed program and the SAC program. 

        #constant *= (2. * np.pi) ** 3 

 

    ### convert poles, zeros and gain in radian to Hertz 

    if tohz: 

        for i, z in enumerate(zeros): 

            if abs(z) > 0.0: 

                zeros[i] /= 2 * np.pi 

        for i, p in enumerate(poles): 

            if abs(p) > 0.0: 

                poles[i] /= 2 * np.pi 

        #constant /= (2. * np.pi) ** 3 

 

    # fill up ObsPy Poles and Zeros AttribDict 

    # In SAC pole-zero files CONSTANT is defined as: 

    # digitizer_gain*seismometer_gain*A0 

 

    tr.stats.paz = AttribDict() 

    tr.stats.paz.seismometer_gain = 1.0 

    tr.stats.paz.digitizer_gain = 1.0 

    tr.stats.paz.poles = poles 

    tr.stats.paz.zeros = zeros 

    # taken from obspy.gse2.paz:145 

    tr.stats.paz.sensitivity = tr.stats.paz.digitizer_gain * \ 

        tr.stats.paz.seismometer_gain 

    tr.stats.paz.gain = constant 

 

 

def attach_resp(tr, resp_file, todisp=False, tovel=False, torad=False, 

               tohz=False): 

    """ 

    Extract key instrument response information from a RESP file, which 

    can be extracted from a dataless SEED volume by, for example, using 

    the script obspy-dataless2resp or the rdseed program. At the moment, 

    you have to determine yourself if the given response is for velocity 

    or displacement and if the values are given in rad or Hz. This is 

    still experimental code (see also documentation for 

    :func:`obspy.sac.sacio.attach_paz`). 

    Attaches to a trace a paz AttribDict containing poles, zeros, and gain. 

 

    :param tr: An ObsPy :class:`~obspy.core.trace.Trace` object 

    :param resp_file: path to RESP-file or file pointer 

    :param todisp: change a velocity transfer function to a displacement 

                   transfer function by adding another zero 

    :param tovel: change a displacement transfer function to a velocity 

                  transfer function by removing one 0,0j zero 

    :param torad: change to radians 

    :param tohz: change to Hertz 

 

    >>> from obspy import Trace 

    >>> tr = Trace() 

    >>> respfile = os.path.join(os.path.dirname(__file__), 'tests', 'data', 

    ...                         'RESP.NZ.CRLZ.10.HHZ') 

    >>> attach_resp(tr, respfile, torad=True, todisp=False) 

    >>> print tr.stats.paz.keys()  # doctest: +NORMALIZE_WHITESPACE 

    ['sensitivity', 'digitizer_gain', 'seismometer_gain', 'zeros', 'gain', 

     't_shift', 'poles'] 

    >>> print tr.stats.paz.poles  # doctest: +SKIP 

    [(-0.15931644664884559+0.15931644664884559j), 

     (-0.15931644664884559-0.15931644664884559j), 

     (-314.15926535897933+202.31856689118268j), 

     (-314.15926535897933-202.31856689118268j)] 

    """ 

    if isinstance(resp_file, str): 

        resp_file = open(resp_file, 'r') 

 

    zeros_pat = r'B053F10-13' 

    poles_pat = r'B053F15-18' 

    a0_pat = r'B053F07' 

    sens_pat = r'B058F04' 

    t_shift_pat = r'B057F08' 

    t_shift = 0.0 

    poles = [] 

    zeros = [] 

    while True: 

        line = resp_file.readline() 

        if not line: 

            break 

        if line.startswith(a0_pat): 

            a0 = float(line.split(':')[1]) 

        if line.startswith(sens_pat): 

            sens = float(line.split(':')[1]) 

        if line.startswith(poles_pat): 

            tmp = line.split() 

            poles.append(complex(float(tmp[2]), float(tmp[3]))) 

        if line.startswith(zeros_pat): 

            tmp = line.split() 

            zeros.append(complex(float(tmp[2]), float(tmp[3]))) 

        if line.startswith(t_shift_pat): 

            t_shift += float(line.split(':')[1]) 

    constant = a0 * sens 

 

    if torad: 

        tmp = [z * 2. * np.pi for z in zeros] 

        zeros = tmp 

        tmp = [p * 2. * np.pi for p in poles] 

        poles = tmp 

 

    if todisp: 

        zeros.append(complex(0, 0j)) 

 

    ### To convert the displacement response to the velocity response, 

    ### division with jw is used. This is equivalent to one less zero 

    ### in the pole-zero representation 

    if tovel: 

        for i, zero in enumerate(list(zeros)): 

            if zero == complex(0, 0j): 

                zeros.pop(i) 

                found_zero = True 

                break 

        if not found_zero: 

            raise Exception("Could not remove (0,0j) zero to change \ 

            displacement response to velocity response") 

 

    ### convert poles, zeros and gain in radian to Hertz 

    if tohz: 

        for i, z in enumerate(zeros): 

            if abs(z) > 0.0: 

                zeros[i] /= 2 * np.pi 

        for i, p in enumerate(poles): 

            if abs(p) > 0.0: 

                poles[i] /= 2 * np.pi 

        constant /= (2. * np.pi) ** 3 

 

    # fill up ObsPy Poles and Zeros AttribDict 

    # In SAC pole-zero files CONSTANT is defined as: 

    # digitizer_gain*seismometer_gain*A0 

 

    tr.stats.paz = AttribDict() 

    tr.stats.paz.seismometer_gain = sens 

    tr.stats.paz.digitizer_gain = 1.0 

    tr.stats.paz.poles = poles 

    tr.stats.paz.zeros = zeros 

    # taken from obspy.gse2.paz:145 

    tr.stats.paz.sensitivity = tr.stats.paz.digitizer_gain * \ 

        tr.stats.paz.seismometer_gain 

    tr.stats.paz.gain = constant 

    tr.stats.paz.t_shift = t_shift 

 

 

if __name__ == "__main__": 

    import doctest 

    doctest.testmod()