version  0.0.1
Defines the C++ API for MsPASS
Public Member Functions | List of all members
mspass::seismic::Seismogram Class Reference

Implemntation of Seismogram for MsPASS. More...

#include <Seismogram.h>

Inheritance diagram for mspass::seismic::Seismogram:
mspass::seismic::CoreSeismogram mspass::utility::ProcessingHistory mspass::seismic::BasicTimeSeries mspass::utility::Metadata mspass::utility::BasicProcessingHistory mspass::utility::BasicMetadata mspass::seismic::SeismogramWGaps

Public Member Functions

 Seismogram ()
 
 Seismogram (const size_t nsamples)
 
 Seismogram (const mspass::seismic::CoreSeismogram &d)
 Construct from lower level CoreSeismogram. More...
 
 Seismogram (const mspass::seismic::CoreSeismogram &d, const std::string alg)
 
 Seismogram (const mspass::seismic::BasicTimeSeries &bts, const mspass::utility::Metadata &md)
 Partial constructor from base classes. More...
 
 Seismogram (const mspass::seismic::BasicTimeSeries &b, const mspass::utility::Metadata &m, const mspass::utility::ProcessingHistory &his, const bool card, const bool ortho, const mspass::utility::dmatrix &tm, const mspass::utility::dmatrix &uin)
 Construct from all pieces. More...
 
 Seismogram (const Metadata &md, const std::string jobname=std::string("test"), const std::string jobid=std::string("UNDEFINED"), const std::string readername=std::string("load3C"), const std::string algid=std::string("0"))
 
 Seismogram (const Metadata &md, bool load_data)
 Construct from Metadata definition that includes data path. More...
 
 Seismogram (const Seismogram &parent)
 
Seismogramoperator= (const Seismogram &parent)
 
void load_history (const mspass::utility::ProcessingHistory &h)
 Load just the ProcessingHistory data from another data source. More...
 
size_t memory_use () const
 
- Public Member Functions inherited from mspass::seismic::CoreSeismogram
 CoreSeismogram ()
 
 CoreSeismogram (const size_t nsamples)
 
 CoreSeismogram (const std::vector< mspass::seismic::CoreTimeSeries > &ts, const unsigned int component_to_clone=0)
 
 CoreSeismogram (const mspass::utility::Metadata &md, const bool load_data=true)
 Construct from Metadata definition that includes data path. More...
 
 CoreSeismogram (const CoreSeismogram &)
 
void set_dt (const double sample_interval)
 Set the sample interval. More...
 
void set_npts (const size_t npts)
 Set the number of samples attribute for data. More...
 
void sync_npts ()
 Sync the number of samples attribute with actual data size. More...
 
void set_t0 (const double t0in)
 Set the data start time. More...
 
CoreSeismogramoperator= (const CoreSeismogram &)
 
CoreSeismogramoperator+= (const CoreSeismogram &d)
 Summation operator. More...
 
const CoreSeismogram operator+ (const CoreSeismogram &other) const
 
CoreSeismogramoperator*= (const double)
 
CoreSeismogramoperator-= (const CoreSeismogram &d)
 Subtraction operator. More...
 
const CoreSeismogram operator- (const CoreSeismogram &other) const
 
std::vector< double > operator[] (const int sample) const
 
std::vector< double > operator[] (const double time) const
 Overloaded version of operator[] for time. More...
 
virtual ~CoreSeismogram ()
 
void rotate_to_standard ()
 
void rotate (mspass::utility::SphericalCoordinate &sc)
 
void rotate (const double nu[3])
 
void rotate (const double phi)
 Rotate horizontals by a simple angle in degrees. More...
 
void transform (const double a[3][3])
 
void free_surface_transformation (const mspass::seismic::SlownessVector u, const double vp0, const double vs0)
 
mspass::utility::dmatrix get_transformation_matrix () const
 
bool set_transformation_matrix (const mspass::utility::dmatrix &A)
 Define the transformaton matrix. More...
 
bool set_transformation_matrix (const double a[3][3])
 Define the transformaton matrix with a C style 3x3 matrix. More...
 
bool set_transformation_matrix (pybind11::object a)
 Define the transformaton matrix with a python object. More...
 
bool cardinal () const
 
bool orthogonal () const
 
double endtime () const noexcept
 
- Public Member Functions inherited from mspass::seismic::BasicTimeSeries
 BasicTimeSeries ()
 
 BasicTimeSeries (const BasicTimeSeries &)
 
virtual ~BasicTimeSeries ()
 Virtual destructor. More...
 
double time (const int i) const
 
int sample_number (double t) const
 
double endtime () const noexcept
 
bool shifted () const
 
double time_reference () const
 
void force_t0_shift (const double t)
 Force a t0 shift value on data. More...
 
virtual void ator (const double tshift)
 
virtual void rtoa ()
 
virtual void shift (const double dt)
 
bool live () const
 
bool dead () const
 
void kill ()
 
void set_live ()
 
double dt () const
 
bool time_is_UTC () const
 
bool time_is_relative () const
 
TimeReferenceType timetype () const
 
double samprate () const
 
size_t npts () const
 
double t0 () const
 
void set_tref (const TimeReferenceType newtref)
 Force the time standard. More...
 
BasicTimeSeriesoperator= (const BasicTimeSeries &parent)
 
- Public Member Functions inherited from mspass::utility::ProcessingHistory
 ProcessingHistory ()
 
 ProcessingHistory (const std::string jobnm, const std::string jid)
 
 ProcessingHistory (const ProcessingHistory &parent)
 
bool is_empty () const
 
bool is_raw () const
 
bool is_origin () const
 
bool is_volatile () const
 
bool is_saved () const
 
size_t number_of_stages () override
 Return number of processing stages that have been applied to this object. More...
 
void set_as_origin (const std::string alg, const std::string algid, const std::string uuid, const AtomicType typ, bool define_as_raw=false)
 
std::string new_ensemble_process (const std::string alg, const std::string algid, const AtomicType typ, const std::vector< ProcessingHistory * > parents, const bool create_newid=true)
 
void add_one_input (const ProcessingHistory &data_to_add)
 Add one datum as an input for current data. More...
 
void add_many_inputs (const std::vector< ProcessingHistory * > &d)
 Define several data objects as inputs. More...
 
void merge (const ProcessingHistory &data_to_add)
 Merge the history nodes from another. More...
 
void accumulate (const std::string alg, const std::string algid, const AtomicType typ, const ProcessingHistory &newinput)
 Method to use with a spark reduce algorithm. More...
 
std::string clean_accumulate_uuids ()
 Clean up inconsistent uuids that can be produced by reduce. More...
 
std::string new_map (const std::string alg, const std::string algid, const AtomicType typ, const ProcessingStatus newstatus=ProcessingStatus::VOLATILE)
 Define this algorithm as a one-to-one map of same type data. More...
 
std::string new_map (const std::string alg, const std::string algid, const AtomicType typ, const ProcessingHistory &data_to_clone, const ProcessingStatus newstatus=ProcessingStatus::VOLATILE)
 Define this algorithm as a one-to-one map. More...
 
std::string map_as_saved (const std::string alg, const std::string algid, const AtomicType typ)
 Prepare the current data for saving. More...
 
void clear ()
 
std::multimap< std::string, mspass::utility::NodeDataget_nodes () const
 
int stage () const
 
ProcessingStatus status () const
 
std::string id () const
 
std::pair< std::string, std::string > created_by () const
 
NodeData current_nodedata () const
 
std::string newid ()
 
int number_inputs () const
 
int number_inputs (const std::string uuidstr) const
 
void set_id (const std::string newid)
 
std::list< mspass::utility::NodeDatainputs (const std::string id_to_find) const
 Return a list of data that define the inputs to a give uuids. More...
 
ProcessingHistoryoperator= (const ProcessingHistory &parent)
 
- Public Member Functions inherited from mspass::utility::BasicProcessingHistory
 BasicProcessingHistory (const std::string jobname, const std::string jobid)
 
 BasicProcessingHistory (const BasicProcessingHistory &parent)
 
std::string jobid () const
 
void set_jobid (const std::string &newjid)
 
std::string jobname () const
 
void set_jobname (const std::string jobname)
 
BasicProcessingHistoryoperator= (const BasicProcessingHistory &parent)
 

Additional Inherited Members

- Public Attributes inherited from mspass::seismic::CoreSeismogram
mspass::utility::dmatrix u
 
- Public Attributes inherited from mspass::utility::ProcessingHistory
ErrorLogger elog
 
- Protected Attributes inherited from mspass::seismic::CoreSeismogram
bool components_are_orthogonal
 
bool components_are_cardinal
 
double tmatrix [3][3]
 
- Protected Attributes inherited from mspass::seismic::BasicTimeSeries
bool mlive
 
double mdt
 
double mt0
 
size_t nsamp
 
TimeReferenceType tref
 
bool t0shift_is_valid
 
double t0shift
 
- Protected Attributes inherited from mspass::utility::ProcessingHistory
std::multimap< std::string, mspass::utility::NodeDatanodes
 
- Protected Attributes inherited from mspass::utility::BasicProcessingHistory
std::string jid
 
std::string jnm
 

Detailed Description

Implemntation of Seismogram for MsPASS.

This is the working version of a three-component seismogram object used in the MsPASS framework. It extends CoreSeismogram by adding ProcessingHistory.

Constructor & Destructor Documentation

◆ Seismogram() [1/9]

mspass::seismic::Seismogram::Seismogram ( )
inline

Default constructor. Only runs subclass default constructors.

Vector (three-component) seismogram data object.
Definition: CoreSeismogram.h:40
Lightweight class to preserve procesing chain of atomic objects.
Definition: ProcessingHistory.h:258

◆ Seismogram() [2/9]

mspass::seismic::Seismogram::Seismogram ( const size_t  nsamples)

Bare bones constructor allocates space and little else.

Sometimes it is helpful to construct a skeleton that can be fleshed out manually. This constructor allocates a 3xnsamples array and sets the matrix to all zeros. It also sets the orientation information to cardinal and defines the data in relative time units. The data are marked dead because the assumption is the caller will fill out commonly needed basic Metadata, load some kind of sample data into the data matrix, and then call the set_live method when that process is completed. That kind of manipulation is most common in preparing simulation or test data where the common tags on real data do not exist and need to be defined manually for the simulatioon or test. The history section is initialized with the default constructor, which currently means it is empty. If a simulation or test requires a history origin the user must load it manaually.

Parameters
nsamplesis the number of 3c vector samples to iniitalize the object with (creates a 3xnsamples matrix).
10  : CoreSeismogram(nsamples),ProcessingHistory()
11 {
12 /* Note this constructor body needs no content. Just a wrapper for CoreSeismogram */
13 }
CoreSeismogram()
Definition: CoreSeismogram.cc:28
ProcessingHistory()
Definition: ProcessingHistory.cc:86

◆ Seismogram() [3/9]

mspass::seismic::Seismogram::Seismogram ( const mspass::seismic::CoreSeismogram d)

Construct from lower level CoreSeismogram.

In MsPASS CoreSeismogram has the primary functions that define the concept of a three-component seismogram. Seismogram implements mspass specific features needed to mesh with the mspass processing system. This constructor clones onlyh the CoreSeismogram components and initializes the ProcessingHistory with the default constructor leaves the history in an empty state. Users should call methods in ProcessingHistory to initiate a valid history chain. Processing can continue if left in that state, but the history chain will have an undefined origin. Job information will also be lost if not initialized (see BasicProcessingHistory)

Parameters
dis the data to be copied to create the new Seismogram object
18 {
19 }

◆ Seismogram() [4/9]

mspass::seismic::Seismogram::Seismogram ( const mspass::seismic::CoreSeismogram d,
const std::string  alg 
)

Contruct from a core seismogram and initialize history as origin.

This constructor is a variant of a similar one built only from a CoreSeismogram. This constuctor is intended to be used mainly on simulation data created by some mechanism (e.g. a python procedure). It is more rigid that the simple one arg constructor as it will create a top level history record. The record will mark the result as an origin with an id set as a uuid created by a random number generator (boost implementation). The jobname and jobid are frozen as "test". Use a different constructor and/or reset job info if more flexibility is desired. Use of this constructor is recommended only for test python programs that do not need to interact with MongoDB.

Parameters
iscore data to be cloned
algis the algorithm name to set for the origin history record.

◆ Seismogram() [5/9]

mspass::seismic::Seismogram::Seismogram ( const mspass::seismic::BasicTimeSeries bts,
const mspass::utility::Metadata md 
)

Partial constructor from base classes.

This is a partial constructor useful sometimes for building a Seismogram object from scratch such as in a simultion. It clones the Metadata and uses attributes in BasicTimeSeries to build a framework that data can be added into. It initialzies the data array to npts stored in the BasicTimeSeries passed to the constructor.

Parameters
btsis the BasicTimeSeries that overrides any md data.
mdis the Metadata cloned to make the partially constructed object.
31  : CoreSeismogram(md,false),ProcessingHistory()
32 {
33  /* the contents of BasicTimeSeries passed will override anything set from
34  Metadata in this section. Note also the very important use of this
35  for the putters to assure proper resolution of the virtual methods */
36  this->kill(); //set dead because buffer has invalid data
37  this->set_t0(bts.t0());
38  this->set_tref(bts.timetype());
39  this->set_npts(bts.npts());
40  /* The handling of these is kind of awkward in the current api. Good to
41  hide it here. Note in current implementation force_t0_shift sets the
42  shift as valid (what shifted tests). I am slightly worried there are
43  cases where the else block could create an inconsistency with active
44  source data if this method were used to convert raw data. Careful if
45  this is used in that context where the time reference is defined
46  implicitly as shot time with no tie to an absolute time reference */
47  if(bts.shifted())
48  {
49  double t0shift;
50  t0shift=bts.time_reference();
51  this->force_t0_shift(t0shift);
52  }
53  else
54  {
55  this->force_t0_shift(0.0);
56  }
57 }
bool shifted() const
Definition: BasicTimeSeries.h:88
size_t npts() const
Definition: BasicTimeSeries.h:183
double t0() const
Definition: BasicTimeSeries.h:186
double time_reference() const
Definition: BasicTimeSeries.cc:86
void force_t0_shift(const double t)
Force a t0 shift value on data.
Definition: BasicTimeSeries.h:115
void set_tref(const TimeReferenceType newtref)
Force the time standard.
Definition: BasicTimeSeries.h:243
void kill()
Definition: BasicTimeSeries.h:151
void set_t0(const double t0in)
Set the data start time.
Definition: CoreSeismogram.cc:1047
void set_npts(const size_t npts)
Set the number of samples attribute for data.
Definition: CoreSeismogram.cc:1067

References mspass::seismic::BasicTimeSeries::force_t0_shift(), mspass::seismic::BasicTimeSeries::kill(), mspass::seismic::BasicTimeSeries::npts(), mspass::seismic::CoreSeismogram::set_npts(), mspass::seismic::CoreSeismogram::set_t0(), mspass::seismic::BasicTimeSeries::set_tref(), mspass::seismic::BasicTimeSeries::shifted(), mspass::seismic::BasicTimeSeries::t0(), and mspass::seismic::BasicTimeSeries::time_reference().

◆ Seismogram() [6/9]

mspass::seismic::Seismogram::Seismogram ( const mspass::seismic::BasicTimeSeries b,
const mspass::utility::Metadata m,
const mspass::utility::ProcessingHistory his,
const bool  card,
const bool  ortho,
const mspass::utility::dmatrix tm,
const mspass::utility::dmatrix uin 
)

Construct from all pieces.

This constructor build a Seismogram object from all the pieces that define this highest level object in mspass. This constructor is planned to be hidden from python programmers and not exposed with pybind11 wrappers because is has some potentially undesirable side effects if not used carefully. The primary purpose of this constructor is for serialization and deserializaton in spark with pickle. The pickle interface is purely in C for this function so again python programmers don't need to see this constructor. The parameter names are obvious because they are associated one for one with the objects they are used to construct. The only exceptions are card and ortho which are booleans used to set the internal booleans components_are_cardinal and components_are_orthogonal respectively (two internals users should not mess with). tm is also a dmatrix representation the tmatrix stored internally as a 2d C array, but we use the dmatrix to mesh with serialization.

69  : ProcessingHistory(his)
70 {
72  this->Metadata::operator=(m);
75  int i,j;
76  for(i=0;i<3;++i)
77  for(j=0;j<3;++j) tmatrix[i][j]=tm(i,j);
78  this->u=uin;
79 }
BasicTimeSeries & operator=(const BasicTimeSeries &parent)
Definition: BasicTimeSeries.cc:62
double tmatrix[3][3]
Definition: CoreSeismogram.h:519
bool components_are_orthogonal
Definition: CoreSeismogram.h:489
mspass::utility::dmatrix u
Definition: CoreSeismogram.h:53
bool components_are_cardinal
Definition: CoreSeismogram.h:511
Metadata & operator=(const Metadata &mdold)
Definition: Metadata.cc:108

References mspass::seismic::CoreSeismogram::components_are_cardinal, mspass::seismic::CoreSeismogram::components_are_orthogonal, mspass::seismic::BasicTimeSeries::operator=(), mspass::utility::Metadata::operator=(), mspass::seismic::CoreSeismogram::tmatrix, and mspass::seismic::CoreSeismogram::u.

◆ Seismogram() [7/9]

mspass::seismic::Seismogram::Seismogram ( const Metadata md,
const std::string  jobname = std::string("test"),
const std::string  jobid = std::string("UNDEFINED"),
const std::string  readername = std::string("load3C"),
const std::string  algid = std::string("0") 
)

Constructor driven by a Metadata object.

The flexibilityof Metadata makes it helpful at times to build a Seismogram object driven by Metadata key:value pairs. This implementation always uses file io to read the actual data and clone the Metadata. ProcessingHistory is initialized to define the result as raw data. The id field is set from Metadata using the special (frozen) key "wf_id" that is assumed to hold the string representation of an ObjectID of the waveform read.

Note all but the Metadata argument are defaulted, but setting the jobname and jobid is strongly recommended to avoid ambiguity - these define a unique run of a particular workflow.

Parameters
mdis the Metadata that is used to drive the constructor.
jobnameis the jobname (defaults to test, but should be changed)
readernameis the algorithm name assigned to the top level history record. Defaults to "load3C"

◆ Seismogram() [8/9]

mspass::seismic::Seismogram::Seismogram ( const Metadata md,
bool  load_data 
)
inline

Construct from Metadata definition that includes data path.

A Metadata object is sufficiently general that it can contain enough information to contruct an object from attributes contained in it. This constuctor uses that approach, with the actual loading of data being an option (on by default). In mspass this is constructor is used to load data with Metadata constructed from MongoDB and then using the path created from two parameters (dir and dfile used as in css3.0 wfdisc) to read data. The API is general but the implementation in mspass is very rigid. It blindly assumes the data being read are binary doubles in the right byte order and ordered in the native order for dmatrix (Fortran order). i.e. the constuctor does a raw fread of ns*3 doubles into the internal array used in the dmatrix implementation.

A second element of the Metadata that is special for MsPASS is the handling of the transformation matrix by this constructor. In MsPASS the transformation matrix is stored as a python object in MongoDB. This constructor aims to fetch that entity with the key 'tmatrix'. To be more robust and simpler to use with data not loaded from mongodb we default tmatrix to assume the data are in standard coordinates. That is, if the key tmatrix is not defined in Metadata passed as arg0, the constructor assumes it should set the transformation matrix to an identity. Use set_transformation_matrix if that assumption is wrong for your data.

Parameters
mdis the Metadata used for the construction.
load_dataif true (default) a file name is constructed from dir+"/"+dfile, the file is openned, fseek is called to foff, data are read with fread, and the file is closed. If false a dmatrix for u is still created of size 3xns, but the matrix is only initialized to all zeros.
Exceptions
Willthrow a MsPASSError if required metadata are missing.

◆ Seismogram() [9/9]

mspass::seismic::Seismogram::Seismogram ( const Seismogram parent)
inline

Standard copy constructor.

Member Function Documentation

◆ load_history()

void mspass::seismic::Seismogram::load_history ( const mspass::utility::ProcessingHistory h)

Load just the ProcessingHistory data from another data source.

Some algorithms don't handle processing history. In those situations it can prove helpful to manage ProcessingHistory separately and load the history data through this method.

Parameters
his the ProcessingHistory data to copy into this Seismogram.
125 {
127 }
ProcessingHistory & operator=(const ProcessingHistory &parent)
Definition: ProcessingHistory.cc:741

References mspass::utility::ProcessingHistory::operator=().

◆ memory_use()

size_t mspass::seismic::Seismogram::memory_use ( ) const

Return an estimate of the memmory use by the data in this object.

Memory consumed by a Seismogram object is needed to implement the sizeof method in python that dask/spark use to manage memory. Without that feature we had memory fault issues. Note the estimate this method returns should not be expected to be exact. The MsPASS implementation or any alternative implementation avoids an exact calculation because it requries an (expensive) traversal of multiple map containers.

129 {
130  size_t memory_estimate;
131  memory_estimate = sizeof(Seismogram);
132  /* data for a seismogram is 3 channels so 3*npts*/
133  memory_estimate += sizeof(double)*3*this->npts();
134  /* We can only estimate the size of the Metadata container.
135  These constants are defined in memory_constants.h */
136  memory_estimate += memory_constants::MD_AVERAGE_SIZE*this->md.size();
137  memory_estimate += memory_constants::KEY_AVERAGE_SIZE*this->changed_or_set.size();
138  /* Similar for history and elog containers */
139  memory_estimate += memory_constants::HISTORYDATA_AVERAGE_SIZE*this->nodes.size();
140  memory_estimate += memory_constants::ELOG_AVERAGE_SIZE*this->elog.size();
141  return memory_estimate;
142 }
Seismogram()
Definition: Seismogram.h:18

References mspass::seismic::BasicTimeSeries::npts(), Seismogram(), and mspass::utility::Metadata::size().

◆ operator=()

Seismogram & mspass::seismic::Seismogram::operator= ( const Seismogram parent)

Standard assignment operator.

116 {
117  if(this!=(&parent))
118  {
119  this->CoreSeismogram::operator=(parent);
120  this->ProcessingHistory::operator=(parent);
121  }
122  return *this;
123 }
CoreSeismogram & operator=(const CoreSeismogram &)
Definition: CoreSeismogram.cc:896

References mspass::seismic::CoreSeismogram::operator=(), and mspass::utility::ProcessingHistory::operator=().


The documentation for this class was generated from the following files: