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

# -*- coding: utf-8 -*- 

""" 

Signal processing functions for RtMemory objects. 

 

For sequential packet processing that requires memory (which includes recursive 

filtering), each processing function (e.g., :mod:`obspy.realtime.signal`) 

needs to manage the initialization and update of 

:class:`~obspy.realtime.rtmemory.RtMemory` object(s), and needs to know when 

and how to get values from this memory. 

 

For example: Boxcar smoothing: For each new data point available past the end 

of the boxcar, the original, un-smoothed data point value at the beginning of 

the boxcar has to be subtracted from the running boxcar sum, this value may be 

in a previous packet, so has to be retrieved from memory see 

:func:`obspy.realtime.signal.boxcar`. 

 

:copyright: 

    The ObsPy Development Team (devs@obspy.org) & Anthony Lomax 

:license: 

    GNU Lesser General Public License, Version 3 

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

""" 

 

import math 

import sys 

import numpy as np 

from obspy.core.trace import Trace, UTCDateTime 

from obspy.realtime.rtmemory import RtMemory 

 

_PI = math.pi 

_TWO_PI = 2.0 * math.pi 

_MIN_FLOAT_VAL = 1.0e-20 

 

 

def scale(trace, factor=1.0, rtmemory_list=None):  # @UnusedVariable 

    """ 

    Scale array data samples by specified factor. 

 

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

    :param trace:  :class:`~obspy.core.trace.Trace` object to append to this 

        RtTrace 

    :type factor: float, optional 

    :param factor: Scale factor (default is 1.0). 

    :type rtmemory_list: list of :class:`~obspy.realtime.rtmemory.RtMemory`, 

        optional 

    :param rtmemory_list: Persistent memory used by this process for specified 

        trace. 

    :rtype: NumPy :class:`numpy.ndarray` 

    :return: Processed trace data from appended Trace object. 

    """ 

    if not isinstance(trace, Trace): 

        msg = "trace parameter must be an obspy.core.trace.Trace object." 

        raise ValueError(msg) 

    trace.data *= factor 

    return trace.data 

 

 

def integrate(trace, rtmemory_list=None): 

    """ 

    Apply simple rectangular integration to array data. 

 

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

    :param trace:  :class:`~obspy.core.trace.Trace` object to append to this 

        RtTrace 

    :type rtmemory_list: list of :class:`~obspy.realtime.rtmemory.RtMemory`, 

        optional 

    :param rtmemory_list: Persistent memory used by this process for specified 

        trace. 

    :rtype: NumPy :class:`numpy.ndarray` 

    :return: Processed trace data from appended Trace object. 

    """ 

    if not isinstance(trace, Trace): 

        msg = "trace parameter must be an obspy.core.trace.Trace object." 

        raise ValueError(msg) 

 

    if not rtmemory_list: 

        rtmemory_list = [RtMemory()] 

 

    sample = trace.data 

    if np.size(sample) < 1: 

        return sample 

 

    delta_time = trace.stats.delta 

 

    rtmemory = rtmemory_list[0] 

 

    # initialize memory object 

    if not rtmemory.initialized: 

        memory_size_input = 0 

        memory_size_output = 1 

        rtmemory.initialize(sample.dtype, memory_size_input, 

                            memory_size_output, 0, 0) 

 

    sum = rtmemory.output[0] 

 

    for i in range(np.size(sample)): 

        sum += sample[i] * delta_time 

        sample[i] = sum 

 

    rtmemory.output[0] = sum 

 

    return sample 

 

 

def differentiate(trace, rtmemory_list=None): 

    """ 

    Apply simple differentiation to array data. 

 

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

    :param trace:  :class:`~obspy.core.trace.Trace` object to append to this 

        RtTrace 

    :type rtmemory_list: list of :class:`~obspy.realtime.rtmemory.RtMemory`, 

        optional 

    :param rtmemory_list: Persistent memory used by this process for specified 

        trace. 

    :rtype: NumPy :class:`numpy.ndarray` 

    :return: Processed trace data from appended Trace object. 

    """ 

    if not isinstance(trace, Trace): 

        msg = "trace parameter must be an obspy.core.trace.Trace object." 

        raise ValueError(msg) 

 

    if not rtmemory_list: 

        rtmemory_list = [RtMemory()] 

 

    sample = trace.data 

    if np.size(sample) < 1: 

        return(sample) 

 

    delta_time = trace.stats.delta 

 

    rtmemory = rtmemory_list[0] 

 

    # initialize memory object 

    if not rtmemory.initialized: 

        memory_size_input = 1 

        memory_size_output = 0 

        rtmemory.initialize(sample.dtype, memory_size_input, 

                            memory_size_output, 0, 0) 

        # avoid large diff value for first output sample 

        rtmemory.input[0] = sample[0] 

 

    previous_sample = rtmemory.input[0] 

 

    for i in range(np.size(sample)): 

        diff = (sample[i] - previous_sample) / delta_time 

        previous_sample = sample[i] 

        sample[i] = diff 

 

    rtmemory.input[0] = previous_sample 

 

    return sample 

 

 

def boxcar(trace, width, rtmemory_list=None): 

    """ 

    Apply boxcar smoothing to data in array sample. 

 

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

    :param trace:  :class:`~obspy.core.trace.Trace` object to append to this 

        RtTrace 

    :type width: int 

    :param width: Width in number of sample points for filter. 

    :type rtmemory_list: list of :class:`~obspy.realtime.rtmemory.RtMemory`, 

        optional 

    :param rtmemory_list: Persistent memory used by this process for specified 

        trace. 

    :rtype: NumPy :class:`numpy.ndarray` 

    :return: Processed trace data from appended Trace object. 

    """ 

    if not isinstance(trace, Trace): 

        msg = "trace parameter must be an obspy.core.trace.Trace object." 

        raise ValueError(msg) 

 

    if not width > 0: 

        msg = "width parameter not specified or < 1." 

        raise ValueError(msg) 

 

    if not rtmemory_list: 

        rtmemory_list = [RtMemory()] 

 

    sample = trace.data 

 

    rtmemory = rtmemory_list[0] 

 

    # initialize memory object 

    if not rtmemory.initialized: 

        memory_size_input = width 

        memory_size_output = 0 

        rtmemory.initialize(sample.dtype, memory_size_input, 

                            memory_size_output, 0, 0) 

 

    # initialize array for time-series results 

    new_sample = np.zeros(np.size(sample), sample.dtype) 

 

    i = 0 

    i1 = i - width 

    i2 = i      # causal boxcar of width width 

    sum = 0.0 

    icount = 0 

    for i in range(np.size(sample)): 

        value = 0.0 

        if (icount == 0):    # first pass, accumulate sum 

            for n in range(i1, i2 + 1): 

                if (n < 0): 

                    value = rtmemory.input[width + n] 

                else: 

                    value = sample[n] 

                sum += value 

                icount = icount + 1 

        else:                # later passes, update sum 

            if ((i1 - 1) < 0): 

                value = rtmemory.input[width + (i1 - 1)] 

            else: 

                value = sample[(i1 - 1)] 

            sum -= value 

            if (i2 < 0): 

                value = rtmemory.input[width + i2] 

            else: 

                value = sample[i2] 

            sum += value 

        if (icount > 0): 

            new_sample[i] = (float)(sum / float(icount)) 

        else: 

            new_sample[i] = 0.0 

        i1 = i1 + 1 

        i2 = i2 + 1 

 

    rtmemory.updateInput(sample) 

 

    return new_sample 

 

 

def tauc(trace, width, rtmemory_list=None): 

    """ 

    Calculate instantaneous period in a fixed window (Tau_c). 

 

    .. seealso:: 

 

        Implements equations 1-3 in [Allen2003]_ except use a fixed width 

        window instead of decay function. 

 

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

    :param trace:  :class:`~obspy.core.trace.Trace` object to append to this 

        RtTrace 

    :type width: int 

    :param width: Width in number of sample points for tauc window. 

    :type rtmemory_list: list of :class:`~obspy.realtime.rtmemory.RtMemory`, 

        optional 

    :param rtmemory_list: Persistent memory used by this process for specified 

        trace. 

    :rtype: NumPy :class:`numpy.ndarray` 

    :return: Processed trace data from appended Trace object. 

    """ 

    if not isinstance(trace, Trace): 

        msg = "trace parameter must be an obspy.core.trace.Trace object." 

        raise ValueError(msg) 

 

    if not width > 0: 

        msg = "tauc: width parameter not specified or < 1." 

        raise ValueError(msg) 

 

    if not rtmemory_list: 

        rtmemory_list = [RtMemory(), RtMemory()] 

 

    sample = trace.data 

    delta_time = trace.stats.delta 

 

    rtmemory = rtmemory_list[0] 

    rtmemory_dval = rtmemory_list[1] 

 

    sample_last = 0.0 

 

    # initialize memory object 

    if not rtmemory.initialized: 

        memory_size_input = width 

        memory_size_output = 1 

        rtmemory.initialize(sample.dtype, memory_size_input, 

                            memory_size_output, 0, 0) 

        sample_last = sample[0] 

    else: 

        sample_last = rtmemory.input[width - 1] 

 

    # initialize memory object 

    if not rtmemory_dval.initialized: 

        memory_size_input = width 

        memory_size_output = 1 

        rtmemory_dval.initialize(sample.dtype, memory_size_input, 

                                 memory_size_output, 0, 0) 

 

    new_sample = np.zeros(np.size(sample), sample.dtype) 

    deriv = np.zeros(np.size(sample), sample.dtype) 

 

    #sample_last = rtmemory.input[width - 1] 

    sample_d = 0.0 

    deriv_d = 0.0 

    xval = rtmemory.output[0] 

    dval = rtmemory_dval.output[0] 

 

    for i in range(np.size(sample)): 

 

        sample_d = sample[i] 

        deriv_d = (sample_d - sample_last) / delta_time 

        indexBegin = i - width 

        if (indexBegin >= 0): 

            xval = xval - (sample[indexBegin]) * (sample[indexBegin]) \ 

                + sample_d * sample_d 

            dval = dval - deriv[indexBegin] * deriv[indexBegin] \ 

                + deriv_d * deriv_d 

        else: 

            index = i 

            xval = xval - rtmemory.input[index] * rtmemory.input[index] \ 

                + sample_d * sample_d 

            dval = dval \ 

                - rtmemory_dval.input[index] * rtmemory_dval.input[index] \ 

                + deriv_d * deriv_d 

        deriv[i] = deriv_d 

        sample_last = sample_d 

        # if (xval > _MIN_FLOAT_VAL &  & dval > _MIN_FLOAT_VAL) { 

        if (dval > _MIN_FLOAT_VAL): 

            new_sample[i] = _TWO_PI * math.sqrt(xval / dval) 

        else: 

            new_sample[i] = 0.0 

 

    # update memory 

    rtmemory.output[0] = xval 

    rtmemory.updateInput(sample) 

    rtmemory_dval.output[0] = dval 

    rtmemory_dval.updateInput(deriv) 

 

    return new_sample 

 

# memory object indices for storing specific values 

_AMP_AT_PICK = 0 

_HAVE_USED_MEMORY = 1 

_FLAG_COMPETE_MWP = 2 

_INT_INT_SUM = 3 

_POLARITY = 4 

_MEMORY_SIZE_OUTPUT = 5 

 

 

def mwpIntegral(trace, max_time, ref_time, mem_time=1.0, gain=1.0, 

                rtmemory_list=None): 

    """ 

    Calculate Mwp integral on a displacement trace. 

 

    .. seealso:: [Tsuboi1999]_ and [Tsuboi1995]_ 

 

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

    :param trace:  :class:`~obspy.core.trace.Trace` object to append to this 

        RtTrace 

    :type max_time: float 

    :param max_time: Maximum time in seconds after ref_time to apply Mwp 

        integration. 

    :type ref_time: :class:`~obspy.core.utcdatetime.UTCDateTime` 

    :param ref_time: Reference date and time of the data sample 

        (e.g. P pick time) at which to begin Mwp integration. 

    :type mem_time: float, optional 

    :param mem_time: Length in seconds of data memory (must be much larger 

        than maximum delay between pick declaration and pick time). Defaults 

        to ``1.0``. 

    :type gain: float, optional 

    :param gain: Nominal gain to convert input displacement trace to meters 

        of ground displacement. Defaults to ``1.0``. 

    :type rtmemory_list: list of :class:`~obspy.realtime.rtmemory.RtMemory`, 

        optional 

    :param rtmemory_list: Persistent memory used by this process for specified 

        trace. 

    :rtype: NumPy :class:`numpy.ndarray` 

    :return: Processed trace data from appended Trace object. 

    """ 

    if not isinstance(trace, Trace): 

        msg = "trace parameter must be an obspy.core.trace.Trace object." 

        raise ValueError(msg) 

 

    if not isinstance(ref_time, UTCDateTime): 

        msg = "ref_time must be an obspy.core.utcdatetime.UTCDateTime object." 

        raise ValueError(msg) 

 

    if not max_time >= 0: 

        msg = "max_time parameter not specified or < 0." 

        raise ValueError(msg) 

 

    if not rtmemory_list: 

        rtmemory_list = [RtMemory()] 

 

    sample = trace.data 

    delta_time = trace.stats.delta 

 

    rtmemory = rtmemory_list[0] 

 

    # initialize memory object 

    if not rtmemory.initialized: 

        memory_size_input = int(0.5 + mem_time * trace.stats.sampling_rate) 

        memory_size_output = _MEMORY_SIZE_OUTPUT 

        rtmemory.initialize(sample.dtype, memory_size_input, 

                            memory_size_output, 0, 0) 

 

    new_sample = np.zeros(np.size(sample), sample.dtype) 

 

    ioffset_pick = int(round( 

                       (ref_time - trace.stats.starttime) 

                       * trace.stats.sampling_rate)) 

    ioffset_mwp_min = ioffset_pick 

 

    # set reference amplitude 

    if ioffset_mwp_min >= 0 and ioffset_mwp_min < trace.data.size: 

        # value in trace data array 

        rtmemory.output[_AMP_AT_PICK] = trace.data[ioffset_mwp_min] 

    elif ioffset_mwp_min >= -(np.size(rtmemory.input)) and ioffset_mwp_min < 0: 

        # value in memory array 

        index = ioffset_mwp_min + np.size(rtmemory.input) 

        rtmemory.output[_AMP_AT_PICK] = rtmemory.input[index] 

    elif ioffset_mwp_min < -(np.size(rtmemory.input)) \ 

        and not rtmemory.output[_HAVE_USED_MEMORY]: 

        msg = "mem_time not large enough to buffer required input data." 

        raise ValueError(msg) 

    if ioffset_mwp_min < 0 and rtmemory.output[_HAVE_USED_MEMORY]: 

        ioffset_mwp_min = 0 

    else: 

        rtmemory.output[_HAVE_USED_MEMORY] = 1 

    # set Mwp end index corresponding to maximum duration 

    mwp_end_index = int(round(max_time / delta_time)) 

    ioffset_mwp_max = mwp_end_index + ioffset_pick 

    if ioffset_mwp_max < trace.data.size: 

        rtmemory.output[_FLAG_COMPETE_MWP] = 1  # will complete 

    if ioffset_mwp_max > trace.data.size: 

        ioffset_mwp_max = trace.data.size 

    # apply double integration, check for extrema 

    mwp_amp_at_pick = rtmemory.output[_AMP_AT_PICK] 

    mwp_int_int_sum = rtmemory.output[_INT_INT_SUM] 

    polarity = rtmemory.output[_POLARITY] 

    amplitude = 0.0 

    for n in range(ioffset_mwp_min, ioffset_mwp_max): 

        if n >= 0: 

            amplitude = trace.data[n] 

        elif n >= -(np.size(rtmemory.input)): 

            # value in memory array 

            index = n + np.size(rtmemory.input) 

            amplitude = rtmemory.input[index] 

        else: 

            msg = "Error: Mwp: attempt to access rtmemory.input array of " + \ 

                "size=%d at invalid index=%d: this should not happen!" % \ 

                (np.size(rtmemory.input), n + np.size(rtmemory.input)) 

            print msg 

            continue  # should never reach here 

        disp_amp = amplitude - mwp_amp_at_pick 

        # check displacement polarity 

        if disp_amp >= 0.0:  # pos 

            # check if past extremum 

            if  polarity < 0:  # passed from neg to pos displacement 

                mwp_int_int_sum *= -1.0 

                mwp_int_int_sum = 0 

            polarity = 1 

        elif disp_amp < 0.0:  # neg 

            # check if past extremum 

            if polarity > 0:  # passed from pos to neg displacement 

                mwp_int_int_sum = 0 

            polarity = -1 

        mwp_int_int_sum += (amplitude - mwp_amp_at_pick) * delta_time / gain 

        new_sample[n] = mwp_int_int_sum 

 

    rtmemory.output[_INT_INT_SUM] = mwp_int_int_sum 

    rtmemory.output[_POLARITY] = polarity 

 

    # update memory 

    rtmemory.updateInput(sample) 

 

    return new_sample 

 

 

MWP_INVALID = -9.9 

# 4.213e19 - Tsuboi 1995, 1999 

MWP_CONST = 4.0 * _PI  # 4 PI 

MWP_CONST *= 3400.0  # rho 

MWP_CONST *= 7900.0 * 7900.0 * 7900.0  # Pvel**3 

MWP_CONST *= 2.0  # FP average radiation pattern 

MWP_CONST *= (10000.0 / 90.0)  # distance deg -> km 

MWP_CONST *= 1000.0  # distance km -> meters 

# http://mail.python.org/pipermail/python-list/2010-February/1235196.html, ff. 

try: 

    FLOAT_MIN = sys.float_info.min 

except AttributeError: 

    FLOAT_MIN = 1.1e-37 

 

 

def calculateMwpMag(peak, epicentral_distance): 

    """ 

    Calculate Mwp magnitude. 

 

    .. seealso:: [Tsuboi1999]_ and [Tsuboi1995]_ 

 

    :type peak: float 

    :param peak: Peak value of integral of displacement seismogram. 

    :type epicentral_distance: float 

    :param epicentral_distance: Great-circle epicentral distance from station 

        in degrees. 

    :rtype: float 

    :returns: Calculated Mwp magnitude. 

    """ 

    moment = MWP_CONST * peak * epicentral_distance 

    mwp_mag = MWP_INVALID 

    if moment > FLOAT_MIN: 

        mwp_mag = (2.0 / 3.0) * (math.log10(moment) - 9.1) 

    return mwp_mag