version  0.0.1
Defines the C++ API for MsPASS
Loading...
Searching...
No Matches
GeneralIterDecon.h
1#ifndef __SIIMPLE_GENERAL_ITER_DECON__
2#define __SIIMPLE_GENERAL_ITER_DECON__
3#include "mspass/algorithms/TimeWindow.h"
4#include "mspass/algorithms/deconvolution/ComplexArray.h"
5#include "mspass/algorithms/deconvolution/FFTDeconOperator.h"
6#include "mspass/algorithms/deconvolution/ScalarDecon.h"
7#include "mspass/seismic/CoreTimeSeries.h"
8#include "mspass/seismic/Seismogram.h"
9#include "mspass/utility/AntelopePf.h"
10#include "mspass/utility/Metadata.h"
11#include "mspass/utility/dmatrix.h"
12#include <vector>
13namespace mspass::algorithms::deconvolution {
14/* Wang's original version allowed XCORR here - dropped for now as this
15code assumes the decon produces a zero phase ideal output while in every
16case the original XCORR iterative method uses the raw wavelet and cross
17correlation. */
18
19enum IterDeconType { WATER_LEVEL, LEAST_SQ, MULTI_TAPER };
21public:
24 int col;
26 double u[3];
28 double amp;
29 /* This is the primary constructor for oneo of thise. In the iterative
30 sequnce of generalized iterative decon we find the column with maximum
31 amplitude and point this constructor at the right column. The
32 parent dmatrix is then updated subtracting the inverse wavelet at
33 this lag. */
37 ThreeCSpike(const ThreeCSpike &parent);
40 ThreeCSpike &operator=(const ThreeCSpike &parent);
41};
67public:
77 void changeparameter(const mspass::utility::Metadata &md) {
78 this->preprocessor->changeparameter(md);
79 };
80 int load(const mspass::seismic::CoreSeismogram &d,
82 int loadnoise(const mspass::seismic::CoreSeismogram &d,
89 int load(const mspass::seismic::CoreSeismogram &d,
92 void process();
95 mspass::seismic::CoreTimeSeries ideal_output() {
96 return this->preprocessor->ideal_output();
97 };
99 return this->preprocessor->actual_output();
100 };
102 return this->preprocessor->inverse_wavelet();
103 };
105 return this->preprocessor->inverse_wavelet(t0parent);
106 };
108
109private:
110 /* These are data at different stages of process. d_all is the
111 largest signal window that is assumed to have been initialized by the
112 load method for this object. d_decon is the
113 result of apply the preprocessor (signal processing) deconvolution to
114 all three components. d_decon is computed from d, which is a windowed
115 version of the input data received by the load method. It have a
116 time duration less than or at least equal to that of d_all.
117 r is the residual, which is accumulated during the iterative method.
118 The time duraction of r is the same as d. It is initalized by
119 convolving the inverse filter with d_all.
120 n is the noise data. It should normally be at least as long a d_all*/
121 mspass::seismic::CoreSeismogram d_all, d_decon, r, n;
122 /* We save the set of data lengths for clarity and a minor bit efficiency.
123 ndwin is the number of samples in d_all and r.
124 nnwin is the number of samples in n.
125 nfft is the size of d_decon */
126 int ndwin, nnwin, nfft;
127 /* this is the time shift in sample applied by preprocessor*/
128 int time_shift;
129 /* Save the TimeWindow objects hat define the extent of d_all, d_decon,
130 and n. Some things need at least some of these downstream */
131 mspass::algorithms::TimeWindow dwin, nwin, fftwin;
134 int noise_component;
138 std::list<ThreeCSpike> spikes;
145 std::vector<double> lag_weights;
146 /* This vector contains the function time shifted and added to lag_weights
147 vector after each iteration. */
148 std::vector<double> wtf;
149 int nwtf; // size of wtf - cached because wtf is inside the deepest loop
150
151 /* This is a pointer to the BasicDeconOperator class used for preprocessing
152 Classic use of inheritance to simplify the api. */
153 ScalarDecon *preprocessor;
154
155 /* This defines a shaping wavelet applied to output. The inverse filter
156 preprocessing algorithm is assumed to have it's own embedded shaping wavelet.
157 i.e. this is the wavelet convolved with the output spike sequence. */
158 ShapingWavelet shapingwavelet;
159 /* This parameter is set in the constructor. It would normally be half the
160 length of the fir representation of the inverse wavelet. (winv_fir below)*/
161 int wavelet_pad;
170 std::vector<double> winv_fir, actual_o_fir;
171 int actual_o_0; // offset from sample zero for zero lag position
172 IterDeconType decon_type;
173 /* This is called by the constructor to create the wtf penalty function */
174 void construct_weight_penalty_function(const mspass::utility::Metadata &md);
180 void update_residual_matrix(ThreeCSpike spk);
187 void update_lag_weights(int col);
193 double compute_resid_linf_floor();
196 bool has_not_converged();
197 /* These are convergence attributes. lw_inf indicates Linf norm of
198 lag_weight array, lw_l2 is L2 metric of lag_weight, resid_inf is Linf
199 norm of residual vector, and resid_l2 is L2 of resid matrix. prev
200 modifier means the size of that quantity in the previous iteration. initial
201 means initial value at the top of the loop.*/
202 double lw_linf_initial, lw_linf_prev;
203 double lw_l2_initial, lw_l2_prev;
204 double resid_linf_initial, resid_linf_prev;
205 double resid_l2_initial, resid_l2_prev;
206 /* These are convergence paramwters for the different tests */
207 int iter_count, iter_max; // actual iteration count and ceiling to break loop
208 /*lw metrics are scaled with range of 0 to 1. l2 gets scaled by number of
209 points and so can use a similar absolute scale. */
210 double lw_linf_floor, lw_l2_floor;
211 /* We use a probability level to define the floor Linf of the residual
212 matrix. For L2 we use the conventional fractional improvment metric. */
213 double resid_linf_prob, resid_linf_floor;
214 double resid_l2_tol;
215
216 /* DEBUG attributes and methods. These should be deleted after testing. */
217 std::vector<double> lw_linf_history, lw_l2_history, resid_l2_history,
218 resid_linf_history;
219 void print_convergence_history(std::ostream &ofs) {
220 ofs << "lw_inf lw_l2 resid_l2, resid_linf" << std::endl;
221 for (int i = 0; i < iter_count; ++i) {
222 if (i < lw_linf_history.size()) {
223 ofs << lw_linf_history[i] << " ";
224 } else {
225 ofs << "xxxx ";
226 }
227 if (i < lw_l2_history.size()) {
228 ofs << lw_l2_history[i] << " ";
229 } else {
230 ofs << "xxxx ";
231 }
232 if (i < resid_linf_history.size()) {
233 ofs << resid_linf_history[i] << " ";
234 } else {
235 ofs << "xxxx ";
236 }
237 if (i < resid_l2_history.size()) {
238 ofs << resid_l2_history[i] << " ";
239 } else {
240 ofs << "xxxx ";
241 }
242 ofs << std::endl;
243 }
244 }
245};
246} // namespace mspass::algorithms::deconvolution
247#endif
Defines a time window.
Definition TimeWindow.h:12
Implements the Generalized Iterative Method of Wang and Pavlis.
Definition GeneralIterDecon.h:66
mspass::seismic::CoreTimeSeries inverse_wavelet()
Return a FIR represention of the inverse filter.
Definition GeneralIterDecon.h:101
mspass::seismic::CoreTimeSeries actual_output()
Return the actual output of the deconvolution operator.
Definition GeneralIterDecon.h:98
mspass::utility::Metadata QCMetrics()
Return appropriate quality measures.
Definition GeneralIterDecon.cc:738
Base class decon operator for single station 3C decon (receiver functions).
Definition ScalarDecon.h:29
virtual mspass::seismic::CoreTimeSeries inverse_wavelet()=0
Return a FIR represention of the inverse filter.
virtual mspass::seismic::CoreTimeSeries actual_output()=0
Return the actual output of the deconvolution operator.
Frequency domain shaping wavelet.
Definition ShapingWavelet.h:21
Definition GeneralIterDecon.h:20
double amp
Definition GeneralIterDecon.h:28
double u[3]
Definition GeneralIterDecon.h:26
ThreeCSpike & operator=(const ThreeCSpike &parent)
Definition ThreeCSpike.cc:27
int col
Definition GeneralIterDecon.h:24
Vector (three-component) seismogram data object.
Definition CoreSeismogram.h:39
Scalar time series data object.
Definition CoreTimeSeries.h:17
C++ object version of a parameter file.
Definition AntelopePf.h:61
Definition Metadata.h:71
Lightweight, simple matrix object.
Definition dmatrix.h:104