#include "ThreeComponentSeismogram.h" using namespace SEISPP; class ThreeComponentSeismogram : public BasicTimeSeries , public Metadata { public: bool components_are_orthogonal; bool components_are_cardinal; double tmatrix[3][3]; dmatrix u; ThreeComponentSeismogram(); ThreeComponentSeismogram(int nsamples); ThreeComponentSeismogram(Metadata& md,bool load_data); ThreeComponentSeismogram(Dbptr db, MetadataList& mdl, AttributeMap& am); ThreeComponentSeismogram(const ThreeComponentSeismogram&); ThreeComponentSeismogram& operator = (const ThreeComponentSeismogram&); ~ThreeComponentSeismogram(); void zero_gaps(); void rotate_to_standard() throw(SeisppError); void rotate(SphericalCoordinate scor); void rotate(double nu[3]); void apply_transformation_matrix(double a[3][3]); void free_surface_transformation(SlownessVector u, double vp0, double vs0); };
A ThreeComponentSeismogram is a special type of time series data. It is special in two ways: (1) it is vector data and not scalar data and (2) that vector size matches our physical universe dimension of 3. Item (1) is encapsulated in this implementation by making this data object a derived class from a BasicTimeSeries object. Item (2) and the fact that the customers of this software are nearly guaranteed to be geophysicists lead to special member functions in this object that might otherwise be better left as external procedures.
A ThreeComponentSeismogram is a familiar concept to all seismologists. This implementation stores such data in a 3 by ns (defined in BasicTimeSeries (3)) matrix referenced with the symbol u to be consistent with the generic symbol used for vector seismic data in the bible of seismology by Aki and Richards.
A ThreeComponentSeismogram inherits three fundamental parameters from a BasicTimeSeries object definition: (1) dt = sample interval, (2) t0 = start time, and (3) ns = number of samples in this data object. It also inherits a couple of auxiliary parameters and data gap processing facilities described in BasicTimeSeries(3). It also inherits a Metadata object (see Metadata(3)) used basically as a variable length header in applications. See the man page for TimeSeries (3) for more details about motivation and use of Metadata for time series data like this.
The key data difference this object adds over a BasicTimeSeries are concepts related to coordinate transformations. This begins with two booleans added for convenience: components_are_orthogonal and components_are_cardinal. The verbose names should make their meaning clear with the statement that "cardinal" means the components are in standard geographical coordinates (+x1=east, +x2=north, and +x3=up). The user should recognize that the implicit assumption is that all transformation matrices are relative to the fixed reference frame of earth geographical coordinates. This has the well known failing for the existing GSN station at the South pole but handling that one special case was not something I chose to address.
The tmatrix stores the current transformation matrix of the data relative to the cardinal coordinate. Hence, if components_are_cardinal is true you should find that tmatrix is an identity.
Several member functions are provided to deal with transformation matrix issues. The rotate_to_standard standard function will return the data to cardinal coordinates using the current tmatrix. (Note if in doubt, call this function. If components_are_cardinal is set true the function will return immediately and do nothing to the data.)
The rotate function will apply a special earth specific transformation to a form of ray coordinate. It comes in two versions specified by a unit vector (double nu[3]) version) or angular definition (SphericalCoordinate version). Both are used to define a coordinate transformation in which: (a) the x3 direction is transformed to be parallel to the direction defined by the unit vector nu (SphericalCoordinate angles defined in scor); x2 (original north) is rotated to be radial from the direction defined by nu (defined by the direction of the projection of nu onto a horizontal plane); and x1 (original east) is transformed to a transverse component relative to nu (defined as a horizontal vector perpendicular to the transformed x2=radial and x3=longitudinal). Note that this yields data in an orthogonal reference frame.
The free_surface_transformation function is similar in some ways to the rotate function but it transforms the data to a coordinate system that is not orthogonal. I compute the free surface transformation matrix defined by Kennett (1991, GJI,vol. 104, p. 153-163). The result is a ray coordinate transformation with x1=transverse, x2=radial, and x3=longitudinal. This transformation is defined by the slowness (ray parameter) in s/km of a plane wave and the surface P (vp0) and S (vs0) velocities (both in units of km/s). Be aware also that not only is this transformation not orthogonal, it also changes the relative amplitudes of the P and S components. That is, it does not preserve original vector amplitudes. This function will throw a SeisppError object if the surface velocities passed to the function are higher than the apparent velocity (1/ray_parameter) as the transformation matrix used will not work for post-critical wavefields.
The apply_transformation_matrix functions more generically. It simply applies the 3x3 transformation matrix to the data.
Note that all coordinate transformations are accumulative. That is, multiple calls to any of the transformation routines will be accumulated in tmatrix. The notable exception is rotate_to_standard which will always return the data to the standard, cardinal coordinate system.
Symbol Metadata name ====Common with parallel TimeSeries constructor====== dt samprate t0 starttime ns nsamp ====Special for ThreeComponentSeismogram====== components_are_cardinal components_are_cardinal tmatrix[0]0] U11 tmatrix[1]0] U21 tmatrix[2]0] U31 tmatrix[0]1] U12 tmatrix[1]1] U22 tmatrix[2]1] U32 tmatrix[0]2] U13 tmatrix[1]2] U23 tmatrix[2]2] U33
Keyword Definition TimeReferenceType tref field of BasicTimeSeries datatype Data type ala CSS3.0 dtype attribute foff File offset in bytes to first sample dir Directory where data will be found dfile File name to read three_component_data_order Matrix order switch
There are two different dbsave functions that save a ThreeComponentSeismogram object's contents to a Datascope Database. The prototypes are:
void dbsave(ThreeComponentSeismogram& ts,Dbptr db,string table, MetadataList& md, AttributeMap& am); void dbsave(ThreeComponentSeismogram& ts,Dbptr db, string table, MetadataList& md, AttributeMap& am, vector<string>chanmap,bool output_as_standard);The best way to understand these functions is that the first is a simplified, overloaded version of the second (more general) function. The basic algorithm is that the three-component matrix is fragmented into three TimeSeries objects and the TimeSeries version of dbsave is called on each of the resulting three TimeSeries objects. The ts, db, table, md, and am arguments are exactly as described in the documentation for TimeSeries(3). The quick version is that the results are written to the Datascope database db in table "table" using the set of names in md with internal to external namespace mapping defined in am. This function has two frozen properties. First, the output channels are ALWAYS called: E, N, and Z. Second, to be consistent with this naming convention rotate_to_standard is always called before output.
The second version of this functions adds more generality. The primary distinction is the chanmap vector. The chanmap variable is assumed to be a three-component, STL vector of strings that define the names that should be assigned to the three components of the output TimeSeries. They are assumed to be in data order so that chanmap[0], for example, corresponds to row 0 of the 3xns data matrix used for the ThreeComponentSeismogram object. When the boolean variable output_as_standard is true the data are rotated to standard coordinate (rotate_to_standard is called) before sending the result to output.
Both of these functions throw SeisppError objects if there are problems on output.
Metadata(3), BasicTimeSeries(3), TimeSeries(3), http://geology.indiana.edu/pavlis/software/seispp/html/index.html
The mix of an object-oriented matrix implementation and the fixed C array tmatrix for the transformation matrix is a small potential confusion. I did this because a transformation matrix has a fixed dimension and so a fixed size C array works cleanly.
If you use the transformation matrix make sure you realize that this is the transformation used to transform the data from standard coordinates to the current state. When components_are_orthogonal is true the inverse transform is the transpose of tmatrix. Otherwise a matrix inversion is necessary, which is what rotate_to_standard does in that situation.