version  0.0.1
Defines the C++ API for MsPASS
CoreSeismogram.h
1 #ifndef _MSPASS_CORESEISMOGRAM_H_
2 #define _MSPASS_CORESEISMOGRAM_H_
3 #include <memory>
4 #include <vector>
5 #include "mspass/utility/Metadata.h"
6 #include "mspass/utility/dmatrix.h"
7 #include "mspass/seismic/BasicTimeSeries.h"
8 #include "mspass/seismic/CoreTimeSeries.h"
9 #include "mspass/utility/SphericalCoordinate.h"
10 #include "mspass/seismic/SlownessVector.h"
11 //#include "mspass/seismic/Ensemble.h"
12 namespace mspass::seismic{
13 
14 /* A Seismogram is viewed as a special collection of Time Series
15 type data that is essentially a special version of a vector time series.
16 It is "special" as the vector is 3D with real components. One could produce a
17 similar inherited type for an n vector time series object.
18 
19 The structure of a vector time series allows the data to be stored in a matrix.
20 Here we use a lightweight matrix object I call dmatrix. This object is
21 contains core concepts that define a seismogram. It can be extended as in
22 MsPASS to add functionality or aliased to Seismogram to simplify the naming as
23 done, for example, with std::basic_string made equivalent to std::string.
24 */
25 
26 
40 {
41 public:
54 
70  CoreSeismogram(const size_t nsamples);
106  CoreSeismogram(const std::vector<mspass::seismic::CoreTimeSeries>& ts,
107  const unsigned int component_to_clone=0);
143  CoreSeismogram(const mspass::utility::Metadata& md,const bool load_data=true);
149  /* These overload virtual methods in BasicTimeSeries. */
162  void set_dt(const double sample_interval);
184  void set_npts(const size_t npts);
191  void sync_npts();
208  void set_t0(const double t0in);
237  const CoreSeismogram operator+(const CoreSeismogram& other) const;
239  CoreSeismogram& operator*=(const double);
264  const CoreSeismogram operator-(const CoreSeismogram& other) const;
287  std::vector<double> operator[](const int sample)const;
299  std::vector<double> operator[](const double time)const;
301  virtual ~CoreSeismogram(){};
315  void rotate_to_standard();
316  // This overloaded pair do the same thing for a vector
317  // specified as a unit vector nu or as spherical coordinate angles
357 
382  void rotate(const double nu[3]);
399  void rotate(const double phi);
408  void transform(const double a[3][3]);
421  void free_surface_transformation(const mspass::seismic::SlownessVector u, const double vp0, const double vs0);
432  {
433  mspass::utility::dmatrix result(3,3);
434  for(int i=0;i<3;++i)
435  for(int j=0;j<3;++j) result(i,j)=tmatrix[i][j];
436  return result;
437  };
466  bool set_transformation_matrix(const double a[3][3]);
476  bool set_transformation_matrix(pybind11::object a);
477 
479  bool cardinal()const {return components_are_cardinal;};
481  bool orthogonal()const {return components_are_orthogonal;};
486  double endtime()const noexcept
487  {
488  return(mt0+mdt*static_cast<double>(u.columns()-1));
489  };
490 
491 protected:
511  bool components_are_cardinal; // true if x1=e, x2=n, x3=up
519  double tmatrix[3][3];
520 private:
521  /* This is used internally when necessary to test if a computed or input
522  matrix is an identity matrix. */
523  bool tmatrix_is_cardinal();
524 };
525 } //end mspass::seismic namespace enscapsulation
526 #endif // End guard
Base class for time series objects.
Definition: BasicTimeSeries.h:35
size_t npts() const
Definition: BasicTimeSeries.h:183
double time(const int i) const
Definition: BasicTimeSeries.h:63
double mdt
Definition: BasicTimeSeries.h:260
double mt0
Definition: BasicTimeSeries.h:264
Vector (three-component) seismogram data object.
Definition: CoreSeismogram.h:40
const CoreSeismogram operator+(const CoreSeismogram &other) const
Definition: CoreSeismogram.cc:1016
std::vector< double > operator[](const int sample) const
Definition: CoreSeismogram.cc:1123
bool cardinal() const
Definition: CoreSeismogram.h:479
virtual ~CoreSeismogram()
Definition: CoreSeismogram.h:301
double tmatrix[3][3]
Definition: CoreSeismogram.h:519
bool orthogonal() const
Definition: CoreSeismogram.h:481
CoreSeismogram & operator+=(const CoreSeismogram &d)
Summation operator.
Definition: CoreSeismogram.cc:926
void set_dt(const double sample_interval)
Set the sample interval.
Definition: CoreSeismogram.cc:1028
void sync_npts()
Sync the number of samples attribute with actual data size.
Definition: CoreSeismogram.cc:1099
bool set_transformation_matrix(pybind11::object a)
Define the transformaton matrix with a python object.
void rotate_to_standard()
Definition: CoreSeismogram.cc:400
void set_t0(const double t0in)
Set the data start time.
Definition: CoreSeismogram.cc:1047
bool components_are_orthogonal
Definition: CoreSeismogram.h:489
CoreSeismogram(const std::vector< mspass::seismic::CoreTimeSeries > &ts, const unsigned int component_to_clone=0)
bool set_transformation_matrix(const mspass::utility::dmatrix &A)
Define the transformaton matrix.
Definition: CoreSeismogram.cc:787
mspass::utility::dmatrix u
Definition: CoreSeismogram.h:53
CoreSeismogram()
Definition: CoreSeismogram.cc:28
CoreSeismogram & operator=(const CoreSeismogram &)
Definition: CoreSeismogram.cc:896
mspass::utility::dmatrix get_transformation_matrix() const
Definition: CoreSeismogram.h:431
void rotate(mspass::utility::SphericalCoordinate &sc)
Definition: CoreSeismogram.cc:544
const CoreSeismogram operator-(const CoreSeismogram &other) const
Definition: CoreSeismogram.cc:1022
CoreSeismogram & operator*=(const double)
Definition: CoreSeismogram.cc:915
void free_surface_transformation(const mspass::seismic::SlownessVector u, const double vp0, const double vs0)
Definition: CoreSeismogram.cc:723
double endtime() const noexcept
Definition: CoreSeismogram.h:486
void set_npts(const size_t npts)
Set the number of samples attribute for data.
Definition: CoreSeismogram.cc:1067
void transform(const double a[3][3])
Definition: CoreSeismogram.cc:673
CoreSeismogram & operator-=(const CoreSeismogram &d)
Subtraction operator.
Definition: CoreSeismogram.cc:972
bool components_are_cardinal
Definition: CoreSeismogram.h:511
Slowness vector object.
Definition: SlownessVector.h:18
Definition: Metadata.h:76
Lightweight, simple matrix object.
Definition: dmatrix.h:102
size_t columns() const
Definition: dmatrix.cc:232
Define metadata keys.
Definition: BasicSpectrum.h:6
Encapsulates spherical coordinates in a data structure.
Definition: SphericalCoordinate.h:14