GPLIB++
RFVelCalc.cpp
Go to the documentation of this file.
1 #include "RFVelCalc.h"
2 #include "WFunc.h"
3 #include <numeric>
4 #include <fstream>
5 #include <cassert>
6 #include "FatalException.h"
7 #include "NumUtil.h"
8 
9 namespace gplib
10  {
11  RFVelCalc::RFVelCalc(const double mysigma, const double myc,
12  const RecCalc::trfmethod themethod) :
13  //we initialize the receiver function calculation with a shift of 5 seconds
14  //this enables us catch the full waveform of the initial correlation peak in most cases
15  //there might be a better way to do this
16  RFCalculator(5, mysigma, myc, false, themethod)
17  {
18  //with normalization absolute velocity estimation does not work
19  //so we make sure the RF is not normalized
20  RFCalculator.SetNormalize(false);
21  }
22 
24  {
25  }
26 
28  RFCalculator(Old.RFCalculator), Velocities(Old.Velocities), Periods(
29  Old.Periods)
30  {
31 
32  }
33 
35  {
36  if (this == &source)
37  return *this;
38  RFCalculator = source.RFCalculator;
39  Velocities = source.Velocities;
40  Periods = source.Periods;
41  return *this;
42  }
43 
44  void RFVelCalc::AbsVelCalc(const double slowness,
45  const SeismicDataComp &RadRec, const SeismicDataComp &VerRec,
46  ttsdata &AppVel)
47  {
48  //we need filtered versions of the two components
49  //so create copies that we can work on
50  ttsdata FilteredR(RadRec.GetData().size()), FilteredZ(
51  VerRec.GetData().size());
52  //make sure all old stuff is deleted before we do another calculation
53  AppVel.clear();
54  Velocities.clear();
55  Periods.clear();
56  //in debug mode we perform some extra checks
57  assert(FilteredR.size() == FilteredZ.size());
58  assert(fcmp(RadRec.GetDt(),VerRec.GetDt(),std::numeric_limits<double>::epsilon())==0);
59  //we have to compensate for the shift of 5 seconds that we use above
60  const double relshift = 5.0
61  / (RadRec.GetDt() * RadRec.GetData().size());
62  //at the moment we calculate apparent velocities at periods independent of the receiver function parameters
63  //we specify the width in seconds
64  for (double width = 0.1; width < 80; width *= 1.1)
65  {
66  //calculate the width of the filter relative to the length of the whole component
67  double sigma = width / (RadRec.GetDt() * RadRec.GetData().size());
68  //window the component with a truncated squared cosine window
69  ApplyWindow(RadRec.GetData().begin(), RadRec.GetData().end(),
70  FilteredR.begin(), TruncCosSq(sigma), relshift);
71  ApplyWindow(VerRec.GetData().begin(), VerRec.GetData().end(),
72  FilteredZ.begin(), TruncCosSq(sigma), relshift);
73  // and sum up the values, this is eq. 5 in Svenningsen
74  double rvalue = std::accumulate(FilteredR.begin(), FilteredR.end(), 0.0);
75  double zvalue = std::accumulate(FilteredZ.begin(), FilteredZ.end(), 0.0);
76 
77  if (zvalue == 0.0)
78  throw FatalException("Zero amplitude for Z-Receiver function !");
79  //calculate the apparent velocity
80  double vsapp = sin(0.5 * atan2(rvalue, zvalue)) / slowness;
81  //store the apparent velocity values and the periods
82  AppVel.push_back(vsapp);
83  Velocities.push_back(vsapp);
84  Periods.push_back(width);
85  }
86  }
87 
88  void RFVelCalc::CalcRFVelFromRec(const double slowness,
89  const SeismicDataComp &RRec, const SeismicDataComp &VerComp,
90  ttsdata &AppVel)
91  {
92  SeismicDataComp ZRecFunc;
93  RFCalculator.CalcRecData(VerComp, VerComp, ZRecFunc);
94  AbsVelCalc(slowness, RRec, ZRecFunc, AppVel);
95  }
96 
97  void RFVelCalc::CalcRFVel(const double slowness,
98  const SeismicDataComp &RadComp, const SeismicDataComp &VerComp,
99  ttsdata &AppVel)
100  {
101  SeismicDataComp RRecFunc, ZRecFunc;
102  RFCalculator.CalcRecData(RadComp, VerComp, RRecFunc);
103  RFCalculator.CalcRecData(VerComp, VerComp, ZRecFunc);
104  AbsVelCalc(slowness, RRecFunc, ZRecFunc, AppVel);
105  }
106 
107  void RFVelCalc::WriteVelocities(const std::string filename)
108  {
109  std::ofstream outfile(filename.c_str());
110  for (size_t i = 0; i < Velocities.size(); ++i)
111  outfile << Periods.at(i) << " " << Velocities.at(i) << std::endl;
112  }
113  }
trfmethod
There are several ways to calculate receiver functions.
Definition: RecCalc.h:18
std::vector< double > & GetData()
Access for data vector, for ease of use and efficiency we return a reference.
void WriteVelocities(const std::string filename)
Write the velocity estiamtes and corresponding periods into an ascii file.
Definition: RFVelCalc.cpp:107
void SetNormalize(const bool what)
Change whether the output receiver function is normalized to a maximum amplitude of 1...
Definition: RecCalc.h:55
void CalcRecData(const SeismicDataComp &RadComp, const SeismicDataComp &VerComp, SeismicDataComp &Receiver)
Calculate Receiver functions from two data components.
Definition: RecCalc.cpp:173
void ApplyWindow(InputIterator inbegin, InputIterator inend, OutputIterator outbegin, WindowFunctype WFunc, double relshift=0.0)
Apply one of the above window functions to a range.
Definition: WFunc.h:104
RFVelCalc(const double mysigma, const double myc, const RecCalc::trfmethod themethod=RecCalc::specdiv)
The constructor takes three parameters for the receiver function estimation.
Definition: RFVelCalc.cpp:11
This class implements the method to calculate absolute S-Wave velocities from Receiver function data ...
Definition: RFVelCalc.h:13
void CalcRFVel(const double slowness, const SeismicDataComp &RadComp, const SeismicDataComp &VerComp, ttsdata &AppVel)
Calculate absolute velocities from the radial and vertical components of the seismogram.
Definition: RFVelCalc.cpp:97
RFVelCalc & operator=(const RFVelCalc &source)
Definition: RFVelCalc.cpp:34
A variable width cosine squared window that is zero outside.
Definition: WFunc.h:77
double GetDt() const
Return dt in s.
void CalcRFVelFromRec(const double slowness, const SeismicDataComp &RRec, const SeismicDataComp &VerComp, ttsdata &AppVel)
Calculate absolute velocities from the radial receiver function and the vertical component of the sei...
Definition: RFVelCalc.cpp:88
virtual ~RFVelCalc()
Definition: RFVelCalc.cpp:23
The basic exception class for all errors that arise in gplib.