version  0.0.1
Defines the C++ API for MsPASS
Loading...
Searching...
No Matches
CoreSeismogram.h
1#ifndef _MSPASS_CORESEISMOGRAM_H_
2#define _MSPASS_CORESEISMOGRAM_H_
3#include "mspass/seismic/BasicTimeSeries.h"
4#include "mspass/seismic/CoreTimeSeries.h"
5#include "mspass/seismic/SlownessVector.h"
6#include "mspass/utility/Metadata.h"
7#include "mspass/utility/SphericalCoordinate.h"
8#include "mspass/utility/dmatrix.h"
9#include <memory>
10#include <vector>
11// #include "mspass/seismic/Ensemble.h"
12namespace mspass::seismic {
13
14/* A Seismogram is viewed as a special collection of Time Series
15type data that is essentially a special version of a vector time series.
16It is "special" as the vector is 3D with real components. One could produce a
17similar inherited type for an n vector time series object.
18
19The structure of a vector time series allows the data to be stored in a matrix.
20Here we use a lightweight matrix object I call dmatrix. This object is
21contains core concepts that define a seismogram. It can be extended as in
22MsPASS to add functionality or aliased to Seismogram to simplify the naming as
23done, for example, with std::basic_string made equivalent to std::string.
24*/
25
40public:
53
69 CoreSeismogram(const size_t nsamples);
105 CoreSeismogram(const std::vector<mspass::seismic::CoreTimeSeries> &ts,
106 const unsigned int component_to_clone = 0);
143 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;
289 std::vector<double> operator[](const int sample) const;
301 std::vector<double> operator[](const double time) const;
303 virtual ~CoreSeismogram() {};
319 void rotate_to_standard();
320 // This overloaded pair do the same thing for a vector
321 // specified as a unit vector nu or as spherical coordinate angles
363
390 void rotate(const double nu[3]);
407 void rotate(const double phi);
416 void transform(const double a[3][3]);
430 const double vp0, const double vs0);
441 mspass::utility::dmatrix result(3, 3);
442 for (int i = 0; i < 3; ++i)
443 for (int j = 0; j < 3; ++j)
444 result(i, j) = tmatrix[i][j];
445 return result;
446 };
475 bool set_transformation_matrix(const double a[3][3]);
487 bool set_transformation_matrix(pybind11::object a);
488
490 bool cardinal() const { return components_are_cardinal; };
492 bool orthogonal() const { return components_are_orthogonal; };
497 double endtime() const noexcept {
498 return (mt0 + mdt * static_cast<double>(u.columns() - 1));
499 };
500
501protected:
521 bool components_are_cardinal; // true if x1=e, x2=n, x3=up
529 double tmatrix[3][3];
530
531private:
532 /* This is used internally when necessary to test if a computed or input
533 matrix is an identity matrix. */
534 bool tmatrix_is_cardinal();
535};
536} // namespace mspass::seismic
537#endif // End guard
Base class for time series objects.
Definition BasicTimeSeries.h:35
size_t npts() const
Definition BasicTimeSeries.h:171
double time(const int i) const
Definition BasicTimeSeries.h:62
double mdt
Definition BasicTimeSeries.h:247
double mt0
Definition BasicTimeSeries.h:251
Vector (three-component) seismogram data object.
Definition CoreSeismogram.h:39
const CoreSeismogram operator+(const CoreSeismogram &other) const
Definition CoreSeismogram.cc:995
std::vector< double > operator[](const int sample) const
Definition CoreSeismogram.cc:1086
bool cardinal() const
Definition CoreSeismogram.h:490
virtual ~CoreSeismogram()
Definition CoreSeismogram.h:303
double tmatrix[3][3]
Definition CoreSeismogram.h:529
bool orthogonal() const
Definition CoreSeismogram.h:492
CoreSeismogram & operator+=(const CoreSeismogram &d)
Summation operator.
Definition CoreSeismogram.cc:914
void set_dt(const double sample_interval)
Set the sample interval.
Definition CoreSeismogram.cc:1006
void sync_npts()
Sync the number of samples attribute with actual data size.
Definition CoreSeismogram.cc:1065
bool set_transformation_matrix(pybind11::object a)
Define the transformaton matrix with a python object.
void rotate_to_standard()
Definition CoreSeismogram.cc:380
void set_t0(const double t0in)
Set the data start time.
Definition CoreSeismogram.cc:1022
bool components_are_orthogonal
Definition CoreSeismogram.h:512
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:775
mspass::utility::dmatrix u
Definition CoreSeismogram.h:52
CoreSeismogram()
Definition CoreSeismogram.cc:26
CoreSeismogram & operator=(const CoreSeismogram &)
Definition CoreSeismogram.cc:888
mspass::utility::dmatrix get_transformation_matrix() const
Definition CoreSeismogram.h:440
void rotate(mspass::utility::SphericalCoordinate &sc)
Definition CoreSeismogram.cc:524
const CoreSeismogram operator-(const CoreSeismogram &other) const
Definition CoreSeismogram.cc:1001
CoreSeismogram & operator*=(const double)
Definition CoreSeismogram.cc:903
void free_surface_transformation(const mspass::seismic::SlownessVector u, const double vp0, const double vs0)
Definition CoreSeismogram.cc:706
double endtime() const noexcept
Definition CoreSeismogram.h:497
void set_npts(const size_t npts)
Set the number of samples attribute for data.
Definition CoreSeismogram.cc:1039
void transform(const double a[3][3])
Definition CoreSeismogram.cc:653
CoreSeismogram & operator-=(const CoreSeismogram &d)
Subtraction operator.
Definition CoreSeismogram.cc:955
bool components_are_cardinal
Definition CoreSeismogram.h:521
Slowness vector object.
Definition SlownessVector.h:16
Definition Metadata.h:71
Lightweight, simple matrix object.
Definition dmatrix.h:104
size_t columns() const
Definition dmatrix.cc:216
Define metadata keys.
Definition BasicSpectrum.h:6
Encapsulates spherical coordinates in a data structure.
Definition SphericalCoordinate.h:15