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

#!/usr/bin/env python 

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

# Filename: invsim.py 

#  Purpose: Python Module for Instrument Correction (Seismology) 

#   Author: Moritz Beyreuther, Yannik Behr 

#    Email: moritz.beyreuther@geophysik.uni-muenchen.de 

# 

# Copyright (C) 2008-2012 Moritz Beyreuther, Yannik Behr 

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

""" 

Python Module for Instrument Correction (Seismology). 

PAZ (Poles and zeros) information must be given in SEED convention, correction 

to m/s. 

 

:copyright: 

    The ObsPy Development Team (devs@obspy.org) 

:license: 

    GNU Lesser General Public License, Version 3 

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

""" 

 

from obspy.core.util.base import NamedTemporaryFile 

from obspy.core.util.decorator import deprecated_keywords 

from obspy.signal.detrend import simple as simpleDetrend 

from obspy.signal.headers import clibevresp 

import ctypes as C 

import math as M 

import numpy as np 

import os 

import scipy.signal 

import util 

import warnings 

 

 

# Sensitivity is 2080 according to: 

# P. Bormann: New Manual of Seismological Observatory Practice 

# IASPEI Chapter 3, page 24 

# (PITSA has 2800) 

WOODANDERSON = {'poles': [-6.283 + 4.7124j, -6.283 - 4.7124j], 

                'zeros': [0 + 0j], 'gain': 1.0, 'sensitivity': 2080} 

 

 

def cosTaper(npts, p=0.1, freqs=None, flimit=None, halfcosine=True, 

             sactaper=False): 

    """ 

    Cosine Taper. 

 

    :type npts: Int 

    :param npts: Number of points of cosine taper. 

    :type p: Float 

    :param p: Decimal percentage of cosine taper (ranging from 0 to 1). Default 

        is 0.1 (10%) which tapers 5% from the beginning and 5% form the end. 

    :rtype: float NumPy ndarray 

    :return: Cosine taper array/vector of length npts. 

    :type freqs: NumPy ndarray 

    :param freqs: Frequencies as, for example, returned by fftfreq 

    :type flimit: List or tuple of floats 

    :param flimit: The list or tuple defines the four corner frequencies 

        (f1, f2, f3, f4) of the cosine taper which is one between f2 and f3 and 

        tapers to zero for f1 < f < f2 and f3 < f < f4. 

    :type halfcosine: Boolean 

    :param halfcosine: If True the taper is a half cosine function. If False it 

        is a quarter cosine function. 

    :type sactaper: Boolean 

    :param sactaper: If set to True the cosine taper already tapers at the 

        corner frequency (SAC behaviour). By default, the taper has a value 

        of 1.0 at the corner frequencies. 

 

    .. rubric:: Example 

 

    >>> tap = cosTaper(100, 1.0) 

    >>> tap2 = 0.5 * (1 + np.cos(np.linspace(np.pi, 2 * np.pi, 50))) 

    >>> np.allclose(tap[0:50], tap2) 

    True 

    >>> npts = 100 

    >>> p = 0.1 

    >>> tap3 = cosTaper(npts, p) 

    >>> ( tap3[npts*p/2.:npts*(1-p/2.)]==np.ones(npts*(1-p)) ).all() 

    True 

    """ 

    if p < 0 or p > 1: 

        msg = "Decimal taper percentage must be between 0 and 1." 

        raise ValueError(msg) 

    if p == 0.0 or p == 1.0: 

        frac = int(npts * p / 2.0) 

    else: 

        frac = int(npts * p / 2.0 + 0.5) 

 

    if freqs is not None and flimit is not None: 

        fl1, fl2, fl3, fl4 = flimit 

        idx1 = np.argmin(abs(freqs - fl1)) 

        idx2 = np.argmin(abs(freqs - fl2)) 

        idx3 = np.argmin(abs(freqs - fl3)) 

        idx4 = np.argmin(abs(freqs - fl4)) 

    else: 

        idx1 = 0 

        idx2 = frac - 1 

        idx3 = npts - frac 

        idx4 = npts - 1 

    if sactaper: 

        # in SAC the second and third 

        # index are already tapered 

        idx2 += 1 

        idx3 -= 1 

 

    # Very small data lengths or small decimal taper percentages can result in 

    # idx1 == idx2 and idx3 == idx4. This breaks the following calculations. 

    if idx1 == idx2: 

        idx2 += 1 

    if idx3 == idx4: 

        idx3 -= 1 

 

    # the taper at idx1 and idx4 equals zero and 

    # at idx2 and idx3 equals one 

    cos_win = np.zeros(npts) 

    if halfcosine: 

        #cos_win[idx1:idx2+1] =  0.5 * (1.0 + np.cos((np.pi * \ 

        #    (idx2 - np.arange(idx1, idx2+1)) / (idx2 - idx1)))) 

        cos_win[idx1:idx2 + 1] = 0.5 * (1.0 - np.cos((np.pi * \ 

            (np.arange(idx1, idx2 + 1) - idx1) / (idx2 - idx1)))) 

        cos_win[idx2 + 1:idx3] = 1.0 

        cos_win[idx3:idx4 + 1] = 0.5 * (1.0 + np.cos((np.pi * \ 

            (idx3 - np.arange(idx3, idx4 + 1)) / (idx4 - idx3)))) 

    else: 

        cos_win[idx1:idx2 + 1] = np.cos(-(np.pi / 2.0 * \ 

               (idx2 - np.arange(idx1, idx2 + 1)) / (idx2 - idx1))) 

        cos_win[idx2 + 1:idx3] = 1.0 

        cos_win[idx3:idx4 + 1] = np.cos((np.pi / 2.0 * \ 

            (idx3 - np.arange(idx3, idx4 + 1)) / (idx4 - idx3))) 

 

    # if indices are identical division by zero 

    # causes NaN values in cos_win 

    if idx1 == idx2: 

        cos_win[idx1] = 0.0 

    if idx3 == idx4: 

        cos_win[idx3] = 0.0 

    return cos_win 

 

 

def c_sac_taper(npts, p=0.1, freqs=None, flimit=None, pitsa=False): 

    twopi = 6.283185307179586 

    dblepi = 0.5 * twopi 

    fl1, fl2, fl3, fl4 = flimit 

    taper = [] 

    for freq in freqs: 

        if freq < fl3 and freq > fl2: 

            taper_v = 1.0 

        if freq >= fl3 and freq <= fl4: 

            taper_v = 0.5 * (1.0 + M.cos(dblepi * (freq - fl3) / (fl4 - fl3))) 

        if freq > fl4 or freq < fl1: 

            taper_v = 0.0 

        if freq >= fl1 and freq <= fl2: 

            taper_v = 0.5 * (1.0 - M.cos(dblepi * (freq - fl1) / (fl2 - fl1))) 

        taper.append(taper_v) 

    return np.array(taper) 

 

 

@deprecated_keywords({'pitsa': None}) 

def evalresp(t_samp, nfft, filename, date, station='*', channel='*', 

             network='*', locid='*', units="VEL", freq=False, 

             debug=False): 

    """ 

    Use the evalresp library to extract instrument response 

    information from a SEED RESP-file. 

 

    :type t_samp: float 

    :param t_samp: Sampling interval in seconds 

    :type nfft: int 

    :param nfft: Number of FFT points of signal which needs correction 

    :type filename: str 

    :param filename: SEED RESP-filename or content of RESP file 

    :type date: UTCDateTime 

    :param date: Date of interest 

    :type station: str 

    :param station: Station id 

    :type channel: str 

    :param channel: Channel id 

    :type network: str 

    :param network: Network id 

    :type locid: str 

    :param locid: Location id 

    :type units: str 

    :param units: Units to return response in. Can be either DIS, VEL or ACC 

    :type debug: bool 

    :param debug: Verbose output to stdout. Disabled by default. 

    :rtype: numpy.ndarray complex128 

    :return: Frequency response from SEED RESP-file of length nfft 

    """ 

    # evalresp needs files with correct line separators depending on OS 

    data = open(filename, 'rb').read() 

    fh = NamedTemporaryFile() 

    tempfile = fh.name 

    fh.write(os.linesep.join(data.splitlines())) 

    fh.close() 

 

    fy = 1 / (t_samp * 2.0) 

    # start at zero to get zero for offset/ DC of fft 

    freqs = np.linspace(0, fy, nfft // 2 + 1) 

    start_stage = C.c_int(-1) 

    stop_stage = C.c_int(0) 

    stdio_flag = C.c_int(0) 

    sta = C.create_string_buffer(station) 

    cha = C.create_string_buffer(channel) 

    net = C.create_string_buffer(network) 

    locid = C.create_string_buffer(locid) 

    unts = C.create_string_buffer(units) 

    if debug: 

        vbs = C.create_string_buffer("-v") 

    else: 

        vbs = C.create_string_buffer("") 

    rtyp = C.create_string_buffer("CS") 

    datime = C.create_string_buffer("%d,%3d" % (date.year, date.julday)) 

    fn = C.create_string_buffer(tempfile) 

    nfreqs = C.c_int(freqs.shape[0]) 

    res = clibevresp.evresp(sta, cha, net, locid, datime, unts, fn, 

                            freqs, nfreqs, rtyp, vbs, start_stage, 

                            stop_stage, stdio_flag, C.c_int(0)) 

    # optimizing performance, see 

    # http://wiki.python.org/moin/PythonSpeed/PerformanceTips 

    nfreqs, rfreqs, rvec = res[0].nfreqs, res[0].freqs, res[0].rvec 

    h = np.empty(nfreqs, dtype='complex128') 

    f = np.empty(nfreqs, dtype='float64') 

    for i in xrange(nfreqs): 

        h[i] = rvec[i].real + rvec[i].imag * 1j 

        f[i] = rfreqs[i] 

    clibevresp.free_response(res) 

    del nfreqs, rfreqs, rvec, res 

    # delete temporary file 

    try: 

        os.remove(tempfile) 

    except: 

        pass 

    if freq: 

        return h, f 

    return h 

 

 

def cornFreq2Paz(fc, damp=0.707): 

    """ 

    Convert corner frequency and damping to poles and zeros. 2 zeros at 

    position (0j, 0j) are given as output  (m/s). 

 

    :param fc: Corner frequency 

    :param damping: Corner frequency 

    :return: Dictionary containing poles, zeros and gain 

    """ 

    poles = [-(damp + M.sqrt(1 - damp ** 2) * 1j) * 2 * np.pi * fc] 

    poles.append(-(damp - M.sqrt(1 - damp ** 2) * 1j) * 2 * np.pi * fc) 

    return {'poles': poles, 'zeros': [0j, 0j], 'gain': 1, 'sensitivity': 1.0} 

 

 

@deprecated_keywords({'pitsa': None}) 

def pazToFreqResp(poles, zeros, scale_fac, t_samp, nfft, freq=False): 

    """ 

    Convert Poles and Zeros (PAZ) to frequency response. The output 

    contains the frequency zero which is the offset of the trace. 

 

    :type poles: List of complex numbers 

    :param poles: The poles of the transfer function 

    :type zeros: List of complex numbers 

    :param zeros: The zeros of the transfer function 

    :type scale_fac: Float 

    :param scale_fac: Gain factor 

    :type t_samp: Float 

    :param t_samp: Sampling interval in seconds 

    :type nfft: Integer 

    :param nfft: Number of FFT points of signal which needs correction 

    :rtype: numpy.ndarray complex128 

    :return: Frequency response of PAZ of length nfft 

 

    .. note:: 

        In order to plot/calculate the phase you need to multiply the 

        complex part by -1. This results from the different definition of 

        the Fourier transform and the phase. The numpy.fft is defined as 

        A(jw) = \int_{-\inf}^{+\inf} a(t) e^{-jwt}; where as the analytic 

        signal is defined A(jw) = | A(jw) | e^{j\phi}. That is in order to 

        calculate the phase the complex conjugate of the signal needs to be 

        taken. E.g. phi = angle(f,conj(h),deg=True) 

        As the range of phi is from -pi to pi you could add 2*pi to the 

        negative values in order to get a plot from [0, 2pi]: 

        where(phi<0,phi+2*pi,phi); plot(f,phi) 

    """ 

    n = nfft // 2 

    b, a = scipy.signal.ltisys.zpk2tf(zeros, poles, scale_fac) 

    # a has to be a list for the scipy.signal.freqs() call later but zpk2tf() 

    # strangely returns it as an integer. 

    if not isinstance(a, np.ndarray) and a == 1.0: 

        a = [1.0] 

    fy = 1 / (t_samp * 2.0) 

    # start at zero to get zero for offset / DC of fft 

    f = np.linspace(0, fy, n + 1) 

    _w, h = scipy.signal.freqs(b, a, f * 2 * np.pi) 

    if freq: 

        return h, f 

    return h 

 

 

def waterlevel(spec, wlev): 

    """ 

    Return the absolute spectral value corresponding 

    to dB wlev in spectrum spec. 

 

    :param spec: The spectrum 

    :param wlev: The water level 

    """ 

    return np.abs(spec).max() * 10.0 ** (-wlev / 20.0) 

 

 

def specInv(spec, wlev): 

    """ 

    Invert Spectrum and shrink values under water-level of max spec 

    amplitude. The water-level is given in db scale. 

 

    :note: In place operations on spec, translated from PITSA spr_sinv.c 

    :param spec: Spectrum as returned by numpy.fft.rfft 

    :param wlev: Water level to use 

    """ 

    # Calculated waterlevel in the scale of spec 

    swamp = waterlevel(spec, wlev) 

 

    # Find length in real fft frequency domain, spec is complex 

    sqrt_len = np.abs(spec) 

    # Set/scale length to swamp, but leave phase untouched 

    # 0 sqrt_len will transform in np.nans when dividing by it 

    idx = np.where((sqrt_len < swamp) & (sqrt_len > 0.0)) 

    spec[idx] *= swamp / sqrt_len[idx] 

    found = len(idx[0]) 

    # Now invert the spectrum for values where sqrt_len is greater than 

    # 0.0, see PITSA spr_sinv.c for details 

    sqrt_len = np.abs(spec)  # Find length of new scaled spec 

    inn = np.where(sqrt_len > 0.0) 

    spec[inn] = 1.0 / spec[inn] 

    # For numerical stability, set all zero length to zero, do not invert 

    spec[sqrt_len == 0.0] = complex(0.0, 0.0) 

    return found 

 

 

def seisSim(data, samp_rate, paz_remove=None, paz_simulate=None, 

            remove_sensitivity=True, simulate_sensitivity=True, 

            water_level=600.0, zero_mean=True, taper=True, 

            taper_fraction=0.05, pre_filt=None, seedresp=None, 

            nfft_pow2=False, pitsasim=True, sacsim=False, shsim=False, 

            **_kwargs): 

    """ 

    Simulate/Correct seismometer. 

 

    :type data: NumPy ndarray 

    :param data: Seismogram, detrend before hand (e.g. zero mean) 

    :type samp_rate: Float 

    :param samp_rate: Sample Rate of Seismogram 

    :type paz_remove: Dictionary, None 

    :param paz_remove: Dictionary containing keys 'poles', 'zeros', 'gain' 

        (A0 normalization factor). poles and zeros must be a list of complex 

        floating point numbers, gain must be of type float. Poles and Zeros are 

        assumed to correct to m/s, SEED convention. Use None for no inverse 

        filtering. 

    :type paz_simulate: Dictionary, None 

    :param paz_simulate: Dictionary containing keys 'poles', 'zeros', 'gain'. 

        Poles and zeros must be a list of complex floating point numbers, gain 

        must be of type float. Or None for no simulation. 

    :type remove_sensitivity: Boolean 

    :param remove_sensitivity: Determines if data is divided by 

        `paz_remove['sensitivity']` to correct for overall sensitivity of 

        recording instrument (seismometer/digitizer) during instrument 

        correction. 

    :type simulate_sensitivity: Boolean 

    :param simulate_sensitivity: Determines if data is multiplied with 

        `paz_simulate['sensitivity']` to simulate overall sensitivity of 

        new instrument (seismometer/digitizer) during instrument simulation. 

    :type water_level: Float 

    :param water_level: Water_Level for spectrum to simulate 

    :type zero_mean: Boolean 

    :param zero_mean: If true the mean of the data is subtracted 

    :type taper: Boolean 

    :param taper: If true a cosine taper is applied. 

    :type taper_fraction: Float 

    :param taper_fraction: Taper fraction of cosine taper to use 

    :type pre_filt: List or tuple of floats 

    :param pre_filt: Apply a bandpass filter to the data trace before 

        deconvolution. The list or tuple defines the four corner frequencies 

        (f1,f2,f3,f4) of a cosine taper which is one between f2 and f3 and 

        tapers to zero for f1 < f < f2 and f3 < f < f4. 

    :type seedresp: Dictionary, None 

    :param seedresp: Dictionary contains keys 'filename', 'date', 'units'. 

        'filename' is the path to a RESP-file generated from a dataless SEED 

        volume; 

        'date' is a `~obspy.core.utcdatetime.UTCDateTime` object for the date 

        that the response function should be extracted for; 

        'units' defines the units of the response function. 

        Can be either 'DIS', 'VEL' or 'ACC'. 

    :type nfft_pow2: Boolean 

    :param nfft_pow2: Number of frequency points to use for FFT. If True, 

        the exact power of two is taken (default in PITSA). If False the 

        data are not zeropadded to the next power of two which makes a 

        slower FFT but is then much faster for e.g. evalresp which scales 

        with the FFT points. 

    :type pitsasim: Boolean 

    :param pitsasim: Choose parameters to match 

        instrument correction as done by PITSA. 

    :type sacsim: Boolean 

    :param sacsim: Choose parameters to match 

        instrument correction as done by SAC. 

    :type shsim: Boolean 

    :param shsim: Choose parameters to match 

        instrument correction as done by Seismic Handler. 

    :return: The corrected data are returned as numpy.ndarray float64 

        array. float64 is chosen to avoid numerical instabilities. 

 

    This function works in the frequency domain, where nfft is the next power 

    of len(data) to avoid wrap around effects during convolution. The inverse 

    of the frequency response of the seismometer (``paz_remove``) is 

    convolved with the spectrum of the data and with the frequency response 

    of the seismometer to simulate (``paz_simulate``). A 5% cosine taper is 

    taken before simulation. The data must be detrended (e.g.) zero mean 

    beforehand. If paz_simulate=None only the instrument correction is done. 

    In the latter case, a broadband filter can be applied to the data trace 

    using pre_filt. This restricts the signal to the valid frequency band and 

    thereby avoids artefacts due to amplification of frequencies outside of the 

    instrument's passband (for a detailed discussion see 

    *Of Poles and Zeros*, F. Scherbaum, Kluwer Academic Publishers). 

 

    .. versionchanged:: 0.5.1 

        The default for `remove_sensitivity` and `simulate_sensitivity` has 

        been changed to ``True``. Old deprecated keyword arguments `paz`, 

        `inst_sim`, `no_inverse_filtering` have been removed. 

    """ 

    # Checking the types 

    if not paz_remove and not paz_simulate and not seedresp: 

        msg = "Neither inverse nor forward instrument simulation specified." 

        raise TypeError(msg) 

 

    for d in [paz_remove, paz_simulate]: 

        if d is None: 

            continue 

        for key in ['poles', 'zeros', 'gain']: 

            if key not in d: 

                raise KeyError("Missing key: %s" % key) 

    # Translated from PITSA: spr_resg.c 

    delta = 1.0 / samp_rate 

    # 

    ndat = len(data) 

    data = data.astype("float64") 

    if zero_mean: 

        data -= data.mean() 

    if taper: 

        if sacsim: 

            data *= cosTaper(ndat, taper_fraction, 

                             sactaper=sacsim, halfcosine=False) 

        else: 

            data *= cosTaper(ndat, taper_fraction) 

    # The number of points for the FFT has to be at least 2 * ndat (in 

    # order to prohibit wrap around effects during convolution) cf. 

    # Numerical Recipes p. 429 calculate next power of 2. 

    if nfft_pow2: 

        nfft = util.nextpow2(2 * ndat) 

    # evalresp scales directly with nfft, therefor taking the next power of 

    # two has a greater negative performance impact than the slow down of a 

    # not power of two in the FFT 

    elif ndat & 0x1:  # check if uneven 

        nfft = 2 * (ndat + 1) 

    else: 

        nfft = 2 * ndat 

    # Transform data in Fourier domain 

    data = np.fft.rfft(data, n=nfft) 

    # Inverse filtering = Instrument correction 

    if paz_remove: 

        freq_response, freqs = pazToFreqResp(paz_remove['poles'], 

                                             paz_remove['zeros'], 

                                             paz_remove['gain'], delta, nfft, 

                                             freq=True) 

    if seedresp: 

        freq_response, freqs = evalresp(delta, nfft, seedresp['filename'], 

                                        seedresp['date'], 

                                        units=seedresp['units'], freq=True) 

        if not remove_sensitivity: 

            msg = "remove_sensitivity is set to False, but since seedresp " + \ 

                  "is selected the overall sensitivity will be corrected " + \ 

                  " for anyway!" 

            warnings.warn(msg) 

    if paz_remove or seedresp: 

        if pre_filt: 

            # make cosine taper 

            fl1, fl2, fl3, fl4 = pre_filt 

            if sacsim: 

                cos_win = c_sac_taper(freqs.size, freqs=freqs, 

                                      flimit=(fl1, fl2, fl3, fl4)) 

            else: 

                cos_win = cosTaper(freqs.size, freqs=freqs, 

                                   flimit=(fl1, fl2, fl3, fl4)) 

            data *= cos_win 

        specInv(freq_response, water_level) 

        data *= freq_response 

        del freq_response 

    # Forward filtering = Instrument simulation 

    if paz_simulate: 

        data *= pazToFreqResp(paz_simulate['poles'], 

                paz_simulate['zeros'], paz_simulate['gain'], delta, nfft) 

 

    data[-1] = abs(data[-1]) + 0.0j 

    # transform data back into the time domain 

    data = np.fft.irfft(data)[0:ndat] 

    if pitsasim: 

        # linear detrend 

        data = simpleDetrend(data) 

    if shsim: 

        # detrend using least squares 

        data = scipy.signal.detrend(data, type="linear") 

    # correct for involved overall sensitivities 

    if paz_remove and remove_sensitivity and not seedresp: 

        data /= paz_remove['sensitivity'] 

    if paz_simulate and simulate_sensitivity: 

        data *= paz_simulate['sensitivity'] 

    return data 

 

 

def paz2AmpValueOfFreqResp(paz, freq): 

    """ 

    Returns Amplitude at one frequency for the given poles and zeros 

 

    :param paz: Given poles and zeros 

    :param freq: Given frequency 

 

    The amplitude of the freq is estimated according to "Of Poles and 

    Zeros", Frank Scherbaum, p 43. 

 

    .. rubric:: Example 

 

    >>> paz = {'poles': [-4.44 + 4.44j, -4.44 - 4.44j], 

    ...        'zeros': [0 + 0j, 0 + 0j], 

    ...        'gain': 0.4} 

    >>> amp = paz2AmpValueOfFreqResp(paz, 1) 

    >>> print(round(amp, 7)) 

    0.2830262 

    """ 

    jw = complex(0, 2 * np.pi * freq)  # angular frequency 

    fac = complex(1, 0) 

    for zero in paz['zeros']:  # numerator 

        fac *= jw - zero 

    for pole in paz['poles']:  # denominator 

        fac /= jw - pole 

    return abs(fac) * paz['gain'] 

 

 

def estimateMagnitude(paz, amplitude, timespan, h_dist): 

    """ 

    Estimates local magnitude from poles and zeros of given instrument, the 

    peak to peak amplitude and the time span from peak to peak. 

    Readings on two components can be used in magnitude estimation by providing 

    lists for ``paz``, ``amplitude`` and ``timespan``. 

 

    :param paz: PAZ of the instrument [m/s] or list of the same 

    :param amplitude: Peak to peak amplitude [counts] or list of the same 

    :param timespan: Timespan of peak to peak amplitude [s] or list of the same 

    :param h_dist: Hypocentral distance [km] 

    :returns: Estimated local magnitude Ml 

 

    .. note:: 

        Magnitude estimation according to Bakun & Joyner, 1984, Eq. (3) page 

        1835. Bakun, W. H. and W. B. Joyner: The Ml scale in central 

        California, Bull. Seismol. Soc. Am., 74, 1827-1843, 1984 

 

    .. rubric:: Example 

 

    >>> paz = {'poles': [-4.444+4.444j, -4.444-4.444j, -1.083+0j], 

    ...        'zeros': [0+0j, 0+0j, 0+0j], 

    ...        'gain': 1.0, 'sensitivity': 671140000.0} 

    >>> mag = estimateMagnitude(paz, 3.34e6, 0.065, 0.255) 

    >>> print(round(mag, 6)) 

    2.132873 

    >>> mag = estimateMagnitude([paz, paz], [3.34e6, 5e6], [0.065, 0.1], 0.255) 

    >>> print(round(mag, 6)) 

    2.347618 

    """ 

    # convert input to lists 

    if not isinstance(paz, list) and not isinstance(paz, tuple): 

        paz = [paz] 

    if not isinstance(amplitude, list) and not isinstance(amplitude, tuple): 

        amplitude = [amplitude] 

    if not isinstance(timespan, list) and not isinstance(timespan, tuple): 

        timespan = [timespan] 

    # convert every input amplitude to Wood Anderson and calculate the mean 

    wa_ampl_mean = 0.0 

    count = 0 

    for paz, amplitude, timespan in zip(paz, amplitude, timespan): 

        wa_ampl_mean += estimateWoodAndersonAmplitude(paz, amplitude, timespan) 

        count += 1 

    wa_ampl_mean /= count 

    # mean of input amplitudes (if more than one) should be used in final 

    # magnitude estimation (usually N and E components) 

    magnitude = np.log10(wa_ampl_mean) + np.log10(h_dist / 100.0) + \ 

                0.00301 * (h_dist - 100.0) + 3.0 

    return magnitude 

 

 

def estimateWoodAndersonAmplitude(paz, amplitude, timespan): 

    """ 

    Convert amplitude in counts measured of instrument with given Poles and 

    Zeros information for use in :func:`estimateMagnitude`. 

    Amplitude should be measured as full peak to peak amplitude, timespan as 

    difference of the two readings. 

 

    :param paz: PAZ of the instrument [m/s] or list of the same 

    :param amplitude: Peak to peak amplitude [counts] or list of the same 

    :param timespan: Timespan of peak to peak amplitude [s] or list of the same 

    :returns: Simulated zero to peak displacement amplitude on Wood Anderson 

        seismometer [mm] for use in local magnitude estimation. 

    """ 

    # analog to pitsa/plt/RCS/plt_wave.c,v, lines 4881-4891 

    freq = 1.0 / (2 * timespan) 

    wa_ampl = amplitude / 2.0  # half peak to peak amplitude 

    wa_ampl /= (paz2AmpValueOfFreqResp(paz, freq) * paz['sensitivity']) 

    wa_ampl *= paz2AmpValueOfFreqResp(WOODANDERSON, freq) * \ 

            WOODANDERSON['sensitivity'] 

    wa_ampl *= 1000  # convert to mm 

    return wa_ampl 

 

 

if __name__ == '__main__': 

    import doctest 

    doctest.testmod(exclude_empty=True)