GPLIB++
ResPkModel.cpp
Go to the documentation of this file.
1 #include "ResPkModel.h"
2 #include "FatalException.h"
3 #include <fstream>
4 #include <iomanip>
5 #include <iostream>
6 #include <algorithm>
7 #include <boost/cast.hpp>
8 using namespace std;
9 
10 namespace gplib
11  {
12  ResPkModel::ResPkModel() :
13  slowness(-1.0), InputWave(PWave)
14  {
15  }
17  {
18  }
19 
20  ResPkModel::ResPkModel(const int nlayers) :
21  slowness(-1.0), InputWave(PWave), SeismicModel(nlayers)
22  {
23 
24  }
25 
27  slowness(-1.0), InputWave(PWave)
28  {
29  this->SeismicModel::operator=(source);
30  }
31 
33  {
34  this->SeismicModel::operator=(source);
35  this->slowness = source.slowness;
36  this->InputWave = source.InputWave;
37  copy(source.Strike.begin(), source.Strike.end(), back_inserter(this->Strike));
38  copy(source.Dip.begin(), source.Dip.end(), back_inserter(this->Dip));
39  }
40 
42  {
43  if (this == &source)
44  return *this;
45  this->SeismicModel::operator=(source);
46  this->slowness = source.slowness;
47  this->InputWave = source.InputWave;
48  Strike.assign(source.Strike.size(), 0);
49  Dip.assign(source.Dip.size(), 0);
50  copy(source.Strike.begin(), source.Strike.end(), Strike.begin());
51  copy(source.Dip.begin(), source.Dip.end(), Dip.begin());
52  return *this;
53  }
54 
55  void ResPkModel::WriteModel(const std::string filename)
56  {
57  ofstream neu; //the file to write the model to
58 
59  neu.open(filename.c_str()); //open file for writing
60  const unsigned int nlayers = GetPVelocity().size();
61  if (Dip.size() != nlayers)
62  {
63  Dip.assign(nlayers, 0.0);
64  }
65  if (Strike.size() != nlayers)
66  {
67  Strike.assign(nlayers, 0.0);
68  }
69  neu << setw(3) << nlayers << " " << filename << endl;
70 
71  for (unsigned int i = 0; i < nlayers; ++i) //write all layers
72  {
73  neu.setf(ios::fixed);
74  neu << setw(3);
75  //we write out everything with a precision of 4 digits to make
76  //sure that we represent the values well
77  //only thickness only has a precision of 2 to allow value > 100
78  //also 10m precision is sufficient for layer thickness
79  //respktn is picky about formatting
80  neu << i + 1 << setw(8) << setprecision(4) << GetPVelocity().at(i);
81  neu << setw(8) << setprecision(4) << GetSVelocity().at(i) << setw(8)
82  << setprecision(4) << GetDensity().at(i);
83  neu << setw(8) << setprecision(2) << GetThickness().at(i) << setw(8)
84  << setprecision(4) << GetQp().at(i);
85  neu << setw(8) << setprecision(4) << GetQs().at(i);
86  neu << setw(8) << setprecision(4) << GetStrike().at(i) << setw(8)
87  << setprecision(4) << GetDip().at(i) << " 0.2500" << endl;
88  }
89  neu.close(); //close file
90  }
91 
92  void ResPkModel::ReadModel(const std::string filename)
93  {
94  ifstream infile(filename.c_str());
95  char dummy[255];
96  double number;
97  int nlayers = -1;
98  infile >> nlayers;
99  infile.getline(dummy, 255);
100  if (nlayers < 1)
101  throw FatalException("Problem reading file: " + filename);
102  Init(nlayers);
103  int i = 0;
104  while (infile.good())
105  {
106  //we throw away the layer number
107  infile >> number;
108  //if the layer number read succeeded we read in the rest
109  if (infile.good())
110  {
111  infile >> SetPVelocity().at(i) >> SetSVelocity().at(i)
112  >> SetDensity().at(i) >> SetThickness().at(i) >> SetQp().at(i)
113  >> SetQs().at(i);
114  infile >> SetStrike().at(i) >> SetDip().at(i) >> number;
115  ++i;
116  }
117  }
118  if (i != nlayers)
119  throw FatalException("Problem reading file: " + filename);
120  }
121 
122  void ResPkModel::WriteRunFile(const std::string &filename)
123  {
124  const unsigned int mintime = 200; //the minimum time the seismogram has to have to be correct
125  ofstream runfile;
126  runfile.open(filename.c_str());
127  runfile << "#!/bin/bash" << endl;
128  runfile << "respknt 2>&1> /dev/null << eof" << endl;
129  runfile << filename + ".mod" << endl;
130  runfile << "n" << endl;
131  if (InputWave == PWave)
132  {
133  runfile << "1" << endl;
134  }
135  else
136  {
137  runfile << "2" << endl;
138  }
139  runfile << GetDt() << endl;
140  runfile
141  << max(boost::numeric_cast<unsigned int>((GetNpts() - 1) * GetDt()), mintime)
142  << endl;
143  runfile << GetSlowness() << endl;
144  runfile << "f" << endl;
145  runfile << "y" << endl;
146  runfile << "eof" << endl;
147  runfile.close();
148  }
149  }
virtual ~ResPkModel()
Definition: ResPkModel.cpp:16
virtual void ReadModel(const std::string filename)
Read the model in its native format from a file.
Definition: ResPkModel.cpp:92
virtual void WriteRunFile(const std::string &filename)
Write out a script that performs a forward calculation for the current model.
Definition: ResPkModel.cpp:122
const trealdata & GetStrike() const
Read-only access to the vector of Strikes in each layer.
Definition: ResPkModel.h:46
trealdata & SetDensity()
Read-write access to the vector of densities in g/cm^3 in each layer.
Definition: SeismicModel.h:96
trealdata & SetStrike()
Read-write access to the vector of Strikes in each layer.
Definition: ResPkModel.h:51
const trealdata & GetDensity() const
Read-only access to the vector of densities in g/cm^3 in each layer.
Definition: SeismicModel.h:91
trealdata & SetQs()
Read-write access to the vector of S quality factors for each layer.
Definition: SeismicModel.h:126
const trealdata & GetQp() const
Read-only access to the vector of P quality factors for each layer.
Definition: SeismicModel.h:111
double GetSlowness() const
Get the slowness in s/km for the synthetic forward calculation.
Definition: ResPkModel.h:36
trealdata & SetThickness()
Read-write access to the vector of thicknesses in km in each layer.
Definition: SeismicModel.h:106
trealdata & SetQp()
Read-write access to the vector of P quality factors for each layer.
Definition: SeismicModel.h:116
trealdata & SetDip()
Read-write access to the vector of Dips in each layer.
Definition: ResPkModel.h:61
trealdata & SetPVelocity()
Read-write access to the vector of P-velocities in km/s in each layer.
Definition: SeismicModel.h:76
The class SeismicModel is the base class for some of the model format for seismic codes...
Definition: SeismicModel.h:17
double GetDt() const
Get the time between two samples in s, this is for synthetic forward calculations.
Definition: SeismicModel.h:61
SeismicModel & operator=(const SeismicModel &source)
const trealdata & GetSVelocity() const
Read-only access to the vector of S-velocities in km/s in each layer.
Definition: SeismicModel.h:81
ResPkModel & operator=(const ResPkModel &source)
Definition: ResPkModel.cpp:41
void Init(const int nlayers)
Init provides a convenient way to allocate memory in all structures for a given number of layers...
trealdata & SetSVelocity()
Read-write access to the vector of S-velocities in km/s in each layer.
Definition: SeismicModel.h:86
unsigned int GetNpts() const
Get the number of points for synthetic seismogram calculation.
Definition: SeismicModel.h:41
virtual void WriteModel(const std::string filename)
Write the model in its native format to a file.
Definition: ResPkModel.cpp:55
const trealdata & GetDip() const
Read-only access to the vector of Dips in each layer.
Definition: ResPkModel.h:56
This class stores and writes model for the respktn 1D seismic code that we use for receiver function ...
Definition: ResPkModel.h:18
const trealdata & GetPVelocity() const
Read-only access to the vector of P-velocities in km/s in each layer.
Definition: SeismicModel.h:71
const trealdata & GetThickness() const
Read-only access to the vector of thicknesses in km in each layer.
Definition: SeismicModel.h:101
The basic exception class for all errors that arise in gplib.
const trealdata & GetQs() const
Read-only access to the vector of S quality factors for each layer.
Definition: SeismicModel.h:121