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

Scalar time series data object. More...

#include <CoreTimeSeries.h>

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

Public Member Functions

 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 Attributes

std::vector< double > s
 

Additional Inherited Members

- Protected Attributes inherited from mspass::seismic::BasicTimeSeries
bool mlive
 
double mdt
 
double mt0
 
size_t nsamp
 
TimeReferenceType tref
 
bool t0shift_is_valid
 
double t0shift
 

Detailed Description

Scalar time series data object.

This data object extends BasicTimeSeries mainly by adding a vector of scalar data. It uses a Metadata object to contain auxiliary parameters that aren't essential to define the data object, but which are necessary for some algorithms.

Author
Gary L. Pavlis

Constructor & Destructor Documentation

◆ CoreTimeSeries() [1/4]

mspass::seismic::CoreTimeSeries::CoreTimeSeries ( )

Default constructor. Initializes object data to zeros and sets the initial STL vector size to 0 length.

16 {
17  s.reserve(0);
18  this->set_dt(1.0);
19  this->set_t0(0.0);
20  this->set_npts(0);
21 }
BasicTimeSeries()
Definition: BasicTimeSeries.cc:42
void set_dt(const double sample_interval)
Set the sample interval.
Definition: CoreTimeSeries.cc:185
std::vector< double > s
Definition: CoreTimeSeries.h:28
void set_npts(const size_t npts)
Set the number of samples attribute for data.
Definition: CoreTimeSeries.cc:222
void set_t0(const double t0in)
Set the data start time.
Definition: CoreTimeSeries.cc:203
Metadata()
Definition: Metadata.h:80

References s, set_dt(), set_npts(), and set_t0().

◆ CoreTimeSeries() [2/4]

mspass::seismic::CoreTimeSeries::CoreTimeSeries ( const size_t  nsin)

Similar to the default constructor but creates a vector of data with nsin samples and initializes all samples to 0.0. This vector can safely be accessed with the vector index operator (i.e. operator []). A corollary is that push_back or push_front applied to this vector will alter it's length so use this only if the size of the data to fill the object is already known.

23 {
24  /* IMPORTANT: this constructor assumes BasicTimeSeries initializes the
25  equivalent of:
26  set_dt(1.0)
27  set_t0(0.0)
28  set_tref(TimeReferenceType::Relative)
29  this->kill() - i.e. marked dead
30 
31  It further assumes current api where set_npts allocates and initializes s
32  to nsin zeros */
33  this->CoreTimeSeries::set_npts(nsin);
34 }

References set_npts().

◆ CoreTimeSeries() [3/4]

mspass::seismic::CoreTimeSeries::CoreTimeSeries ( const BasicTimeSeries bts,
const Metadata md 
)

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.

58  : BasicTimeSeries(bd), Metadata(md)
59 {
60  /* this assumes set_npts initializes the vector containers, s, to zeros
61  AND that BasicTimeSeries constructor initializes ns (npts) to the value
62  desired. */
63  this->CoreTimeSeries::set_npts(this->nsamp);
64 }
size_t nsamp
Definition: BasicTimeSeries.h:268

References mspass::seismic::BasicTimeSeries::nsamp, and set_npts().

◆ CoreTimeSeries() [4/4]

mspass::seismic::CoreTimeSeries::CoreTimeSeries ( const CoreTimeSeries tsi)

Standard copy constructor.

37  :
38  BasicTimeSeries(tsi),
39  Metadata(tsi)
40 {
41  if(mlive)
42  {
43  s=tsi.s;
44  }
45  else if(tsi.s.size()>0)
46  {
47  /* This is needed to preserve the contents of data vector when something
48  marks the data dead, but one wants to restore it later. Classic example
49  is an interactive trace editor. Found mysterious errors can occur
50  without this features. */
51  s=tsi.s;
52  }
53  /* Do nothing if the parent s is empty as std::vector will be properly
54  initialized*/
55 }
bool mlive
Definition: BasicTimeSeries.h:256

References mspass::seismic::BasicTimeSeries::mlive, and s.

Member Function Documentation

◆ operator*=()

CoreTimeSeries & mspass::seismic::CoreTimeSeries::operator*= ( const double  scale)

Multiply data by a scalar.

181 {
182  dscal(this->npts(),scale,&(this->s[0]),1);
183  return *this;
184 }
size_t npts() const
Definition: BasicTimeSeries.h:183

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

◆ operator+()

const CoreTimeSeries mspass::seismic::CoreTimeSeries::operator+ ( const CoreTimeSeries other) const

Addition operator.

This operator is implemented in a standard way utilizing operator+=. For data with irregular start and end times that has an important consequence; the operator is not communative. i.e given x an y z=x+y will not yield the same result as z=y+x.

169 {
170  CoreTimeSeries result(*this);
171  result += other;
172  return result;
173 }
CoreTimeSeries()
Definition: CoreTimeSeries.cc:15

◆ operator+=()

CoreTimeSeries & mspass::seismic::CoreTimeSeries::operator+= ( const CoreTimeSeries d)

Summation operator.

Summing data from signals of irregular length requires handling potential mismatches in size and overlap. This behaves the way a += operator should logically behave in that situation. That is, because the lhs is where the sum is being accumulated, the size is always controlled by the left hand side of the operator. Any portions of the right hand side that are outside the t0 to endtime() of the left hand side are silently discarded. If the start time of the right hand side is greater than t0 or the endtime is less than endtime of the lhs there will be discontinuties in the sum there the ends of the rhs are inside the range of the lhs.

Parameters
dis other signal to add to this.
Exceptions
MsPASSErrorcan be thrown if lhs and rhs do not have matching time standards.
79 {
80  int i,iend,jend;
81  size_t i0;
82  size_t j,j0;
83  // Sun's compiler complains about const objects without this.
84  CoreTimeSeries& d=const_cast<CoreTimeSeries&>(data);
85  // Silently do nothing if d is marked dead
86  if(!d.mlive) return(*this);
87  // Silently do nothing if d does not overlap with data to contain sum
88  if( (d.endtime()<mt0)
89  || (d.mt0>(this->endtime())) ) return(*this);
90  if(d.tref!=(this->tref))
91  throw MsPASSError("CoreTimeSeries += operator cannot handle data with inconsistent time base\n",
92  ErrorSeverity::Invalid);
93  /* this defines the range of left and right hand sides to be summed */
94  i=d.sample_number(this->mt0);
95  if(i<0)
96  {
97  j0=this->sample_number(d.t0());
98  i0=0;
99  }
100  else
101  {
102  j0=0;
103  i0=i;
104  }
105  iend=d.sample_number(this->endtime());
106  jend=this->sample_number(d.endtime());
107  if(iend>=(d.npts()))
108  {
109  iend=d.npts()-1;
110  }
111  if(jend>=this->npts())
112  {
113  jend=this->npts()-1;
114  }
115  //cout << "i0="<<i0<<" j0="<<j0<<" iend="<<iend<<" jend="<<jend<<endl;
116  /* Now do the actual sum using the computed ranges */
117  for(i=i0,j=j0; i<=iend && j<=jend; ++i,++j)
118  this->s[j]+=d.s[i];
119  return(*this);
120 }
int sample_number(double t) const
Definition: BasicTimeSeries.h:74
double mt0
Definition: BasicTimeSeries.h:264
double endtime() const noexcept
Definition: BasicTimeSeries.h:82
Base class for error object thrown by MsPASS library routines.
Definition: MsPASSError.h:40

References mspass::seismic::BasicTimeSeries::endtime(), mspass::seismic::BasicTimeSeries::mlive, mspass::seismic::BasicTimeSeries::mt0, mspass::seismic::BasicTimeSeries::npts(), s, mspass::seismic::BasicTimeSeries::sample_number(), mspass::seismic::BasicTimeSeries::t0(), and mspass::seismic::BasicTimeSeries::tref.

◆ operator-()

const CoreTimeSeries mspass::seismic::CoreTimeSeries::operator- ( const CoreTimeSeries other) const

Subtraction operator.

This operator is implemented in a standard way utilizing operator-=. For data with irregular start and end times that has an important consequence; the operator is not communative. i.e given x an y z=x-y will not yield the same result as z=-(y-x).

175 {
176  CoreTimeSeries result(*this);
177  result -= other;
178  return result;
179 }

◆ operator-=()

CoreTimeSeries & mspass::seismic::CoreTimeSeries::operator-= ( const CoreTimeSeries d)

Subtraction operator.

Differencing data from signals of irregular length requires handling potential mismatches in size and overlap. This behaves the way a -= operator should logically behave in that situation. That is, because the lhs is where the sum is being accumulated, the size is always controlled by the left hand side of the operator. Any portions of the right hand side that are outside the t0 to endtime() of the left hand side are silently discarded. If the start time of the right hand side is greater than t0 or the endtime is less than endtime of the lhs there will be discontinuties in the sum there the ends of the rhs are inside the range of the lhs.

Parameters
dis other signal to subract from this.
Exceptions
MsPASSErrorcan be thrown if lhs and rhs do not have matching time standards.
126 {
127  int i,iend,jend;
128  size_t i0;
129  size_t j,j0;
130  // Sun's compiler complains about const objects without this.
131  CoreTimeSeries& d=const_cast<CoreTimeSeries&>(data);
132  // Silently do nothing if d is marked dead
133  if(!d.mlive) return(*this);
134  // Silently do nothing if d does not overlap with data to contain sum
135  if( (d.endtime()<mt0)
136  || (d.mt0>(this->endtime())) ) return(*this);
137  if(d.tref!=(this->tref))
138  throw MsPASSError("CoreTimeSeries += operator cannot handle data with inconsistent time base\n",
139  ErrorSeverity::Invalid);
140  /* this defines the range of left and right hand sides to be summed */
141  i=d.sample_number(this->mt0);
142  if(i<0)
143  {
144  j0=this->sample_number(d.t0());
145  i0=0;
146  }
147  else
148  {
149  j0=0;
150  i0=i;
151  }
152  iend=d.sample_number(this->endtime());
153  jend=this->sample_number(d.endtime());
154  if(iend>=(d.npts()))
155  {
156  iend=d.npts()-1;
157  }
158  if(jend>=this->npts())
159  {
160  jend=this->npts()-1;
161  }
162  //cout << "i0="<<i0<<" j0="<<j0<<" iend="<<iend<<" jend="<<jend<<endl;
163  /* Now do the actual sum using the computed ranges */
164  for(i=i0,j=j0; i<=iend && j<=jend; ++i,++j)
165  this->s[j]-=d.s[i];
166  return(*this);
167 }

References mspass::seismic::BasicTimeSeries::endtime(), mspass::seismic::BasicTimeSeries::mlive, mspass::seismic::BasicTimeSeries::mt0, mspass::seismic::BasicTimeSeries::npts(), s, mspass::seismic::BasicTimeSeries::sample_number(), mspass::seismic::BasicTimeSeries::t0(), and mspass::seismic::BasicTimeSeries::tref.

◆ operator=()

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

Standard assignment operator.

67 {
68  if(this!=&tsi)
69  {
70  this->BasicTimeSeries::operator=(tsi);
71  this->Metadata::operator=(tsi);
72  s=tsi.s;
73  }
74  return(*this);
75 }
BasicTimeSeries & operator=(const BasicTimeSeries &parent)
Definition: BasicTimeSeries.cc:62
Metadata & operator=(const Metadata &mdold)
Definition: Metadata.cc:108

References mspass::seismic::BasicTimeSeries::operator=(), mspass::utility::Metadata::operator=(), and s.

◆ operator[]()

double mspass::seismic::CoreTimeSeries::operator[] ( size_t const  sample) const

Extract a sample from data vector with range checking. Because the data vector is public in this interface this operator is simply an alterative interface to this->s[sample].

Exceptions
SeisppErrorexception if the requested sample is outside the range of the data. Note this includes an implicit "outside" defined when the contents are marked dead.
Parameters
sampleis the integer sample number of data desired.
272 {
273  if(!mlive)
274  throw MsPASSError(string("CoreTimeSeries operator[]: attempting to access data marked as dead"),
275  ErrorSeverity::Invalid);
276  if( (i<0) || (i>=s.size()) )
277  {
278  throw MsPASSError(
279  string("CoreTimeSeries operator[]: request for sample outside range of data"),
280  ErrorSeverity::Invalid);
281  }
282  return(s[i]);
283 }

References mspass::seismic::BasicTimeSeries::mlive, and s.

◆ set_dt()

void mspass::seismic::CoreTimeSeries::set_dt ( const double  sample_interval)
virtual

Set the sample interval.

This method is complicated by the need to sync the changed value with Metadata. That is further complicated by the need to support aliases for the keys used to defined dt in Metadata. That is handled by first setting the internal dt value and then going through a fixed list of valid alias keys for dt. Any that exist are changed. If none were previously defined the unique name (see documentation) is added to Metadata.

Parameters
sample_intervalis the new data sample interval to be used.

Reimplemented from mspass::seismic::BasicTimeSeries.

186 {
187  this->BasicTimeSeries::set_dt(sample_interval);
188  /* This is the unique name - we always set it.
189  Feb 2021 - changed to used const string value set in keywords.h*/
190  this->put(SEISMICMD_dt,sample_interval);
191  /* these are hard coded aliases for sample_interval */
192  std::set<string> aliases;
193  std::set<string>::iterator aptr;
194  aliases.insert("dt");
195  for(aptr=aliases.begin();aptr!=aliases.end();++aptr)
196  {
197  if(this->is_defined(*aptr))
198  {
199  this->put(*aptr,sample_interval);
200  }
201  }
202 }
virtual void set_dt(const double sample_interval)
Set the sample interval.
Definition: BasicTimeSeries.h:198
bool is_defined(const std::string key) const noexcept
Definition: Metadata.cc:73
const std::string SEISMICMD_dt("delta")

References mspass::utility::Metadata::is_defined(), mspass::seismic::SEISMICMD_dt(), and mspass::seismic::BasicTimeSeries::set_dt().

◆ set_npts()

void mspass::seismic::CoreTimeSeries::set_npts ( const size_t  npts)
virtual

Set the number of samples attribute for data.

This method is complicated by the need to sync the changed value with Metadata. That is further complicated by the need to support aliases for the keys used to defined npts in Metadata. That is handled by first setting the internal npts value (actually ns) and then going through a fixed list of valid alias keys for npts. Any that exist are changed. If none were previously defined the unique name (see documentation) is added to Metadata.

This attribute has an additional complication compared to other setter that are overrides from BasicTimeSeries. That is, the number of points define the data buffer size to hold the sample data. To guarantee the buffer size and the internal remain consistent this method clears any existing content of the vector s and initializes npts points to 0.0. Note this means if one is using this to assemble a data object in pieces you MUST call this method before loading any data or it will be cleared and you will mysteriously find the data are all zeros.

Parameters
nptsis the new number of points to set.

Reimplemented from mspass::seismic::BasicTimeSeries.

223 {
225  /* This is the unique name - we always set it. Cast is necessary to
226  avoid type mismatch in python for unsigned.
227  Changed Feb 2021 to use key defined in in keywords.h*/
228  this->put(SEISMICMD_npts,(long int)npts);
229  /* these are hard coded aliases for sample_interval */
230  std::set<string> aliases;
231  std::set<string>::iterator aptr;
232  aliases.insert("nsamp");
233  aliases.insert("wfdisc.nsamp");
234  for(aptr=aliases.begin();aptr!=aliases.end();++aptr)
235  {
236  if(this->is_defined(*aptr))
237  {
238  this->put(*aptr,(long int)npts);
239  }
240  }
241  /* this method has the further complication that npts sets the size of the
242  data buffer. We clear it an initialize it to 0 to be consistent with
243  how constructors handle this. */
244  std::vector<double>().swap(this->s); // Clear the memory allocation of s
245  this->s.reserve(npts);
246  for(size_t i=0;i<npts;++i)this->s.push_back(0.0);
247 }
virtual void set_npts(const size_t npts)
Set the number of samples attribute for data.
Definition: BasicTimeSeries.h:213
const std::string SEISMICMD_npts("npts")

References mspass::utility::Metadata::is_defined(), mspass::seismic::BasicTimeSeries::npts(), s, mspass::seismic::SEISMICMD_npts(), and mspass::seismic::BasicTimeSeries::set_npts().

◆ set_t0()

void mspass::seismic::CoreTimeSeries::set_t0 ( const double  t0in)
virtual

Set the data start time.

This method is complicated by the need to sync the changed value with Metadata. That is further complicated by the need to support aliases for the keys used to defined npts in Metadata. That is handled by first setting the internal t0 value and then going through a fixed list of valid alias keys for it. Any that exist are changed. If none were previously defined the unique name (see documentation) is added to Metadata.

This is a dangerous method to use on real data as it can mess up the time if not handled correctly. It should be used only when that sharp knife is needed such as in assembling data outside of constructors in a test program.

Parameters
t0inis the new data sample interval to be used.

Reimplemented from mspass::seismic::BasicTimeSeries.

204 {
205  this->BasicTimeSeries::set_t0(t0in);
206  /* This is the unique name - we always set it.
207  Changed Feb 2021 to use const string value defined in keywords.h*/
208  this->put(SEISMICMD_t0,t0in);
209  /* these are hard coded aliases for sample_interval */
210  std::set<string> aliases;
211  std::set<string>::iterator aptr;
212  aliases.insert("t0");
213  aliases.insert("time");
214  for(aptr=aliases.begin();aptr!=aliases.end();++aptr)
215  {
216  if(this->is_defined(*aptr))
217  {
218  this->put(*aptr,t0in);
219  }
220  }
221 }
virtual void set_t0(const double t0in)
Set the data start time.
Definition: BasicTimeSeries.h:228
const std::string SEISMICMD_t0("starttime")

References mspass::utility::Metadata::is_defined(), mspass::seismic::SEISMICMD_t0(), and mspass::seismic::BasicTimeSeries::set_t0().

◆ sync_npts()

void mspass::seismic::CoreTimeSeries::sync_npts ( )

Sync the number of samples attribute with actual data size.

This method syncs the npts attribute with the actual size of the vector s. It also syncs aliases in the same way as the set_npts method.

249 {
250  if(nsamp != this->s.size()) {
251  this->BasicTimeSeries::set_npts(this->s.size());
252  /* This is the unique name - we always set it. The weird
253  cast is necessary to avoid type mismatch with unsigned
254  Changed Feb 2021 to use key defined in keywords.h*/
255  this->put(SEISMICMD_npts,(long int)nsamp);
256  /* these are hard coded aliases for sample_interval */
257  std::set<string> aliases;
258  std::set<string>::iterator aptr;
259  aliases.insert("nsamp");
260  aliases.insert("wfdisc.nsamp");
261  for(aptr=aliases.begin();aptr!=aliases.end();++aptr)
262  {
263  if(this->is_defined(*aptr))
264  {
265  this->put(*aptr,(long int)nsamp);
266  }
267  }
268  }
269 }

References mspass::utility::Metadata::is_defined(), mspass::seismic::BasicTimeSeries::nsamp, s, mspass::seismic::SEISMICMD_npts(), and mspass::seismic::BasicTimeSeries::set_npts().

Member Data Documentation

◆ s

std::vector<double> mspass::seismic::CoreTimeSeries::s

Actual data stored as an STL vector container. Note the STL guarantees the data elements of vector container are contiguous in memory like FORTRAN vectors. As a result things like the BLAS can be used with data object by using a syntax like this: if d is a CoreTimeSeries object, the address of the first sample of the data is &(d.s[0]).


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