Package llnl.gnem.core.util
Class SeriesMath
java.lang.Object
llnl.gnem.core.util.SeriesMath
A class that provides a number of static methods useful for operating on
series data.
- Author:
- Doug Dodge
-
Nested Class Summary
-
Constructor Summary
-
Method Summary
Modifier and TypeMethodDescriptionstatic double[]
abs
(double[] data) convert all the elements of the data series to their absolute valuestatic float[]
abs
(float[] data) convert all the elements of the data series to their absolute valuestatic double[]
add
(double[] data, double value) Add a scalar value to all elements of the array and return a new arraystatic double[]
add
(double[] data1, double[] data2) Element by element addition of two data series of equal lengthstatic float[]
add
(float[] data, double value) static float[]
add
(float[] data, float value) static void
addScalar
(double[] data, double value) Add a scalar to the series.static void
addScalar
(float[] data, double value) Add a scalar to the series.static float[]
arrayMultiply
(float[] array1, float[] array2) static double[]
autocorrelate
(double[] x) Calculate the autocorrelation of a data seriesstatic double
autocorrelate
(double[] x, int shift) Calculate the autocorrelation at a specific offsetstatic double[]
autocorrelate
(float[] x) Calculate the autocorrelation of a data seriesstatic double
autocorrelate
(float[] x, int shift) Calculate the autocorrelation at a specific offsetstatic int
countGlitches
(float[] data, int windowLength, double threshold) double
covariance
(double[] vec1, double[] vec2) Calculates the covariance of two equally sized vectors.static double[]
crosscorrelate
(double[] x, double[] y) Get the cross-correlation series of two data sets for shift = 0:size result[shift] = sum( x(t) * y(t+shift))static double
crosscorrelate
(double[] x, double[] y, int shift) Get the value of the cross-correlation series at a specific offset result[shift] = sum( x(t) * y(t+shift) )static double[]
crosscorrelate
(float[] x, float[] y) Get the cross-correlation series of two data sets for shift = 0:size result[shift] = sum( x(t) * y(t+shift))static double
crosscorrelate
(float[] x, float[] y, int shift) Get the value of the cross-correlation series at a specific offset result[shift] = sum( x(t) * y(t+shift) )static double[]
cut
(double[] data, double minVal, double maxVal) Produce a cut version of the input series.static double[]
decimate
(double[] data, int decimationfactor) A double value version of the decimation methodstatic float[]
decimate
(float[] data, int decimationfactor) From the SeriesMath method given a data series it is decimated so that only every Nth point is retained This should be equivalent to the more versatile interpolate methods where N is the decimationfactorstatic void
differentiate
(float[] data, double samprate) static Double
dotProduct
(double[] data1, double[] data2) dot product of two equal length Series dot(A,B) = SUM(AiBi) = A1B1 + A2B2+...AnBnstatic Double
dotProduct
(float[] data1, float[] data2) dot product of two equal length Series dot(A,B) = SUM(AiBi) = A1B1 + A2B2+...AnBnstatic float[]
doubleToFloat
(double[] data) static float[]
elementAddition
(float[] data1, float[] data2) adds data1's values to data2'sstatic double[]
elementDivision
(double[] data1, double[] data2) divide the individual elements of one data series by those of anotherstatic float[]
elementDivision
(float[] data1, float[] data2) divide the individual elements of one data series by those of anotherstatic double[]
elementMultiplication
(double[] data1, double[] data2) Element-by-element multiplication, C = A.*Bstatic float[]
elementMultiplication
(float[] data1, float[] data2) Element-by-element multiplication, C = A.*Bstatic float[]
elementSubtraction
(float[] data1, float[] data2) Subtracts from data1 data2's valuesstatic float[]
envelope
(float[] data) Calculate the envelope of real valued data using a Hilbert transform E(t) = sqrt( data^2 + Hilbert(data) ^2 )static int[]
Finds the index of values meeting specified conditions.static int[]
For 2 arrays of equal length find the elements of array 1 that are ("GT", "LT", "==" ...) the elements of array 2static int[]
Finds the index of values meeting specified conditions.static int[]
Finds the index of values meeting specified conditions.static int[]
Finds the index of values meeting specified conditions.static DiscontinuityCollection
findDiscontinuities
(float[] data, double sampRate, int winLength, double factor) findIndicesOf
(double[] data, double value) Find the indices in the data series that match a specific valuestatic double[]
floatToDouble
(float[] data) static double
getDefiniteIntegral
(float[] data, double delta, int idx1, int idx2) Compute definite integral using extended trapezoidal rule.static double
getDefiniteIntegral
(SeriesMath.Function f, double delta, int idx1, int idx2) static double
getEnergy
(float[] x) static float
getExtremum
(float[] data) static float[]
getFirstDifference
(float[] data) static Collection<Glitch>
getGlitches
(float[] data, int windowLength, double threshold) static double
getKurtosis
(float[] data) static double
getMax
(double[] data) Gets the maximum value of the series.static float
getMax
(float[] data) Gets the maximum value of the series.static float
getMax
(float[] data, int idx1, int idx2) static double
getMax
(Collection<Double> data) getMaxIndex
(double[] data) Gets the maximum value of the series and the index it occurs at.getMaxIndex
(float[] data) Gets the maximum value of the series and the index it occurs at.static double
getMean
(double[] data) Gets the mean value of the series.static double
getMean
(float[] data) Gets the mean value of the series.static double
getMean
(Collection<? extends Number> c) static double
getMedian
(double[] values) static double
getMedian
(float[] data) Gets the median value of the time series.static double
static double
getMin
(double[] data) Gets the minimum value of the series.static float
getMin
(float[] data) Gets the minimum value of the series.static float
getMin
(float[] data, int idx1, int idx2) static double
getMin
(Collection<Double> data) getMinIndex
(double[] data) Gets the maximum value of the series and the index it occurs at.getMinIndex
(float[] data) Gets the minimum value of the series and the index it occurs at.static double
getNorm
(double[] data) Treat the series as a vector and calculate the norm (aka the direction or the normal) of the vectorstatic double
getPeakToPeakAmplitude
(float[] data) static double
getPeakToPeakAmplitude
(float[] data, double sampleInterval, double period) Gets the peakToPeakAmplitude attribute of an input timeseries at the specified period.static double
getRange
(float[] u) static double
getRMS
(double[] data) Gets the RMS value of the data in the input float array.static double
getRMS
(float[] data) Gets the RMS value of the data in the input float array.static double
getRMS
(Collection<? extends Number> c) static double
getSnr
(float[] data, double sampleRate, double dataStartTime, double pickEpochTime, double prePickWindowLength, double postPickWindowLength) static double
getSnr
(float[] data, int pickIndex, int preWinStartIndex, int postWinEndIndex) static double
getStandardDeviation
(Collection<? extends Number> c) static double
getStd
(double[] vec, int idx1, int idx2) static double
getStDev
(double[] data) Gets the standard deviation of the time series.static double
getStDev
(float[] data) Gets the standard deviation of the time series.static double[]
getSubSection
(double[] v, int minIdx, int maxIdx) static float[]
getSubSection
(float[] v, int minIdx, int maxIdx) static int[]
getSubSection
(int[] v, int minIdx, int maxIdx) static double
getSum
(double[] data) Gets the sum of the values in the input arraystatic double
getSum
(float[] data) Gets the sum of the values in the input arraystatic double
getSum
(Collection<? extends Number> c) static double
getSumOfSquares
(float[] data) static double[]
getUnitVector
(double[] data) Treat the series as a vector and calculate the unit vector (aka the direction vector or the normalized vector) based on the norm unit_vector[ii] = data[ii]/norm;static double
getVariance
(double[] data) Gets the variance of the time series.static double
getVariance
(float[] data) Gets the variance of the time series.static double
getVariance
(Collection<? extends Number> c) static double
getWindowPeakToPeak
(float[] data, int idx, int Nsamps) static boolean
hasFlatSegments
(float[] data, int segmentlength) Check whether the data series has flat segments - where every element of the segment is identicalstatic boolean
hasVariance
(float[] data, int start, int end) static float[]
hilbert
(float[] data) Calculate the Discrete Hilbert Transform of real valued data H(f=0) == data(f=0) H(f=length/2) == data(f=length/2) H(f) = ( 2 * data(f) ) : t = 1 : length/2 - 1 H(f) = 0.0 : t = length/2 +1 : length H(t) = ifft(H(f))static void
Integrate
(float[] data, double samprate) Integrate a series in the time domainstatic void
IntegrateInFreqDomain
(org.apache.commons.math3.complex.Complex[] spectrum, double dt) Integrate a complex spectrum by dividing by j*omega.static double
interpolate
(double[] x, double[] y, double xinterp) Interpolate a single y-value given x- and y-arrays.static double[]
interpolate
(double[] x, double[] y, double[] xinterp) Produce a vector of interpolated values.static double
interpolate
(double xStart, double dx, double[] y, double xinterp) Interpolate a single y-valuestatic double[]
interpolate
(double xStart, double dx, double[] y, double[] xinterp) Interpolate an array of y-valuesstatic float[]
interpolate
(double xStart, double oldDx, float[] y, double newDx) static float
interpolate
(float[] x, float[] y, float xinterp) Interpolate a single y-value given x- and y-arrays.static float[]
interpolate
(float[] x, float[] y, float[] xinterp) Produce a vector of interpolated values.static float[]
interpolate
(float xStart, float dx, float[] y, float[] xinterp) Interpolate an array of y-valuesstatic float
interpolateValue
(double xStart, double dx, float[] y, double xinterp) Interpolate a single y-valuestatic boolean
isConstant
(float[] data) static float[]
log10
(float[] data) Create an array for which each point is the log10 value of every point in the original input array Note values less than 0.static float[]
matlabDiff
(float[] data) static double[]
meanSmooth
(double[] data, int halfWidth) Smooth the input series using a sliding window of width halfWidth.static float[]
meanSmooth
(float[] data, int halfWidth) A float version of the MeanSmooth methodstatic double[]
multiply
(double[] data, double value) static void
multiplyScalar
(double[] data, double value) Multiply the series values by a scalar constant.static void
multiplyScalar
(float[] data, double value) Multiply the series values by a scalar constant.static <T extends Number>
T[]A method of padding your array with some elementsstatic double
quickMedian
(double[] values) static double
quickMedian
(float[] values) static double
quickMedian
(NumericalList values) static void
removeGlitches
(float[] data, double threshhold) A simple glitch removal method if the (data[j] - median) falls above a defined threshhold replace data[j] with the median valuestatic int
removeGlitches
(float[] data, int windowLength, double threshold) static int
removeGlitches
(float[] data, int windowLength, double MinGlitchAmp, double Thresh2) static void
removeMean
(float[] data) Remove the mean of the series.static void
removeMedian
(float[] data) static void
removeTrend
(float[] data) Remove a linear trend from a seriesstatic void
replaceValue
(float[] data, float value, float replacement) static void
reverseArray
(double[] data) Reverse the data in an array.static void
reverseArray
(float[] data) Reverse the data in an array.static void
rotate
(float[] signal1, float[] signal2, double angle, boolean npinput, boolean npoutput) Perform a clockwise rotation of a pair of signals (Based on the SAC C code rotate.c) This method replaces the original signals with their rotated counterparts.static void
setMaximumRange
(float[] data, double maxRange) static void
setPtoPAmplitude
(float[] data, double newRange) static int[]
sign
(float[] x) static float[]
signedSqrt
(float[] data) Create an array for which each point is the signed sqrt value of every point in the original input array Note values less than 0.static float[]
signedSquare
(float[] data) Create an array for which each point is the signed square value of every point in the original input array Note values less than 0.static void
signum
(double[] data) static void
signum
(float[] data) static double
slope
(float[] data) static float[]
sqrt
(float[] data) Create an array for which each point is the sqrt value of every point in the original input array Note values less than 0.static double[]
square
(double[] data) Create an array for which each point is the squared value of every point in the original input arraystatic float[]
square
(float[] data) Create an array for which each point is the squared value of every point in the original input arraystatic double[]
subtract
(double[] data1, double[] data2) Element by element addition of two data series of equal lengthstatic float[]
subtract
(float[] data1, float[] data2) Element by element addition of two data series of equal lengthstatic void
taper
(float[] data, double TaperPercent) Apply a cosine taper to an array of floatsstatic double[]
testdebugSmooth
(double[] data, int halfWidth) static float[]
testdebugSmooth
(float[] data, int halfWidth) static void
triangleTaper
(float[] data, double taperPercent) static float[]
whiten
(float[] x, double[] a) a whitening filter y(t) = x(t) + SUM(k = 1:p){a(k)*x(t-k)}
-
Constructor Details
-
SeriesMath
public SeriesMath()
-
-
Method Details
-
add
public static double[] add(double[] data1, double[] data2) Element by element addition of two data series of equal length- Parameters:
data1
- - the first data seriesdata2
- - the second data series- Returns:
- result[ii] = data1[ii] + data2[ii]
-
add
public static double[] add(double[] data, double value) Add a scalar value to all elements of the array and return a new array- Parameters:
data
-value
-- Returns:
-
add
public static float[] add(float[] data, float value) -
add
public static float[] add(float[] data, double value) -
addScalar
public static void addScalar(float[] data, double value) Add a scalar to the series.- Parameters:
value
- The scalar value to be added to the seriesdata
- Array containing the series values
-
addScalar
public static void addScalar(double[] data, double value) Add a scalar to the series.- Parameters:
value
- The scalar value to be added to the seriesdata
- Array containing the series values
-
differentiate
public static void differentiate(float[] data, double samprate) -
doubleToFloat
public static float[] doubleToFloat(double[] data) -
floatToDouble
public static double[] floatToDouble(float[] data) -
Integrate
public static void Integrate(float[] data, double samprate) Integrate a series in the time domain- Parameters:
data
- Array containing the series values assumed to be evenly sampled.samprate
- The sample rate
-
IntegrateInFreqDomain
public static void IntegrateInFreqDomain(org.apache.commons.math3.complex.Complex[] spectrum, double dt) Integrate a complex spectrum by dividing by j*omega. Leaves the value for the zero frequency unchanged.- Parameters:
spectrum
- The complex spectrum to be integrated.dt
- The sample interval in the time domain.
-
meanSmooth
public static float[] meanSmooth(float[] data, int halfWidth) A float version of the MeanSmooth method- Parameters:
data
- The data to be smoothed.halfWidth
- The half-width of the smoothing window.- Returns:
- The smoothed data.
-
meanSmooth
public static double[] meanSmooth(double[] data, int halfWidth) Smooth the input series using a sliding window of width halfWidth. At each point, the mean is computed for the set of points centered on the point and of width 2 * W -1. The value of the point is then replaced by the mean. At the ends of the series, shorter windows are used as required. The end points are left unchanged. The input series is left unmodified.- Parameters:
data
- The data to be smoothed.halfWidth
- The half-width of the smoothing window.- Returns:
- The smoothed data series.
-
multiply
public static double[] multiply(double[] data, double value) -
multiplyScalar
public static void multiplyScalar(float[] data, double value) Multiply the series values by a scalar constant.- Parameters:
value
- The scalar value with which to multiply the series valuesdata
- Array containing the series values
-
multiplyScalar
public static void multiplyScalar(double[] data, double value) Multiply the series values by a scalar constant.- Parameters:
value
- The scalar value with which to multiply the series valuesdata
- Array containing the series values
-
removeMean
public static void removeMean(float[] data) Remove the mean of the series.- Parameters:
data
- Array containing the series values
-
removeTrend
public static void removeTrend(float[] data) Remove a linear trend from a series- Parameters:
data
- Array containing the series values
-
reverseArray
public static void reverseArray(float[] data) Reverse the data in an array.- Parameters:
data
- Description of the Parameter
-
reverseArray
public static void reverseArray(double[] data) Reverse the data in an array.- Parameters:
data
- Description of the Parameter
-
slope
public static double slope(float[] data) -
subtract
public static double[] subtract(double[] data1, double[] data2) Element by element addition of two data series of equal length- Parameters:
data1
- - the first data seriesdata2
- - the second data series- Returns:
- result[ii] = data1[ii] - data2[ii]
-
subtract
public static float[] subtract(float[] data1, float[] data2) Element by element addition of two data series of equal length- Parameters:
data1
- - the first data seriesdata2
- - the second data series- Returns:
- result[ii] = data1[ii] - data2[ii]
-
taper
public static void taper(float[] data, double TaperPercent) Apply a cosine taper to an array of floats- Parameters:
TaperPercent
- The (one-sided) percent of the time series to which a taper will be applied. The value ranges from 0 (no taper) to 50 ( The taper extends half the length of the series ). Since the taper is symmetric, a 50% taper means that all but the center value of the series will be scaled by some value less than 1.0.data
- Array containing the series values
-
abs
public static float[] abs(float[] data) convert all the elements of the data series to their absolute value- Parameters:
data
-- Returns:
-
abs
public static double[] abs(double[] data) convert all the elements of the data series to their absolute value- Parameters:
data
-- Returns:
-
arrayMultiply
public static float[] arrayMultiply(float[] array1, float[] array2) -
autocorrelate
public static double[] autocorrelate(float[] x) Calculate the autocorrelation of a data series- Parameters:
x
- - the data series- Returns:
- The autocorrelation of the input sequence.
-
autocorrelate
public static double autocorrelate(float[] x, int shift) Calculate the autocorrelation at a specific offset- Parameters:
x
- - the data seriesshift
- - the specific point desired Note autocorrelate(data,0) gives the autocorrelation value at 0 shift also known as the "energy" of the trace- Returns:
- The value of the autocorrelation at shift.
-
autocorrelate
public static double[] autocorrelate(double[] x) Calculate the autocorrelation of a data series- Parameters:
x
- - the data series
-
autocorrelate
public static double autocorrelate(double[] x, int shift) Calculate the autocorrelation at a specific offset- Parameters:
x
- - the data seriesshift
- - the specific point desired Note autocorrelate(data,0) gives the autocorrelation value at 0 shift also known as the "energy" of the trace- Returns:
- autocorrelation.
-
countGlitches
public static int countGlitches(float[] data, int windowLength, double threshold) -
crosscorrelate
public static double[] crosscorrelate(float[] x, float[] y) Get the cross-correlation series of two data sets for shift = 0:size result[shift] = sum( x(t) * y(t+shift))- Parameters:
x
- the first data seriesy
- the second data series- Returns:
- the cross-correlation series
-
crosscorrelate
public static double crosscorrelate(float[] x, float[] y, int shift) Get the value of the cross-correlation series at a specific offset result[shift] = sum( x(t) * y(t+shift) )- Parameters:
x
- the first data seriesy
- the second data seriesshift
- - the offset between the data series- Returns:
- the cross correlation at a specific point Note the energy of a trace is the autocorrelation at zero shift
-
crosscorrelate
public static double[] crosscorrelate(double[] x, double[] y) Get the cross-correlation series of two data sets for shift = 0:size result[shift] = sum( x(t) * y(t+shift))- Parameters:
x
- the first data seriesy
- the second data series- Returns:
- the cross-correlation series
-
crosscorrelate
public static double crosscorrelate(double[] x, double[] y, int shift) Get the value of the cross-correlation series at a specific offset result[shift] = sum( x(t) * y(t+shift) )- Parameters:
x
- the first data seriesy
- the second data seriesshift
- - the offset between the data series- Returns:
- the cross correlation at a specific point Note the energy of a trace is the autocorrelation at zero shift
-
cut
public static double[] cut(double[] data, double minVal, double maxVal) Produce a cut version of the input series. The cut is made such that no value in the series is less than minVal and no value in the series is greater than maxVal. Assumes that the input series is sorted from least to greatest. Note that if the minVal and maxVal are such that all values in the input series are retained, then the input is simply passed out as the result.- Parameters:
data
- The data to be cut.minVal
- The value bounding the cut series on the low side.maxVal
- The value bounding the series on the high side.- Returns:
- The cut series.
-
decimate
public static float[] decimate(float[] data, int decimationfactor) From the SeriesMath method given a data series it is decimated so that only every Nth point is retained This should be equivalent to the more versatile interpolate methods where N is the decimationfactor- Parameters:
data
- the original data series arraydecimationfactor
- should be an integer greater than 2- Returns:
- the decimated series
-
decimate
public static double[] decimate(double[] data, int decimationfactor) A double value version of the decimation method- Parameters:
data
-decimationfactor
-- Returns:
- the decimated series
-
dotProduct
dot product of two equal length Series dot(A,B) = SUM(AiBi) = A1B1 + A2B2+...AnBn- Parameters:
data1
-data2
-- Returns:
-
dotProduct
dot product of two equal length Series dot(A,B) = SUM(AiBi) = A1B1 + A2B2+...AnBn- Parameters:
data1
-data2
-- Returns:
-
elementAddition
public static float[] elementAddition(float[] data1, float[] data2) adds data1's values to data2's- Parameters:
data1
-data2
-- Returns:
-
elementDivision
public static float[] elementDivision(float[] data1, float[] data2) divide the individual elements of one data series by those of another- Parameters:
data1
- - the first data seriesdata2
- - the denominator data series must have the same length as the first- Returns:
- a new double[] series in which result[ii] = data1[ii]/data2[ii] Akin to array left division in Matlab A.\B
-
elementDivision
public static double[] elementDivision(double[] data1, double[] data2) divide the individual elements of one data series by those of another- Parameters:
data1
- - the first data seriesdata2
- - the denominator data series must have the same length as the first- Returns:
- a new double[] series in which result[ii] = data1[ii]/data2[ii] Akin to array left division in Matlab A.\B
-
elementMultiplication
public static float[] elementMultiplication(float[] data1, float[] data2) Element-by-element multiplication, C = A.*B- Parameters:
data1
-data2
-- Returns:
- result[ii] = data1[ii] * data2[ii]
-
elementMultiplication
public static double[] elementMultiplication(double[] data1, double[] data2) Element-by-element multiplication, C = A.*B- Parameters:
data1
-data2
-- Returns:
- result[ii] = data1[ii] * data2[ii]
-
elementSubtraction
public static float[] elementSubtraction(float[] data1, float[] data2) Subtracts from data1 data2's values- Parameters:
data1
-data2
-- Returns:
-
envelope
public static float[] envelope(float[] data) Calculate the envelope of real valued data using a Hilbert transform E(t) = sqrt( data^2 + Hilbert(data) ^2 )- Parameters:
data
- : the data series- Returns:
- the envelope function for the data
-
find
Finds the index of values meeting specified conditions. Each element of array is tested against value using specified logical ,condition (==,LTEQ,GTEQ,LT,GT) are used in conjuntion with value usage e.g. find( longarray, "greater than equal to", 7); -
find
Finds the index of values meeting specified conditions. Each element of array is tested against value using specified logical ,condition (==,LTEQ,GTEQ,LT,GT) are used in conjuntion with value -
find
Finds the index of values meeting specified conditions. Each element of array is tested against value using specified logical ,condition (==,LTEQ,GTEQ,LT,GT) are used in conjuntion with value -
find
Finds the index of values meeting specified conditions. Each element of array is tested against value using specified logical ,condition (==,LTEQ,GTEQ,LT,GT) are used in conjuntion with value- Parameters:
array
-condition
-value
-- Returns:
-
find
For 2 arrays of equal length find the elements of array 1 that are ("GT", "LT", "==" ...) the elements of array 2- Parameters:
array1
- an array of doublescondition
-array2
- a second array of doubles- Returns:
- an array containing the indices for each element that passes the conditional test e.g. array1 = { 1.1, 2.2, 3.3}, array2 = {4.4, 2.2, 3.3} find(array1, "==" array2) will return result = {1, 2}
-
findDiscontinuities
public static DiscontinuityCollection findDiscontinuities(float[] data, double sampRate, int winLength, double factor) -
findIndicesOf
Find the indices in the data series that match a specific value- Parameters:
data
- - an array of doublesvalue
- - the specific double value you are looking for- Returns:
- - a Vector of all the indices which match that value
-
getDefiniteIntegral
public static double getDefiniteIntegral(float[] data, double delta, int idx1, int idx2) Compute definite integral using extended trapezoidal rule.- Parameters:
data
- The evenly-spaced function to integrate.delta
- The sample interval.idx1
- The starting index over which to integrate.idx2
- The ending index over which to integrate.- Returns:
- The value of the definite integral.
-
getDefiniteIntegral
-
getEnergy
public static double getEnergy(float[] x) -
getExtremum
public static float getExtremum(float[] data) - Parameters:
data
-- Returns:
- the largest absolute value in the series
-
getFirstDifference
public static float[] getFirstDifference(float[] data) -
getGlitches
-
getKurtosis
public static double getKurtosis(float[] data) -
getMax
public static float getMax(float[] data, int idx1, int idx2) -
getMax
public static float getMax(float[] data) Gets the maximum value of the series.- Parameters:
data
- Array containing the series values- Returns:
- The max value of the series
-
getMax
public static double getMax(double[] data) Gets the maximum value of the series.- Parameters:
data
- Array containing the series values- Returns:
- The max value of the series
-
getMax
-
getMaxIndex
Gets the maximum value of the series and the index it occurs at.- Parameters:
data
- Array containing the series values- Returns:
- The (index, max value) of the series
-
getMaxIndex
Gets the maximum value of the series and the index it occurs at.- Parameters:
data
- Array containing the series values- Returns:
- The (index, max value) of the series
-
getMean
public static double getMean(float[] data) Gets the mean value of the series.- Parameters:
data
- Array containing the series values- Returns:
- The mean value of the series
-
getMean
-
getMean
public static double getMean(double[] data) Gets the mean value of the series.- Parameters:
data
- Array containing the series values- Returns:
- The mean value of the series
-
getMedian
public static double getMedian(float[] data) Gets the median value of the time series.- Parameters:
data
- Array containing the series values- Returns:
- The median value of the series
-
getMedian
public static double getMedian(double[] values) -
getMedian
-
getMin
-
getMin
public static float getMin(float[] data, int idx1, int idx2) -
getMin
public static float getMin(float[] data) Gets the minimum value of the series.- Parameters:
data
- Array containing the series values- Returns:
- The min value of the series
-
getMin
public static double getMin(double[] data) Gets the minimum value of the series.- Parameters:
data
- Array containing the series values- Returns:
- The min value of the series
-
getMinIndex
Gets the maximum value of the series and the index it occurs at.- Parameters:
data
- Array containing the series values- Returns:
- The (index, min value) of the series
-
getMinIndex
Gets the minimum value of the series and the index it occurs at.- Parameters:
data
- Array containing the series values- Returns:
- The (index, max value) of the series
-
getNorm
public static double getNorm(double[] data) Treat the series as a vector and calculate the norm (aka the direction or the normal) of the vector- Parameters:
data
-- Returns:
-
getPeakToPeakAmplitude
public static double getPeakToPeakAmplitude(float[] data) -
getPeakToPeakAmplitude
public static double getPeakToPeakAmplitude(float[] data, double sampleInterval, double period) Gets the peakToPeakAmplitude attribute of an input timeseries at the specified period. Advances a sliding window of length period through the time series a point at a time. At each position, the Peak-To-Peak range of the current window is computed. The maximum value of all these values is returned.- Parameters:
period
- The period in seconds at which to compute the value.data
- An array of floats whose Peak-To-Peak value is to be determined.sampleInterval
- The sample interval of the data in the data array.- Returns:
- The maximum Peak-To-Peak value for the array.
-
getRMS
public static double getRMS(float[] data) Gets the RMS value of the data in the input float array.- Parameters:
data
- The input float array.- Returns:
- The RMS value
-
getRMS
public static double getRMS(double[] data) Gets the RMS value of the data in the input float array.- Parameters:
data
- The input float array.- Returns:
- The RMS value
-
getRMS
-
getRange
public static double getRange(float[] u) -
getSnr
public static double getSnr(float[] data, double sampleRate, double dataStartTime, double pickEpochTime, double prePickWindowLength, double postPickWindowLength) -
getSnr
public static double getSnr(float[] data, int pickIndex, int preWinStartIndex, int postWinEndIndex) -
getStDev
public static double getStDev(float[] data) Gets the standard deviation of the time series.- Parameters:
data
- Array containing the series values- Returns:
- The standard deviation value
-
getStDev
public static double getStDev(double[] data) Gets the standard deviation of the time series.- Parameters:
data
- Array containing the series values- Returns:
- The standard deviation value
-
getStandardDeviation
-
getStd
public static double getStd(double[] vec, int idx1, int idx2) -
getSubSection
public static double[] getSubSection(double[] v, int minIdx, int maxIdx) -
getSubSection
public static float[] getSubSection(float[] v, int minIdx, int maxIdx) -
getSubSection
public static int[] getSubSection(int[] v, int minIdx, int maxIdx) -
getSum
public static double getSum(float[] data) Gets the sum of the values in the input array- Parameters:
data
- Array containing the series values- Returns:
- The sum of the values
-
getSum
-
getSum
public static double getSum(double[] data) Gets the sum of the values in the input array- Parameters:
data
- Array containing the series values- Returns:
- The sum of the values
-
getSumOfSquares
public static double getSumOfSquares(float[] data) -
getUnitVector
public static double[] getUnitVector(double[] data) Treat the series as a vector and calculate the unit vector (aka the direction vector or the normalized vector) based on the norm unit_vector[ii] = data[ii]/norm;- Parameters:
data
- is an Array containing the series values
-
getVariance
public static double getVariance(float[] data) Gets the variance of the time series.- Parameters:
data
- Array containing the series values- Returns:
- The variance value
-
getVariance
-
getVariance
public static double getVariance(double[] data) Gets the variance of the time series.- Parameters:
data
- Array containing the series values- Returns:
- The variance value
-
getWindowPeakToPeak
public static double getWindowPeakToPeak(float[] data, int idx, int Nsamps) -
hasVariance
public static boolean hasVariance(float[] data, int start, int end) -
hilbert
public static float[] hilbert(float[] data) Calculate the Discrete Hilbert Transform of real valued data H(f=0) == data(f=0) H(f=length/2) == data(f=length/2) H(f) = ( 2 * data(f) ) : t = 1 : length/2 - 1 H(f) = 0.0 : t = length/2 +1 : length H(t) = ifft(H(f))- Parameters:
data
- : the data series- Returns:
- the complex valued Hilbert transformation (H(t))
-
interpolate
public static double interpolate(double[] x, double[] y, double xinterp) Interpolate a single y-value given x- and y-arrays.- Parameters:
x
- The array of x-values, assumed to be monitonically increasing but not-necessarily evenly-spaced..y
- The array of y-values corresponding to the samples in the x-array.xinterp
- The x-value at which an interpolated y-value is to be computed. xinterp is assumed to be in the range spanned by x.- Returns:
- The interpolated y-value.
-
interpolate
public static float interpolate(float[] x, float[] y, float xinterp) Interpolate a single y-value given x- and y-arrays.- Parameters:
x
- The array of x-values, assumed to be monitonically increasing but not-necessarily evenly-spaced..y
- The array of y-values corresponding to the samples in the x-array.xinterp
- The x-value at which an interpolated y-value is to be computed. xinterp is assumed to be in the range spanned by x.- Returns:
- The interpolated y-value.
-
interpolate
public static double interpolate(double xStart, double dx, double[] y, double xinterp) Interpolate a single y-value- Parameters:
xStart
- The x-value corresponding to the first y-value.dx
- The sample interval of the data.y
- The array of y-values.xinterp
- The x-value at which an interpolated y-value is to be computed. xinterp is assumed to be in the range spanned by x, e.g. xStart LTEQ xinterp LTEQ (y.length-1) * dx. @return The interpolated y-value.
-
interpolate
public static double[] interpolate(double[] x, double[] y, double[] xinterp) Produce a vector of interpolated values.- Parameters:
x
- The array of x-values, assumed to be monotonically increasing but not-necessarily evenly-spaced..y
- The array of y-values corresponding to the samples in the x-array.xinterp
- The array of x-values for which interpolated y-values are to be computed. All values in xinterp are assumed to be within the range spanned by x.- Returns:
- The array of interpolated y-values.
-
interpolate
public static float[] interpolate(float[] x, float[] y, float[] xinterp) Produce a vector of interpolated values.- Parameters:
x
- The array of x-values, assumed to be monotonically increasing but not-necessarily evenly-spaced..y
- The array of y-values corresponding to the samples in the x-array.xinterp
- The array of x-values for which interpolated y-values are to be computed. All values in xinterp are assumed to be within the range spanned by x.- Returns:
- The array of interpolated y-values.
-
interpolate
public static double[] interpolate(double xStart, double dx, double[] y, double[] xinterp) Interpolate an array of y-values- Parameters:
xStart
- The x-value corresponding to the first y-value.dx
- The sample interval of the data.y
- The array of y-values.xinterp
- An array of x-values for which y-values are to be interpolated. All values in xinterp must be in the range defined by xStart, dx, and y.length.- Returns:
- The interpolated y-value.
-
interpolate
public static float[] interpolate(float xStart, float dx, float[] y, float[] xinterp) Interpolate an array of y-values- Parameters:
xStart
- The x-value corresponding to the first y-value.dx
- The sample interval of the data.y
- The array of y-values.xinterp
- An array of x-values for which y-values are to be interpolated. All values in xinterp must be in the range defined by xStart, dx, and y.length.- Returns:
- The interpolated y-value.
-
interpolate
public static float[] interpolate(double xStart, double oldDx, float[] y, double newDx) -
interpolateValue
public static float interpolateValue(double xStart, double dx, float[] y, double xinterp) Interpolate a single y-value- Parameters:
xStart
- The x-value corresponding to the first y-value.dx
- The sample interval of the data.y
- The array of y-values.xinterp
- The x-value at which an interpolated y-value is to be computed. xinterp is assumed to be in the range spanned by x, e.g. xStart LTEQ xinterp LTEQ (y.length-1) * dx. @return The interpolated y-value.
-
isConstant
public static boolean isConstant(float[] data) -
log10
public static float[] log10(float[] data) Create an array for which each point is the log10 value of every point in the original input array Note values less than 0. are replaced by 0.- Parameters:
data
-- Returns:
- the log10 array
-
matlabDiff
public static float[] matlabDiff(float[] data) -
padRight
public static <T extends Number> T[] padRight(T[] original, T padElement, int newLength, Class<T> javaSucks) A method of padding your array with some elements- Type Parameters:
T
-- Parameters:
original
-padElement
- what to fill the padded elements withnewLength
- the length to pad to.javaSucks
- because it throws away the generic type information forcing us to pass it this kludgy way!- Returns:
-
quickMedian
public static double quickMedian(float[] values) -
quickMedian
public static double quickMedian(double[] values) -
quickMedian
-
removeGlitches
public static int removeGlitches(float[] data, int windowLength, double threshold) -
removeGlitches
public static int removeGlitches(float[] data, int windowLength, double MinGlitchAmp, double Thresh2) -
removeGlitches
public static void removeGlitches(float[] data, double threshhold) A simple glitch removal method if the (data[j] - median) falls above a defined threshhold replace data[j] with the median value- Parameters:
data
- The data array to be deglitched.threshhold
- The threshold in terms of signal variance.
-
hasFlatSegments
public static boolean hasFlatSegments(float[] data, int segmentlength) Check whether the data series has flat segments - where every element of the segment is identical- Parameters:
data
-segmentlength
- an integer number of elements defining the shortest segment- Returns:
- true if there are any single-valued segments
-
removeMedian
public static void removeMedian(float[] data) -
replaceValue
public static void replaceValue(float[] data, float value, float replacement) -
rotate
public static void rotate(float[] signal1, float[] signal2, double angle, boolean npinput, boolean npoutput) Perform a clockwise rotation of a pair of signals (Based on the SAC C code rotate.c) This method replaces the original signals with their rotated counterparts.- Parameters:
signal1
- - first input signalsignal2
- - second input signalangle
- - angle of rotation (clockwise from direction of signal1) in degreesnpinput
- - TRUE if the input signals have "normal" polaritynpoutput
- - TRUE if the output signals have "normal" polarity "normal" polarity is such that the second component leads the first component by 90 degrees in a clockwise rotation.
-
setMaximumRange
public static void setMaximumRange(float[] data, double maxRange) -
setPtoPAmplitude
public static void setPtoPAmplitude(float[] data, double newRange) -
sign
public static int[] sign(float[] x) -
signedSqrt
public static float[] signedSqrt(float[] data) Create an array for which each point is the signed sqrt value of every point in the original input array Note values less than 0. are replaced by -1* sqrt(abs(data))- Parameters:
data
-- Returns:
- the signed sqrt array
-
signedSquare
public static float[] signedSquare(float[] data) Create an array for which each point is the signed square value of every point in the original input array Note values less than 0. are replaced by -1*data*data- Parameters:
data
-- Returns:
- the signed squared value array
-
signum
public static void signum(float[] data) -
signum
public static void signum(double[] data) -
sqrt
public static float[] sqrt(float[] data) Create an array for which each point is the sqrt value of every point in the original input array Note values less than 0. are replaced by 0.- Parameters:
data
-- Returns:
- the sqrt array
-
square
public static float[] square(float[] data) Create an array for which each point is the squared value of every point in the original input array- Parameters:
data
-- Returns:
- the squared value array
-
square
public static double[] square(double[] data) Create an array for which each point is the squared value of every point in the original input array- Parameters:
data
-- Returns:
- the squared value array
-
testdebugSmooth
public static float[] testdebugSmooth(float[] data, int halfWidth) -
testdebugSmooth
public static double[] testdebugSmooth(double[] data, int halfWidth) - Parameters:
data
- The data to be smoothed.halfWidth
- The half-width of the smoothing window.- Returns:
- The smoothed data series.
-
triangleTaper
public static void triangleTaper(float[] data, double taperPercent) -
whiten
public static float[] whiten(float[] x, double[] a) a whitening filter y(t) = x(t) + SUM(k = 1:p){a(k)*x(t-k)}- Parameters:
x
- the initial data seta
- a vector of whitening coefficients- Returns:
- the whitened data
-
covariance
public double covariance(double[] vec1, double[] vec2) Calculates the covariance of two equally sized vectors. cov(v1,v2) = 1/(n-1)*sum(xi-xavg)*sum(yi-yavg)- Parameters:
vec1
-vec2
-- Returns:
- the covariance between the provided vectors. Value will be in the range [-1,1]
-