def mean(trace, win=1.0, rtmemory_list=None): """ Calculate recursive mean. win is a window length """ if not isinstance(trace, Trace): msg = "Trace parameter must be an obspy.core.trace.Trace object." raise ValueError(msg) # if this is the first appended trace, the rtmemory_list will be None if not rtmemory_list: rtmemory_list = [RtMemory()] # deal with case of empty trace sample = trace.data dt = trace.stats.delta if np.size(sample) < 1: return sample # get simple info from trace npts = len(sample) # prepare the output array mu1 = np.empty(npts, sample.dtype) # prepare the rt memory rtmemory_mu1 = rtmemory_list[0] if not rtmemory_mu1.initialized: memory_size_input = 1 memory_size_output = 0 rtmemory_mu1.initialize(sample.dtype, memory_size_input,\ memory_size_output, 0, 0) # initialize from memory mu1_last = rtmemory_mu1.input[0] C1 = dt / float(win) a1 = 1 - C1 # do recursive mean for i in xrange(npts): mu1[i] = a1 * mu1_last + C1 * sample[i] mu1_last = mu1[i] # save to memory rtmemory_mu1.input[0] = mu1_last return mu1
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 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 kurtosis(trace, win=3.0, rtmemory_list=None): """ Apply recursive kurtosis calculation on data. Recursive kurtosis is computed using the [ChassandeMottin2002]_ formulation adjusted to give the kurtosis of a gaussian distribution = 0.0. :type trace: :class:`~obspy.core.trace.Trace` :param trace: :class:`~obspy.core.trace.Trace` object to append to this RtTrace :type win: float, optional :param win: window length in seconds for the kurtosis (default is 3.0 s) :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 this is the first appended trace, the rtmemory_list will be None if not rtmemory_list: rtmemory_list = [RtMemory(), RtMemory(), RtMemory()] # deal with case of empty trace sample = trace.data if np.size(sample) < 1: return sample # get simple info from trace npts = len(sample) dt = trace.stats.delta # set some constants for the kurtosis calulation C1 = dt / float(win) a1 = 1.0 - C1 C2 = (1.0 - a1 * a1) / 2.0 bias = -3 * C1 - 3.0 # prepare the output array kappa4 = np.empty(npts, sample.dtype) # initialize the real-time memory needed to store # the recursive kurtosis coefficients until the # next bloc of data is added rtmemory_mu1 = rtmemory_list[0] rtmemory_mu2 = rtmemory_list[1] rtmemory_k4_bar = rtmemory_list[2] # there are three memory objects, one for each "last" coefficient # that needs carrying over # initialize mu1_last to 0 if not rtmemory_mu1.initialized: memory_size_input = 1 memory_size_output = 0 rtmemory_mu1.initialize(sample.dtype, memory_size_input, memory_size_output, 0, 0) # initialize mu2_last (sigma) to 1 if not rtmemory_mu2.initialized: memory_size_input = 1 memory_size_output = 0 rtmemory_mu2.initialize(sample.dtype, memory_size_input, memory_size_output, 1, 0) # initialize k4_bar_last to 0 if not rtmemory_k4_bar.initialized: memory_size_input = 1 memory_size_output = 0 rtmemory_k4_bar.initialize(sample.dtype, memory_size_input, memory_size_output, 0, 0) mu1_last = rtmemory_mu1.input[0] mu2_last = rtmemory_mu2.input[0] k4_bar_last = rtmemory_k4_bar.input[0] # do recursive kurtosis for i in xrange(npts): mu1 = a1 * mu1_last + C1 * sample[i] dx2 = (sample[i] - mu1_last) * (sample[i] - mu1_last) mu2 = a1 * mu2_last + C2 * dx2 dx2 = dx2 / mu2_last k4_bar = (1 + C1 - 2 * C1 * dx2) * k4_bar_last + C1 * dx2 * dx2 kappa4[i] = k4_bar + bias mu1_last = mu1 mu2_last = mu2 k4_bar_last = k4_bar rtmemory_mu1.input[0] = mu1_last rtmemory_mu2.input[0] = mu2_last rtmemory_k4_bar.input[0] = k4_bar_last return kappa4
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
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
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 registerRtProcess(self, process, **options): """ Adds real-time processing algorithm to processing list of this RtTrace. Processing function must be one of: %s. % REALTIME_PROCESS_FUNCTIONS.keys() or a non-recursive, time-domain np or obspy function which takes a single array as an argument and returns an array :type process: str or function :param process: Specifies which processing function is added, e.g. ``"boxcar"`` or ``np.abs``` (functions without brackets). See :mod:`obspy.realtime.signal` for all predefined processing functions. :type options: dict, optional :param options: Required keyword arguments to be passed the respective processing function, e.g. ``width=100`` for ``'boxcar'`` process. See :mod:`obspy.realtime.signal` for all options. :rtype: int :return: Length of processing list after registering new processing function. """ # create process_name either from string or function name process_name = ("%s" % process).lower() # set processing entry for this process entry = False rtmemory_list = None if hasattr(process, '__call__'): # direct function call entry = (process, options, None) elif process_name in REALTIME_PROCESS_FUNCTIONS: # predefined function num = REALTIME_PROCESS_FUNCTIONS[process_name][1] if num: # make sure we have num new RtMemory instances rtmemory_list = [RtMemory() for _i in range(num)] entry = (process_name, options, rtmemory_list) else: # check if process name is contained within a predefined function, # e.g. 'int' for 'integrate' for key in REALTIME_PROCESS_FUNCTIONS: if not key.startswith(process_name): continue process_name = key num = REALTIME_PROCESS_FUNCTIONS[process_name][1] if num: # make sure we have num new RtMemory instances rtmemory_list = [RtMemory() for _i in range(num)] entry = (process_name, options, rtmemory_list) break if not entry: raise NotImplementedError("Can't register process %s" % (process)) # add process entry self.processing.append(entry) # add processing information to the stats dictionary proc_info = "realtime_process:%s:%s" % (process_name, options) self._addProcessingInfo(proc_info) return len(self.processing)
def append(self, trace, gap_overlap_check=False, verbose=False): """ Appends a Trace object to this RtTrace. Registered real-time processing will be applied to copy of appended Trace object before it is appended. This RtTrace will be truncated from the beginning to RtTrace.max_length, if specified. Sampling rate, data type and trace.id of both traces must match. :type trace: :class:`~obspy.core.trace.Trace` :param trace: :class:`~obspy.core.trace.Trace` object to append to this RtTrace :type gap_overlap_check: bool, optional :param gap_overlap_check: Action to take when there is a gap or overlap between the end of this RtTrace and start of appended Trace: If True, raise TypeError. If False, all trace processing memory will be re-initialized to prevent false signal in processed trace. (default is ``True``). :type verbose: bool, optional :param verbose: Print additional information to stdout :return: NumPy :class:`np.ndarray` object containing processed trace data from appended Trace object. """ if not isinstance(trace, Trace): # only add Trace objects raise TypeError("Only obspy.core.trace.Trace objects are allowed") # sanity checks if self.have_appended_data: # check id if self.getId() != trace.getId(): raise TypeError("Trace ID differs:", self.getId(), trace.getId()) # check sample rate if self.stats.sampling_rate != trace.stats.sampling_rate: raise TypeError("Sampling rate differs:", self.stats.sampling_rate, trace.stats.sampling_rate) # check calibration factor if self.stats.calib != trace.stats.calib: raise TypeError("Calibration factor differs:", self.stats.calib, trace.stats.calib) # check data type if self.data.dtype != trace.data.dtype: raise TypeError("Data type differs:", self.data.dtype, trace.data.dtype) # TODO: IMPORTANT? Should improve check for gaps and overlaps # and handle more elegantly # check times gap_or_overlap = False if self.have_appended_data: #delta = int(math.floor(\ # round((rt.stats.starttime - lt.stats.endtime) * sr, 5) )) - 1 diff = trace.stats.starttime - self.stats.endtime delta = diff * self.stats.sampling_rate - 1.0 if verbose: msg = "%s: Overlap/gap of (%g) samples in data: (%s) (%s) " + \ "diff=%gs dt=%gs" print msg % (self.__class__.__name__, delta, self.stats.endtime, trace.stats.starttime, diff, self.stats.delta) if delta < -0.1: msg = "Overlap of (%g) samples in data: (%s) (%s) diff=%gs" + \ " dt=%gs" msg = msg % (-delta, self.stats.endtime, trace.stats.starttime, diff, self.stats.delta) if gap_overlap_check: raise TypeError(msg) gap_or_overlap = True if delta > 0.1: msg = "Gap of (%g) samples in data: (%s) (%s) diff=%gs" + \ " dt=%gs" msg = msg % (delta, self.stats.endtime, trace.stats.starttime, diff, self.stats.delta) if gap_overlap_check: raise TypeError(msg) gap_or_overlap = True if gap_or_overlap: msg += " - Trace processing memory will be re-initialized." warnings.warn(msg, UserWarning) else: # correct start time to pin absolute trace timing to start of # appended trace, this prevents slow drift of nominal trace # timing from absolute time when nominal sample rate differs # from true sample rate self.stats.starttime = \ self.stats.starttime + diff - self.stats.delta if verbose: print "%s: self.stats.starttime adjusted by: %gs" \ % (self.__class__.__name__, diff - self.stats.delta) # first apply all registered processing to Trace for proc in self.processing: process_name, options, rtmemory_list = proc # if gap or overlap, clear memory if gap_or_overlap and rtmemory_list != None: for n in range(len(rtmemory_list)): rtmemory_list[n] = RtMemory() # apply processing trace = trace.copy() dtype = trace.data.dtype if hasattr(process_name, '__call__'): # check if direct function call trace.data = process_name(trace.data, **options) else: # got predefined function func = REALTIME_PROCESS_FUNCTIONS[process_name.lower()][0] options['rtmemory_list'] = rtmemory_list trace.data = func(trace, **options) # assure dtype is not changed trace.data = np.require(trace.data, dtype=dtype) # if first data, set stats if not self.have_appended_data: self.data = np.array(trace.data) self.stats = Stats(header=trace.stats) self.have_appended_data = True return trace # handle all following data sets # fix Trace.__add__ parameters # TODO: IMPORTANT? Should check for gaps and overlaps and handle # more elegantly sum_trace = Trace.__add__(self, trace, method=0, interpolation_samples=0, fill_value='latest', sanity_checks=True) # Trace.__add__ returns new Trace, so update to this RtTrace self.data = sum_trace.data # left trim if data length exceeds max_length if self.max_length != None: max_samples = int(self.max_length * self.stats.sampling_rate + 0.5) if np.size(self.data) > max_samples: starttime = self.stats.starttime + (np.size(self.data) - \ max_samples) / self.stats.sampling_rate self._ltrim(starttime, pad=False, nearest_sample=True, fill_value=None) return trace
def sw_kurtosis(trace, win=3.0, rtmemory_list=None): """ Compute kurtosis using a sliding window method. Calls scipy.stats.kurtosis :type trace: :class:`~obspy.core.trace.Trace` :param trace: :class:`~obspy.core.trace.Trace` object to append to this RtTrace :type win: float :param win: sliding window length in seconds :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()] # deal with case of empty trace sample = trace.data if np.size(sample) < 1: return sample # get info from trace dt = trace.stats.delta npts = int(np.round(win / float(dt))) # set up temporary array for kurtosis and output array k_array = np.empty((len(sample), npts), dtype=sample.dtype) x_out = np.empty(len(sample)) rtmemory = rtmemory_list[0] mem_size = npts - 1 if not rtmemory.initialized: sample_start = sample[0:mem_size] first_trace = True memory_size_input = mem_size memory_size_output = 0 rtmemory.initialize(sample.dtype, memory_size_input, memory_size_output) rtmemory.input = sample_start[::-1] # make an array of the right dimension x = np.empty(len(sample) + mem_size, dtype=sample.dtype) # fill it up partly with the memory, partly with the new data x[0:mem_size] = rtmemory.input[:] x[mem_size:] = sample[:] # do the sliding window kurtosis windows = _window(x, npts) i = 0 for w in windows: k_array[i, 0:npts] = w i = i + 1 xout = ss.kurtosis(k_array, axis=1) # put new data into memory for next trace rtmemory.updateInput(sample) return xout
def convolve(trace, conv_signal=None, rtmemory_list=None): """ Convolve data with a (complex) signal (conv_signal). Note that the signal length should be odd. For a signal of (2N+1) points, the resulting output trace will be time shifted by -N*dt where dt is the sampling rate of the trace. You may want to consider shifting the trace starttime by the same amount before appending to this RtTrace :: tshift = N*dt tr.stats.starttime -= tshift someRtTrace.append(tr) :type trace: :class:`~obspy.core.trace.Trace` :param trace: :class:`~obspy.core.trace.Trace` object to append to this RtTrace :type conv_signal: :class:`numpy.ndarray`, optional :param conv_signal: signal with which to perform convolution :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 conv_signal == None: return trace if not rtmemory_list: rtmemory_list = [RtMemory()] # deal with case of empty trace sample = trace.data if np.size(sample) < 1: return sample flen = len(conv_signal) flen2 = (flen - 1) / 2 mem_size = 3 * flen2 + 1 rtmemory = rtmemory_list[0] if not rtmemory.initialized: first_trace = True memory_size_input = mem_size memory_size_output = 0 rtmemory.initialize(sample.dtype, memory_size_input,\ memory_size_output, 0, 0) # make an array of the right dimension x = np.empty(len(sample) + mem_size, dtype=complex) # fill it up partly with the memory, partly with the new data x[0:mem_size] = rtmemory.input[:] x[mem_size:] = sample[:] # do the convolution x_new = np.real(np.convolve(x, conv_signal, 'same')) i_start = mem_size - flen2 i_end = i_start + len(sample) # put new data into memory for next trace rtmemory.updateInput(sample) return x_new[i_start:i_end]
def dx2(trace, win=1.0, rtmemory_list=None): """ Calculate recursive variance. C is a scaling constant """ if not isinstance(trace, Trace): msg = "Trace parameter must be an obspy.core.trace.Trace object." raise ValueError(msg) # if this is the first appended trace, the rtmemory_list will be None if not rtmemory_list: rtmemory_list = [RtMemory(), RtMemory()] # deal with case of empty trace # are going to need double precision here sample = trace.data.astype('float64') if np.size(sample) < 1: return sample # get simple info from trace npts = len(sample) dt = trace.stats.delta # prepare the output array dx2 = np.empty(npts, sample.dtype) # prepare the rt memory rtmemory_mu1 = rtmemory_list[0] rtmemory_mu2 = rtmemory_list[1] if not rtmemory_mu1.initialized: memory_size_input = 1 memory_size_output = 0 rtmemory_mu1.initialize(sample.dtype, memory_size_input,\ memory_size_output, 0, 0) if not rtmemory_mu2.initialized: memory_size_input = 1 memory_size_output = 0 rtmemory_mu2.initialize(sample.dtype, memory_size_input,\ memory_size_output, sample[0]*sample[0], 0) # initialize from memory mu1_last = rtmemory_mu1.input[0] mu2_last = rtmemory_mu2.input[0] C1 = dt / float(win) a1 = 1.0 - C1 C2 = (1.0 - a1 * a1) / 2.0 # do recursive mean for i in xrange(npts): mu1 = a1 * mu1_last + C1 * sample[i] mu1_last = mu1 dx2[i] = (sample[i] - mu1_last) * (sample[i] - mu1_last) mu2 = a1 * mu2_last + C2 * dx2[i] dx2[i] = dx2[i] / mu2_last mu2_last = mu2 # save to memory rtmemory_mu1.input[0] = mu1_last rtmemory_mu2.input[0] = mu2_last return dx2