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

Implemntation of TimeSeries for MsPASS. More...

#include <TimeSeries.h>

Inheritance diagram for mspass::seismic::TimeSeries:
mspass::seismic::CoreTimeSeries mspass::utility::ProcessingHistory mspass::seismic::BasicTimeSeries mspass::utility::Metadata mspass::utility::BasicProcessingHistory mspass::utility::BasicMetadata mspass::seismic::TimeSeriesWGaps

Public Member Functions

 TimeSeries ()
 
 TimeSeries (const size_t nsamples)
 
 TimeSeries (const BasicTimeSeries &bts, const Metadata &md)
 
 TimeSeries (const Metadata &md)
 
 TimeSeries (const mspass::seismic::CoreTimeSeries &d)
 Construct from lower level CoreTimeSeries. More...
 
 TimeSeries (const mspass::seismic::CoreTimeSeries &d, const std::string alg)
 
 TimeSeries (const mspass::seismic::BasicTimeSeries &b, const mspass::utility::Metadata &m, const mspass::utility::ProcessingHistory &mcts, const std::vector< double > &d)
 
 TimeSeries (const TimeSeries &parent)
 
TimeSeriesoperator= (const TimeSeries &parent)
 
TimeSeriesoperator+= (const TimeSeries &d)
 
TimeSeriesoperator*= (const double scale)
 
TimeSeriesoperator-= (const TimeSeries &d)
 
void load_history (const mspass::utility::ProcessingHistory &h)
 
size_t memory_use () const
 
- Public Member Functions inherited from mspass::seismic::CoreTimeSeries
 CoreTimeSeries ()
 
 CoreTimeSeries (const size_t nsin)
 
 CoreTimeSeries (const BasicTimeSeries &bts, const Metadata &md)
 
 CoreTimeSeries (const CoreTimeSeries &)
 
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...
 
CoreTimeSeriesoperator= (const CoreTimeSeries &parent)
 
CoreTimeSeriesoperator+= (const CoreTimeSeries &d)
 Summation operator. More...
 
const CoreTimeSeries operator+ (const CoreTimeSeries &other) const
 
CoreTimeSeriesoperator*= (const double)
 
CoreTimeSeriesoperator-= (const CoreTimeSeries &d)
 Subtraction operator. More...
 
const CoreTimeSeries operator- (const CoreTimeSeries &other) const
 
double operator[] (size_t const sample) const
 
- 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::CoreTimeSeries
std::vector< double > s
 
- Public Attributes inherited from mspass::utility::ProcessingHistory
ErrorLogger elog
 
- 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 TimeSeries for MsPASS.

This is the working version of a three-component TimeSeries object used in the MsPASS framework. It extends CoreTimeSeries by adding common MsPASS components (ProcessingHistory). It may evolve with additional special features.

Constructor & Destructor Documentation

◆ TimeSeries() [1/8]

mspass::seismic::TimeSeries::TimeSeries ( )
inline

Default constructor. Only runs subclass default constructors.

20  {};
Scalar time series data object.
Definition: CoreTimeSeries.h:18
Lightweight class to preserve procesing chain of atomic objects.
Definition: ProcessingHistory.h:258

◆ TimeSeries() [2/8]

mspass::seismic::TimeSeries::TimeSeries ( const size_t  nsamples)
inline

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 an nsamples vector and initializes it to all zeros. 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 vector, 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 samples needed for storing the data vector.

◆ TimeSeries() [3/8]

mspass::seismic::TimeSeries::TimeSeries ( const BasicTimeSeries bts,
const Metadata md 
)
inline

Partially construct from components.

There are times one wants to use the Metadata area as a template to flesh out a CoreTimeSeries as what might be called skin and bones: skin is Metadata and bones as BasicTimeSeries data. This constructor initializes those two base classes but does not fully a valid data vector. It only attempts to fetch the number of points expected for the data vector using the npts metadata (integer) key (i.e. it sets npts to md.get_int("npts")). It then creates the data vector of that length and initialzies it to all zeros.

This constructor is largely a wrapper for the CoreTimeSeries version with the same signature but it also initalizes the ProcessingHistory to null (calling the default constructor)

◆ TimeSeries() [4/8]

mspass::seismic::TimeSeries::TimeSeries ( const Metadata md)

Partially construct from Metadata alone.

This constructor is useful for interaction with MongoDB where the Metadata container is constructed directly from MongoDB documents the database uses for storage. A key restrition of this constructor is that BasicTimeSeries attributes and the size of the internal array buffer, which is set by npts, are extracted from Metadata with keys fixed in the C++ code. The following keys are required or this contructor will throw a MwPASSWrror:

dt - sample sample_interval t0 - data startttime npts - number of samples (s array will be initialized to this many zeros) time_standard - UTC or Relative (anything but UTC is taken as relative)

It will also handle but not require the two attributes used in the mspass schema to handle shifting from absolute to relative time. These keys are t0_shift - sets amount t0 has been shiftee when originally utc but set to relative with the shift method. Ignored if time is UTC or it is not defined.

30 {
31  mlive=false;
32  try {
33  this->Metadata::operator=(md);
34  /* Names used are from mspass defintions as of Jan 2020.
35  We don't need to call the set methods for these attributes as they
36  would add the overhead of setting delta, startime, and npts to the
37  same value passed. */
38  this->mdt = this->get_double(SEISMICMD_dt);
39  this->mt0 = this->get_double(SEISMICMD_t0);
41  {
42  string tstd = this->get_string(SEISMICMD_time_standard);
43  if(tstd == "UTC")
45  else if(tstd == "Relative")
47  else
48  {
50  this->elog.log_error("TimeSeries Metadata constructor",
51  SEISMICMD_time_standard+" attribute is not defined - set to Relative",
52  ErrorSeverity::Complaint);
53  }
54  }
55  if(this->time_is_relative())
56  {
57  /* It is not an error if a t0 shift is not defined and we are
58  in relative time. That is the norm for active source data. */
60  {
61  double t0shift=this->get_double(SEISMICMD_t0_shift);
62  this->force_t0_shift(t0shift);
63  }
64  }
65  /* this default construct is needed to handle miniseed data in
66  MsPASS. It perhaps should generate a log message a information
67  but for now we do this silently. since the constructor returns
68  a result marked dead in all cases a default of 0 is sensible.*/
69  long int ns;
71  {
72  ns = md.get_long(SEISMICMD_npts);
73  }
74  else
75  {
76  ns = 0;
77  }
78  /* this CoreTimeSeries method sets the npts attribute and
79  initializes the s buffer to all zeros */
80  this->set_npts(ns);
81  }catch(...) {throw;};
82 }
bool mlive
Definition: BasicTimeSeries.h:256
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
double mdt
Definition: BasicTimeSeries.h:260
double mt0
Definition: BasicTimeSeries.h:264
void set_npts(const size_t npts)
Set the number of samples attribute for data.
Definition: CoreTimeSeries.cc:222
int log_error(const mspass::utility::MsPASSError &merr)
Definition: ErrorLogger.cc:81
bool is_defined(const std::string key) const noexcept
Definition: Metadata.cc:73
long get_long(const std::string key) const
Definition: Metadata.h:184
std::string get_string(const std::string key) const override
Definition: Metadata.h:213
Metadata & operator=(const Metadata &mdold)
Definition: Metadata.cc:108
double get_double(const std::string key) const override
Definition: Metadata.h:135
ProcessingHistory()
Definition: ProcessingHistory.cc:86
const std::string SEISMICMD_t0_shift("starttime_shift")
const std::string SEISMICMD_t0("starttime")
const std::string SEISMICMD_dt("delta")
const std::string SEISMICMD_time_standard("time_standard")
const std::string SEISMICMD_npts("npts")

References mspass::seismic::BasicTimeSeries::force_t0_shift(), mspass::utility::Metadata::get_double(), mspass::utility::Metadata::get_long(), mspass::utility::Metadata::get_string(), mspass::utility::Metadata::is_defined(), mspass::utility::ErrorLogger::log_error(), mspass::seismic::BasicTimeSeries::mdt, mspass::seismic::BasicTimeSeries::mlive, mspass::seismic::BasicTimeSeries::mt0, mspass::utility::Metadata::operator=(), mspass::seismic::Relative, mspass::seismic::SEISMICMD_dt(), mspass::seismic::SEISMICMD_npts(), mspass::seismic::SEISMICMD_t0(), mspass::seismic::SEISMICMD_t0_shift(), mspass::seismic::SEISMICMD_time_standard(), mspass::seismic::CoreTimeSeries::set_npts(), mspass::seismic::BasicTimeSeries::set_tref(), and mspass::seismic::UTC.

◆ TimeSeries() [5/8]

mspass::seismic::TimeSeries::TimeSeries ( const mspass::seismic::CoreTimeSeries d)
inline

Construct from lower level CoreTimeSeries.

In MsPASS CoreTimeSeries has the primary functions that define the concept of a a single channel seismogram. TimeSeries implements mspass specific features needed to mesh with the mspass processing system. This constructor clones only the CoreTimeSeries components and initializes the ProcessingHistory with the default constructor leaving 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 TimeSeries object. Note that means a deep copy wherein the data vector is copied.

◆ TimeSeries() [6/8]

mspass::seismic::TimeSeries::TimeSeries ( const mspass::seismic::CoreTimeSeries d,
const std::string  alg 
)

Contruct from a core time series and initialize history as origin.

This constructor is a variant of a similar one built only from a CoreTimeSeries. 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.
12 {
13  /* Not sure this is a good idea, but will give each instance
14  created by this constructor a uuid.*/
15  string id=this->newid();
16  this->ProcessingHistory::set_as_origin(alg,id,id,AtomicType::SEISMOGRAM,false);
17  this->ProcessingHistory::set_jobname(string("test"));
18  this->ProcessingHistory::set_jobid(string("test"));
19 }
CoreTimeSeries()
Definition: CoreTimeSeries.cc:15
void set_as_origin(const std::string alg, const std::string algid, const std::string uuid, const AtomicType typ, bool define_as_raw=false)
Definition: ProcessingHistory.cc:161
std::string newid()
Definition: ProcessingHistory.cc:700

References mspass::utility::ProcessingHistory::newid(), and mspass::utility::ProcessingHistory::set_as_origin().

◆ TimeSeries() [7/8]

mspass::seismic::TimeSeries::TimeSeries ( const mspass::seismic::BasicTimeSeries b,
const mspass::utility::Metadata m,
const mspass::utility::ProcessingHistory mcts,
const std::vector< double > &  d 
)

Special constructor for pickle interface.

The pickle interface required by spark presented problems for MsPASS. The complicated data objects of TimeSeries and TimeSeries have to be serialized in pieces. This constructor is only used in the function called indirectly by pickle.load. It essentially contains a TimeSeries dismembered into the pieces that can be the serialized independently. The parameters are each associated with one of those required pieces and are simply copied to build a valid TimeSeries object in the pickle.load function

◆ TimeSeries() [8/8]

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

Standard copy constructor.

Member Function Documentation

◆ memory_use()

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

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

Memory consumed by a TimeSeries 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.

97 {
98  size_t memory_estimate;
99  memory_estimate = sizeof(TimeSeries);
100  memory_estimate += sizeof(double)*this->npts();
101  /* We can only estimate the size of the Metadata container.
102  These constants are defined in memory_constants.h */
103  memory_estimate += memory_constants::MD_AVERAGE_SIZE*this->md.size();
104  memory_estimate += memory_constants::KEY_AVERAGE_SIZE*this->changed_or_set.size();
105  /* Similar for history and elog containers */
106  memory_estimate += memory_constants::HISTORYDATA_AVERAGE_SIZE*this->nodes.size();
107  memory_estimate += memory_constants::ELOG_AVERAGE_SIZE*this->elog.size();
108  return memory_estimate;
109 }
size_t npts() const
Definition: BasicTimeSeries.h:183
TimeSeries()
Definition: TimeSeries.h:18

References mspass::seismic::BasicTimeSeries::npts(), and TimeSeries().

◆ operator=()

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

Standard assignment operator.

84 {
85  if(this!=(&parent))
86  {
87  this->CoreTimeSeries::operator=(parent);
88  this->ProcessingHistory::operator=(parent);
89  }
90  return *this;
91 }
CoreTimeSeries & operator=(const CoreTimeSeries &parent)
Definition: CoreTimeSeries.cc:66
ProcessingHistory & operator=(const ProcessingHistory &parent)
Definition: ProcessingHistory.cc:741

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


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