GPLIB++
RecCalc.h
Go to the documentation of this file.
1 #ifndef CRECFUNC_H
2 #define CRECFUNC_H
3 #include "SeismicDataComp.h"
4 #include "ResPkModel.h"
5 #include "TsSpectrum.h"
6 #include <string>
7 
8 namespace gplib
9  {
10  /** \addtogroup seistools Seismic data analysis and modeling */
11  /* @{ */
12 
13  //! This class is used to calculate receiver functions from seismic data
14  class RecCalc
15  {
16  public:
17  //! There are several ways to calculate receiver functions
18  enum trfmethod
19  {
21  };
22  private:
23  //! The method we want to use to calculate the receiver function
24  trfmethod method;
25  //! For spectral division we need an object to calculate spectra
26  TsSpectrum Spectrum;
27  //! The waterlevel to fill in spectral holes
28  double c;
29  //! The width of the gaussian filter to smooth the receiver functions
30  double sigma;
31  //! Shift the data to move the initial correlation peak
32  int shift;
33  //! Should the receiver function be normalized to maximum amplitude
34  bool normalize;
35  //! The components used for the receiver function calculation
36  SeismicDataComp RadComp;
37  SeismicDataComp VerComp;
38  //! Calculate the receiver function with the spectral division waterlevel method
39  void SpectralDivision(const SeismicDataComp &RComp,
40  const SeismicDataComp &VComp, SeismicDataComp &Receiver);
41  //! Calculate the receiver function by iterative convolution
42  void IterativeDeconvolution(const SeismicDataComp &RComp,
43  const SeismicDataComp &VComp, SeismicDataComp &Receiver);
44  public:
45  //! Get the radial component, mostly needed for synthetic data
47  {
48  return RadComp;
49  }
51  {
52  return VerComp;
53  }
54  //! Change whether the output receiver function is normalized to a maximum amplitude of 1
55  void SetNormalize(const bool what)
56  {
57  normalize = what;
58  }
59  //! The three Synth*Parallel methods provide alternative acces to the steps in CalcRecSynth for safe parallel execution
60  void SynthPreParallel(const std::string &filename, ResPkModel &Model,
61  SeismicDataComp &Receiver, const bool cleanfiles = false);
62  //! All operations that are safe to execute in parallel
63  void SynthSafeParallel(const std::string &filename, ResPkModel &Model,
64  SeismicDataComp &Receiver, const bool cleanfiles = false);
65  //! Operations of the synthetic receiver function calculation that are not safe in parallel and hafe to be executed after the parallel part
66  void SynthPostParallel(const std::string &filename, ResPkModel &Model,
67  SeismicDataComp &Receiver, const bool cleanfiles = false);
68  //! Calculate Receiver functions from two data components
69  void CalcRecData(const SeismicDataComp &RadComp,
70  const SeismicDataComp &VerComp, SeismicDataComp &Receiver);
71  //! Calculate synthetic receiver funtions from a model
72  void CalcRecSynth(const std::string &filename, ResPkModel &Model,
73  SeismicDataComp &Receiver, const bool cleanfiles = false);
74  RecCalc& operator=(const RecCalc& source);
75  //! The constructor takes the essential parameters that shouldn't change during different calculations
76  RecCalc(int myshift, double mysigma, double myc,
77  bool multicalc = false, trfmethod themethod = iterdecon);
78  RecCalc(const RecCalc &Old);
79  virtual ~RecCalc();
80  };
81  /* @} */
82  }
83 #endif // CRECFUNC_H
trfmethod
There are several ways to calculate receiver functions.
Definition: RecCalc.h:18
virtual ~RecCalc()
Definition: RecCalc.cpp:58
This class is used to calculate receiver functions from seismic data.
Definition: RecCalc.h:14
void SynthSafeParallel(const std::string &filename, ResPkModel &Model, SeismicDataComp &Receiver, const bool cleanfiles=false)
All operations that are safe to execute in parallel.
Definition: RecCalc.cpp:210
void CalcRecSynth(const std::string &filename, ResPkModel &Model, SeismicDataComp &Receiver, const bool cleanfiles=false)
Calculate synthetic receiver funtions from a model.
Definition: RecCalc.cpp:195
void SetNormalize(const bool what)
Change whether the output receiver function is normalized to a maximum amplitude of 1...
Definition: RecCalc.h:55
void SynthPostParallel(const std::string &filename, ResPkModel &Model, SeismicDataComp &Receiver, const bool cleanfiles=false)
Operations of the synthetic receiver function calculation that are not safe in parallel and hafe to b...
Definition: RecCalc.cpp:219
void CalcRecData(const SeismicDataComp &RadComp, const SeismicDataComp &VerComp, SeismicDataComp &Receiver)
Calculate Receiver functions from two data components.
Definition: RecCalc.cpp:173
void SynthPreParallel(const std::string &filename, ResPkModel &Model, SeismicDataComp &Receiver, const bool cleanfiles=false)
The three Synth*Parallel methods provide alternative acces to the steps in CalcRecSynth for safe para...
Definition: RecCalc.cpp:203
The class CTsSpectrum is used to calculate spectra from (real) time series data.
Definition: TsSpectrum.h:21
RecCalc & operator=(const RecCalc &source)
Definition: RecCalc.cpp:43
const SeismicDataComp & GetVerComp()
Definition: RecCalc.h:50
RecCalc(int myshift, double mysigma, double myc, bool multicalc=false, trfmethod themethod=iterdecon)
The constructor takes the essential parameters that shouldn't change during different calculations...
Definition: RecCalc.cpp:27
const SeismicDataComp & GetRadComp()
Get the radial component, mostly needed for synthetic data.
Definition: RecCalc.h:46
This class stores and writes model for the respktn 1D seismic code that we use for receiver function ...
Definition: ResPkModel.h:18