version  0.0.1
Defines the C++ API for MsPASS
Loading...
Searching...
No Matches
CNRDeconEngine.h
1#ifndef __CNR_DECON_ENGINE_H__
2#define __CNR_DECON_ENGINE_H__
3#include <boost/archive/text_iarchive.hpp>
4#include <boost/archive/text_oarchive.hpp>
5#include <boost/serialization/base_object.hpp>
6#include <boost/serialization/shared_ptr.hpp>
7
8#include "mspass/algorithms/Taper.h"
9#include "mspass/algorithms/deconvolution/FFTDeconOperator.h"
10#include "mspass/algorithms/deconvolution/MTPowerSpectrumEngine.h"
11#include "mspass/algorithms/deconvolution/ShapingWavelet.h"
12#include "mspass/seismic/PowerSpectrum.h"
13#include "mspass/seismic/Seismogram.h"
14#include "mspass/seismic/TimeSeries.h"
15#include "mspass/utility/AntelopePf.h"
16
17namespace mspass::algorithms::deconvolution {
18/* This enum is file scope to intentionally exclude it from python wrappers.
19It is used internally to define the algorithm the processor is to run.
20I (glp) chose that approach over the inheritance approach used in the scalar
21methods as an artistic choice. It is a matter of opinion which approach
22is better. This makes one symbol do multiple things with changes done in
23the parameter setup as opposed to having to select the right symbolic name
24to construct. Anyway, this enum defines algorithms that can be chosen for
25processing.
26*/
27enum class CNR3C_algorithms { generalized_water_level, colored_noise_damping };
29public:
31 /* design note - delete when finished.
32
33 The constructor uses the pf to initialize operator properties.
34 Important in that is the multitaper engine that has a significant
35 initialization cost. the initialize_inverse_operator is a messy
36 way to define load data and then compute the inverse stored
37 internally. Doing it that way, however, allows the same code
38 to be used for both single station and array decon. In single
39 station use every datum has to call initiallize_inverse_operator
40 while with an array decon the operator is define once for an
41 ensemble and applied to all members. THE ASSUMPTION is all the
42 grungy work to assure all that gets done correction is handled
43 in python.
44 */
63 CNRDeconEngine(const CNRDeconEngine &parent);
64 void
65 initialize_inverse_operator(const mspass::seismic::TimeSeries &wavelet,
66 const mspass::seismic::TimeSeries &noise_data);
67 void initialize_inverse_operator(
68 const mspass::seismic::TimeSeries &wavelet,
69 const mspass::seismic::PowerSpectrum &noise_spectrum);
70 virtual ~CNRDeconEngine() {};
72 process(const mspass::seismic::Seismogram &d,
73 const mspass::seismic::PowerSpectrum &psnoise, const double fl,
74 const double fh);
75 double get_operator_dt() const { return this->operator_dt; }
91 mspass::seismic::TimeSeries ideal_output();
93 actual_output(const mspass::seismic::TimeSeries &wavelet);
95 inverse_wavelet(const mspass::seismic::TimeSeries &wavelet,
96 const double t0shift);
97 mspass::utility::Metadata QCMetrics();
98 CNRDeconEngine &operator=(const CNRDeconEngine &parent);
99
100private:
101 CNR3C_algorithms algorithm;
102 /* For the colored noise damping algorithm the damper is frequency dependent.
103 The same issue in water level that requires a floor on the water level
104 applies to damping. We use noise_floor to create a lower bound on
105 damper values. Note the damping constant at each frequency is
106 damp*noise except where noise is below noise_floor defined relative to
107 maximum noise value where it is set to n_peak*noise_floor*damp. */
108 double damp;
109 double noise_floor;
110 /* SNR bandbwidth estimates count frequencies with snr above this value */
111 double band_snr_floor;
112 double operator_dt; // Data must match this sample interval
113 int shaping_wavelet_number_poles;
115 /* Expected time window size in samples. When signal lengths
116 match this value the slepian tapers are not recomputed. When there
117 is a mismatch it will change. That means this can change dynamically
118 when run on multiple data objects. */
119 int winlength;
122
123 /* This algorithm uses a mix of damping and water level. Above this floor,
124 which acts a bit like a water level, no regularization is done. If
125 snr is less than this value we regularize with damp*noise_amplitude.
126 Note the noise_floor parameter puts a lower bound on the frequency dependent
127 regularization. If noise amplitude (not power) is less than noise_floor
128 the floor is set like a water level as noise_max*noise_level.*/
129 double snr_regularization_floor;
130 /* These are QC metrics computed by process method. Saved to allow them
131 to be use in QCmetrics method. */
132 double regularization_bandwidth_fraction;
133 double peak_snr[3];
134 double signal_bandwidth_fraction[3];
136 /* This is the lag from sample 0 for the time defines as 0 for the
137 wavelet used to compute the inverse. It is needed to resolve time
138 in the actual_output method.*/
139 int winv_t0_lag;
140 /*** Private methods *****/
141 void update_shaping_wavelet(const double fl, const double fh);
142 /* These are two algorithms for computing inverse operator in the frequency
143 * domain*/
144 void compute_winv(const mspass::seismic::TimeSeries &wavelet,
145 const mspass::seismic::PowerSpectrum &psnoise);
146 void compute_gwl_inverse(const mspass::seismic::TimeSeries &wavelet,
147 const mspass::seismic::PowerSpectrum &psnoise);
148 void compute_gdamp_inverse(const mspass::seismic::TimeSeries &wavelet,
149 const mspass::seismic::PowerSpectrum &psnoise);
150 friend boost::serialization::access;
151 template <class Archive>
152 void serialize(Archive &ar, const unsigned int version) {
153 // std::cout <<"Entered serialize function"<<std::endl;
154 ar &boost::serialization::base_object<FFTDeconOperator>(*this);
155 // std::cout << "Serializing first group of simple parameters"<<std::endl;
156 ar & algorithm;
157 ar & damp;
158 ar & noise_floor;
159 ar & band_snr_floor;
160 ar & operator_dt;
161 ar & shaping_wavelet_number_poles;
162 // std::cout << "Serializing shapingwavelet"<<std::endl;
163 ar & shapingwavelet;
164 ar & winlength;
165 // std::cout<<"Serializing power spectrum engine objects"<<std::endl;
166 ar & signal_engine;
167 ar & noise_engine;
168 ar & snr_regularization_floor;
169 // std::cout << "Serializin winv vector"<<std::endl;
170 ar & winv;
171 ar & winv_t0_lag;
172 // std::cout<<"Serializing final block of parameters"<<std::endl;
173 ar & regularization_bandwidth_fraction;
174 /* These fixed length arrays caused probems - seg faults.
175 * Apparently boost doesn't handle that corectly. There may
176 * be a more concise way to do this but this should always work. */
177 // std::cout << "Entering block for 3 component arrays"<<std::endl;
178 for (auto k = 0; k < 3; ++k) {
179 // std::cout << "k="<<k<<std::endl;
180 ar &peak_snr[k];
181 // std::cout << "bandwidth_fraction"<<std::endl;
182 ar &signal_bandwidth_fraction[k];
183 }
184 // std::cout << "Exiting serialize function"<<std::endl;
185 }
186};
187} // namespace mspass::algorithms::deconvolution
188#endif
mspass::seismic::PowerSpectrum compute_noise_spectrum(const mspass::seismic::TimeSeries &d2use)
Definition CNRDeconEngine.cc:217
Object to hold components needed in all fft based decon algorithms.
Definition FFTDeconOperator.h:20
Multittaper power spectral estimator.
Definition MTPowerSpectrumEngine.h:29
Frequency domain shaping wavelet.
Definition ShapingWavelet.h:21
Definition PowerSpectrum.h:11
Implemntation of Seismogram for MsPASS.
Definition Seismogram.h:14
Implemntation of TimeSeries for MsPASS.
Definition TimeSeries.h:14
C++ object version of a parameter file.
Definition AntelopePf.h:61
Definition Metadata.h:71