version  0.0.1
Defines the C++ API for MsPASS
Public Member Functions | List of all members
mspass::algorithms::deconvolution::WaterLevelDecon Class Reference
Inheritance diagram for mspass::algorithms::deconvolution::WaterLevelDecon:
mspass::algorithms::deconvolution::FFTDeconOperator mspass::algorithms::deconvolution::ScalarDecon mspass::algorithms::deconvolution::BasicDeconOperator

Public Member Functions

 WaterLevelDecon (const WaterLevelDecon &parent)
 
 WaterLevelDecon (const mspass::utility::Metadata &md)
 
 WaterLevelDecon (const mspass::utility::Metadata &md, const std::vector< double > &wavelet, const std::vector< double > &data)
 
void changeparameter (const mspass::utility::Metadata &md)
 
void process ()
 
mspass::seismic::CoreTimeSeries actual_output ()
 Return the actual output of the deconvolution operator. More...
 
mspass::seismic::CoreTimeSeries inverse_wavelet (const double t0parent=0.0)
 Return a FIR respresentation of the inverse filter. More...
 
mspass::seismic::CoreTimeSeries inverse_wavelet ()
 Return default FIR represesentation of the inverse filter. More...
 
mspass::utility::Metadata QCMetrics ()
 Return appropriate quality measures. More...
 
- Public Member Functions inherited from mspass::algorithms::deconvolution::FFTDeconOperator
 FFTDeconOperator (const mspass::utility::Metadata &md)
 
 FFTDeconOperator (const FFTDeconOperator &parent)
 
FFTDeconOperatoroperator= (const FFTDeconOperator &parent)
 
void changeparameter (const mspass::utility::Metadata &md)
 
void change_size (const int nfft_new)
 
void change_shift (const int shift)
 
int get_size ()
 
int get_shift ()
 
int operator_size ()
 
int operator_shift ()
 
double df (const double dt)
 
mspass::seismic::CoreTimeSeries FourierInverse (const ComplexArray &winv, const ComplexArray &sw, const double dt, const double t0parent)
 Return inverse wavelet for Fourier methods. More...
 
- Public Member Functions inherited from mspass::algorithms::deconvolution::ScalarDecon
 ScalarDecon (const mspass::utility::Metadata &md)
 
 ScalarDecon (const std::vector< double > &d, const std::vector< double > &w)
 
 ScalarDecon (const ScalarDecon &parent)
 
int load (const std::vector< double > &wavelet, const std::vector< double > &data)
 Load all data required for decon. More...
 
int loaddata (const std::vector< double > &data)
 
int loadwavelet (const std::vector< double > &wavelet)
 
ScalarDeconoperator= (const ScalarDecon &parent)
 
std::vector< double > getresult ()
 
void change_shaping_wavelet (const ShapingWavelet &nsw)
 
mspass::seismic::CoreTimeSeries ideal_output ()
 

Additional Inherited Members

- Protected Attributes inherited from mspass::algorithms::deconvolution::FFTDeconOperator
int nfft
 
int sample_shift
 
gsl_fft_complex_wavetable * wavetable
 
gsl_fft_complex_workspace * workspace
 
ComplexArray winv
 
- Protected Attributes inherited from mspass::algorithms::deconvolution::ScalarDecon
std::vector< double > data
 
std::vector< double > wavelet
 
std::vector< double > result
 
ShapingWavelet shapingwavelet
 

Member Function Documentation

◆ actual_output()

CoreTimeSeries mspass::algorithms::deconvolution::WaterLevelDecon::actual_output ( )
virtual

Return the actual output of the deconvolution operator.

The actual output is defined as w^-1*w and is compable to resolution kernels in linear inverse theory. Although not required we would normally expect this function to be peaked at 0. Offsets from 0 would imply a bias.

Implements mspass::algorithms::deconvolution::ScalarDecon.

138 {
139  try {
140  ComplexArray W(nfft,&(wavelet[0]));
141  gsl_fft_complex_forward(W.ptr(),1,nfft,wavetable,workspace);
142  ComplexArray ao_fft;
143  ao_fft=winv*W;
144  /* We always apply the shaping wavelet - this perhaps should be optional
145  but probably better done with a none option for the shaping wavelet */
146  ao_fft=(*shapingwavelet.wavelet())*ao_fft;
147  gsl_fft_complex_inverse(ao_fft.ptr(),1,nfft,wavetable,workspace);
148  vector<double> ao;
149  ao.reserve(nfft);
150  for(unsigned int k=0; k<ao_fft.size(); ++k) ao.push_back(ao_fft[k].real());
151  /* We always shift this wavelet to the center of the data vector.
152  We handle the time through the CoreTimeSeries object. */
153  int i0=nfft/2;
154  ao=circular_shift(ao,i0);
155  CoreTimeSeries result(nfft);
156  /* Getting dt from here is unquestionably a flaw in the api, but will
157  retain for now. Perhaps should a copy of dt in the ScalarDecon object. */
158  double dt=this->shapingwavelet.sample_interval();
159  /* t0 is time of sample zero - hence normally negative*/
160  /*Old API
161  result.t0=dt*(-(double)i0);
162  result.dt=dt;
163  result.live=true;
164  result.tref=TimeReferenceType::Relative;
165  result.ns=nfft;
166  */
167 
168  result.set_t0(-dt*((double)i0));
169  result.set_dt(dt);
170  result.set_live();
171  result.set_tref(TimeReferenceType::Relative);
172  result.set_npts(nfft);
173  for(int k=0;k<nfft;++k)result.s[k]=ao[k];;
174  return result;
175  } catch(...) {
176  throw;
177  };
178 }
ComplexArray * wavelet()
Definition: ShapingWavelet.h:74
Scalar time series data object.
Definition: CoreTimeSeries.h:18

References mspass::seismic::CoreTimeSeries::s, mspass::seismic::CoreTimeSeries::set_dt(), mspass::seismic::BasicTimeSeries::set_live(), mspass::seismic::CoreTimeSeries::set_npts(), mspass::seismic::CoreTimeSeries::set_t0(), and mspass::seismic::BasicTimeSeries::set_tref().

◆ inverse_wavelet() [1/2]

CoreTimeSeries mspass::algorithms::deconvolution::WaterLevelDecon::inverse_wavelet ( )
virtual

Return default FIR represesentation of the inverse filter.

This is an overloaded version of the parameterized method. It is equivalent to this->inverse_wavelet(0.0,0.0);

Implements mspass::algorithms::deconvolution::ScalarDecon.

193 {
194  try{
195  return this->inverse_wavelet(0.0);
196  }catch(...){throw;};
197 }
mspass::seismic::CoreTimeSeries inverse_wavelet()
Return default FIR represesentation of the inverse filter.
Definition: WaterLevelDecon.cc:192

◆ inverse_wavelet() [2/2]

CoreTimeSeries mspass::algorithms::deconvolution::WaterLevelDecon::inverse_wavelet ( const double  t0parent = 0.0)
virtual

Return a FIR respresentation of the inverse filter.

An inverse filter has an impulse response. For some wavelets this can be respresented by a FIR filter with finite numbers of coefficients. Since this is a Fourier method the best we can do is return the inverse fft of the regularized operator. The output usually needs to be phase shifted to be most useful. For typical seismic source wavelets that are approximately minimum phase the shift can be small, but for zero phase input it should be approximately half the window size. This method also has an optional argument for t0parent. Because this processor was written to be agnostic about a time standard it implicitly assumes time 0 is sample 0 of the input waveforms. If the original data have a nonzero start time this should be passed as t0parent or the output will contain a time shift of t0parent. Note that tshift and t0parent do very different things. tshift is used to apply circular phase shift to the output (e.g. a shift of 10 samples causes the last 10 samples in the wavelet to be wrapped to the first 10 samples). t0parent only changes the time standard so the output has t0 -= parent.t0.

Output wavelet is always circular shifted with 0 lag at center.

Parameters
t0parent- time zero of parent seismograms (see above).
Exceptions
SeisppErroris thrown if the tshift value is more than half the length of the data (nfft*dt). Reason is justified above.

Implements mspass::algorithms::deconvolution::ScalarDecon.

181 {
182  try {
183  /* Getting dt from here is unquestionably a flaw in the api, but will
184  * retain for now. Perhaps should a copy of dt in the ScalarDecon object. */
185  double dt=this->shapingwavelet.sample_interval();
186  return (this->FFTDeconOperator::FourierInverse(this->winv,
187  *shapingwavelet.wavelet(),dt,t0parent));
188  } catch(...) {
189  throw;
190  };
191 }
mspass::seismic::CoreTimeSeries FourierInverse(const ComplexArray &winv, const ComplexArray &sw, const double dt, const double t0parent)
Return inverse wavelet for Fourier methods.
Definition: FFTDeconOperator.cc:119

◆ QCMetrics()

Metadata mspass::algorithms::deconvolution::WaterLevelDecon::QCMetrics ( )
virtual

Return appropriate quality measures.

Each operator commonly has different was to measure the quality of the result. This method should return these in a generic Metadata object.

Implements mspass::algorithms::deconvolution::ScalarDecon.

199 {
200  try {
201  Metadata md;
202  md.put("underwater_fraction",regularization_fraction);
203  return md;
204  } catch(...) {
205  throw;
206  };
207 }
Definition: Metadata.h:76

The documentation for this class was generated from the following files: