h5geo 0.4.0
C++17 and python API to work with geo-data (seismic, wells, maps, other in process) based on HDF5. Aimed at geoscientists and developers.
Loading...
Searching...
No Matches
H5Seis Class Referenceabstract

Provides API to work with seismic. More...

#include <h5seis.h>

Inheritance diagram for H5Seis:
H5BaseObject H5Base H5BaseImpl< H5Seis > H5BaseObjectImpl< H5Seis > H5SeisImpl

Public Member Functions

virtual bool readSEGYTextHeader (const std::string &segy, h5geo::TextEncoding encoding=static_cast< h5geo::TextEncoding >(0))=0
 Read text header from SEGY file.
 
virtual bool readSEGYBinHeader (const std::string &segy, h5geo::Endian endian=static_cast< h5geo::Endian >(0))=0
 Read binary header from SEGY file.
 
virtual bool readSEGYTraces (const std::vector< std::string > &segyFiles, std::vector< h5geo::SegyFormat > formats=std::vector< h5geo::SegyFormat >(), std::vector< h5geo::Endian > endians=std::vector< h5geo::Endian >(), std::vector< std::vector< std::string > > trcHdrNamesArr=std::vector< std::vector< std::string > >(), size_t trcBuffer=10000, std::function< void(double)> progressCallback=nullptr)=0
 Read trace headers and trace data from SEGY file.
 
virtual bool readSEGYTracesMMap (const std::vector< std::string > &segyFiles, std::vector< h5geo::SegyFormat > formats=std::vector< h5geo::SegyFormat >(), std::vector< h5geo::Endian > endians=std::vector< h5geo::Endian >(), std::vector< std::vector< std::string > > trcHdrNamesArr=std::vector< std::vector< std::string > >(), size_t trcBuffer=10000, int nThreads=-1, std::function< void(double)> progressCallback=nullptr)=0
 Read trace headers and trace data from SEGY file using Memory Mapping.
 
virtual bool writeTextHeader (const char(&txtHdr)[40][80])=0
 Write text header.
 
virtual bool writeTextHeader (const std::vector< std::string > &txtHdr)=0
 Write text header.
 
virtual bool writeBinHeader (const double(&binHdr)[30])=0
 Write binary header.
 
virtual bool writeBinHeader (const std::vector< double > &binHdrVec)=0
 Write text header.
 
virtual bool writeBinHeader (const Eigen::Ref< const Eigen::VectorXd > &binHdrVec)=0
 Write text header.
 
virtual bool writeBinHeader (const std::string &hdrName, const double &value, const std::string &unitsFrom="", const std::string &unitsTo="")=0
 Write text header.
 
virtual bool writeTrace (Eigen::Ref< Eigen::MatrixXf > TRACE, const size_t &fromTrc=0, const size_t &fromSampInd=0, const std::string &dataUnits="")=0
 Write block of traces starting from trace fromTrc and from sample fromSampInd
 
virtual bool writeTrace (Eigen::Ref< Eigen::MatrixXf > TRACE, const Eigen::Ref< const Eigen::VectorX< size_t > > &trcInd, const size_t &fromSampInd=0, const std::string &dataUnits="")=0
 Write traces using indexes.
 
virtual bool writeTraceHeader (const Eigen::Ref< const Eigen::MatrixXd > &HDR, const size_t &fromTrc=0, const size_t &fromHdrInd=0)=0
 Write block of trace headers starting from trace fromTrc and from header index fromHdrInd
 
virtual bool writeTraceHeader (const std::string &hdrName, Eigen::Ref< Eigen::MatrixXd > hdr, const size_t &fromTrc=0, const std::string &unitsFrom="", const std::string &unitsTo="")=0
 Write trace header by name.
 
virtual bool writeTraceHeader (const std::string &hdrName, Eigen::Ref< Eigen::MatrixXd > hdr, const Eigen::Ref< const Eigen::VectorX< size_t > > &trcInd, const std::string &unitsFrom="", const std::string &unitsTo="")=0
 Write trace header by name and trace indexes.
 
virtual bool writeXYTraceHeaders (const std::vector< std::string > &xyHdrNames, Eigen::Ref< Eigen::MatrixX2d > &xy, const size_t &fromTrc=0, const std::string &lengthUnits="", bool doCoordTransform=false)=0
 Write XY trace headers (two columns in Eigen column-major matrix)
 
virtual bool writeXYTraceHeaders (const std::vector< std::string > &xyHdrNames, Eigen::Ref< Eigen::MatrixX2d > &xy, const Eigen::Ref< const Eigen::VectorX< size_t > > &trcInd, const std::string &lengthUnits="", bool doCoordTransform=false)=0
 Write XY trace headers (two columns in Eigen column-major matrix)
 
virtual bool setNTrc (size_t nTrc)=0
 Resize trace and trace header DataSets.
 
virtual bool setNSamp (size_t nSamp)=0
 Resize trace DataSet.
 
virtual std::vector< std::string > getTextHeader ()=0
 Get text header.
 
virtual std::map< std::string, double > getBinHeader ()=0
 Get binary header.
 
virtual double getBinHeader (const std::string &hdrName, const std::string &unitsFrom="", const std::string &unitsTo="")=0
 Get text header by name.
 
virtual Eigen::MatrixXf getTrace (const size_t &fromTrc, size_t nTrc=1, const size_t &fromSampInd=0, size_t nSamp=std::numeric_limits< size_t >::max(), const std::string &dataUnits="")=0
 Get block of traces.
 
virtual Eigen::MatrixXf getTrace (const Eigen::Ref< const Eigen::VectorX< size_t > > &trcInd, const size_t &fromSampInd=0, size_t nSamp=std::numeric_limits< size_t >::max(), const std::string &dataUnits="")=0
 Get traces by indexes.
 
virtual Eigen::MatrixXd getTraceHeader (const size_t &fromTrc, size_t nTrc=1, const size_t &fromHdr=0, size_t nHdr=std::numeric_limits< size_t >::max(), const std::vector< std::string > &unitsFrom=std::vector< std::string >(), const std::vector< std::string > &unitsTo=std::vector< std::string >())=0
 Get block of trace headers.
 
virtual Eigen::VectorXd getTraceHeader (const std::string &hdrName, const size_t &fromTrc=0, size_t nTrc=1, const std::string &unitsFrom="", const std::string &unitsTo="")=0
 Get block of trace header by name.
 
virtual Eigen::MatrixXd getTraceHeader (const std::vector< size_t > &trcInd, const std::vector< size_t > &trcHdrInd, const std::vector< std::string > &unitsFrom=std::vector< std::string >(), const std::vector< std::string > &unitsTo=std::vector< std::string >())=0
 Get trace headers by indexes.
 
virtual Eigen::MatrixXd getTraceHeader (const Eigen::Ref< const Eigen::VectorX< size_t > > &trcInd, const Eigen::Ref< const Eigen::VectorX< size_t > > &trcHdrInd, const std::vector< std::string > &unitsFrom=std::vector< std::string >(), const std::vector< std::string > &unitsTo=std::vector< std::string >())=0
 Get trace headers by indexes.
 
virtual Eigen::MatrixXd getTraceHeader (const std::vector< std::string > &hdrNames, const std::vector< size_t > &trcInd, const std::vector< std::string > &unitsFrom=std::vector< std::string >(), const std::vector< std::string > &unitsTo=std::vector< std::string >())=0
 Get trace header by indexes and name.
 
virtual Eigen::MatrixXd getTraceHeader (const std::vector< std::string > &hdrNames, const Eigen::Ref< const Eigen::VectorX< size_t > > &trcInd, const std::vector< std::string > &unitsFrom=std::vector< std::string >(), const std::vector< std::string > &unitsTo=std::vector< std::string >())=0
 Get trace headers by indexes and names.
 
virtual Eigen::MatrixXd getXYTraceHeaders (const std::vector< std::string > &xyHdrNames, const size_t &fromTrc=0, size_t nTrc=std::numeric_limits< size_t >::max(), const std::string &lengthUnits="", bool doCoordTransform=false)=0
 Get XY trace headers (two columns in Eigen column-major matrix)
 
virtual Eigen::MatrixXd getXYTraceHeaders (const std::vector< std::string > &xyHdrNames, const Eigen::Ref< const Eigen::VectorX< size_t > > &trcInd, const std::string &lengthUnits="", bool doCoordTransform=false)=0
 Get XY trace headers (two columns in Eigen column-major matrix)
 
virtual Eigen::VectorX< size_t > getSortedData (Eigen::MatrixXf &TRACE, Eigen::MatrixXd &HDR, const std::vector< std::string > &keyList, const std::vector< double > &minList, const std::vector< double > &maxList, size_t pStep=1, size_t fromSampInd=0, size_t nSamp=std::numeric_limits< size_t >::max(), const std::string &dataUnits="", const std::string &lengthUnits="", bool doCoordTransform=false)=0
 Get sorted data based on precalculated primary sort keys.
 
virtual ptrdiff_t getBinHeaderIndex (const std::string &hdrName)=0
 Get index (position within 1D dataset) for a given binary header.
 
virtual ptrdiff_t getTraceHeaderIndex (const std::string &hdrName)=0
 Get index (row/col within 2D dataset) for a given trace header.
 
virtual Eigen::VectorXd getSamples (const size_t &trcInd, const std::string &units="")=0
 Get vector of equally spaced samples in specified units.
 
virtual double getFirstSample (const size_t &trcInd, const std::string &units="")=0
 Get first sample.
 
virtual double getLastSample (const size_t &trcInd, const std::string &units="")=0
 Get last sample.
 
virtual double getSampRate (const std::string &units="")=0
 Get sampling rate.
 
virtual size_t getNSamp ()=0
 Get number of samples.
 
virtual size_t getNTrc ()=0
 Get number of traces.
 
virtual size_t getNTrcHdr ()=0
 Get number trace headers (usually 78)
 
virtual size_t getNBinHdr ()=0
 Get number of binary headers.
 
virtual size_t getNTextHdrRows ()=0
 Get number of text header rows.
 
virtual Eigen::VectorX< size_t > getPKeyIndexes (const std::string &pKey, double pMin, double pMax, size_t pStep=1)=0
 Get trace indexes for given PKey
 
virtual Eigen::VectorXd getPKeyValues (const std::string &pKey, const std::string &unitsFrom="", const std::string &unitsTo="")=0
 Get PKey unique values.
 
virtual size_t getPKeySize (const std::string &pKey)=0
 Get number of unique values for a given PKey
 
virtual size_t getPKeySize (const std::string &pKey, double pMin, double pMax, size_t pStep=1)=0
 Get number of unique values for a given PKey
 
virtual size_t getPKeyTraceSize (const std::string &pKey, double pMin, double pMax, size_t pStep=1)=0
 Get number of traces to be selected for a given PKey
 
virtual std::vector< std::string > getPKeyNames ()=0
 Get names of prepared PKeys (names of prepared sortings PKeys)
 
virtual std::map< std::string, double > getTraceHeaderMin ()=0
 Get trace header minimal values.
 
virtual std::map< std::string, double > getTraceHeaderMax ()=0
 Get trace header maximal values.
 
virtual double getTraceHeaderMin (const std::string &hdrName, const std::string &unitsFrom="", const std::string &unitsTo="")=0
 Get trace header minimal value for a given trace header.
 
virtual double getTraceHeaderMax (const std::string &hdrName, const std::string &unitsFrom="", const std::string &unitsTo="")=0
 Get trace header maximal value for a given trace header.
 
virtual H5SeisParam getParam ()=0
 Get parameters that were used to create current seis.
 
virtual bool checkTraceLimits (const size_t &fromTrc, size_t &nTrc)=0
 Check fromTrc, nTrc (passed by reference) and diminish nTrc to fit in data limits (if fromTrc is inside limit)
 
virtual bool checkTraceHeaderLimits (const size_t &fromHdr, size_t &nHdr)=0
 Check fromHdr and nHdr (passed by reference) and diminish nHdr to fit in data limits (if fromTrc is inside limit)
 
virtual bool checkSampleLimits (const size_t &fromSampInd, size_t &nSamp)=0
 Check fromSampInd and nSamp (passed by reference) and diminish nSamp to fit in data limits (if fromSampInd is inside limit)
 
virtual bool generatePRESTKGeometry (double src_x0, double src_dx, size_t src_nx, double src_y0, double src_dy, size_t src_ny, double src_z, double rec_x0, double rec_dx, size_t rec_nx, double rec_y0, double rec_dy, size_t rec_ny, double rec_z, double orientation, bool moveRec, const std::string &lengthUnits="", bool doCoordTransform=false)=0
 Convenient function to prepare geometry for PRESTACK data.
 
virtual bool generateSTKGeometry (double x0, double dx, size_t nx, double y0, double dy, size_t ny, double z, double orientation, const std::string &lengthUnits="", bool doCoordTransform=false)=0
 Convenient function to prepare geometry for STACK data.
 
virtual bool setDomain (const h5geo::Domain &domain)=0
 Set domain for the seismic (TVD, TVDSS, TWT, OWT)
 
virtual bool setDataType (const h5geo::SeisDataType &seisType)=0
 Set datatype for the seismic (STACK or PRESTACK)
 
virtual bool setSurveyType (const h5geo::SurveyType &surveyType)=0
 Set survey type for the seismic (TWO_D or THREE_D)
 
virtual bool setSRD (double val, const std::string &lengthUnits="")=0
 Set Seismic Reference Datum.
 
virtual bool setSampRate (double val, const std::string &units="")=0
 Set sampling rate.
 
virtual bool setFirstSample (double val, const std::string &units="")=0
 Set first sample.
 
virtual h5geo::Domain getDomain ()=0
 Get domain (TVD, TVDSS, TWT, OWT)
 
virtual h5geo::SeisDataType getDataType ()=0
 Get datatype for the seismic (STACK or PRESTACK)
 
virtual h5geo::SurveyType getSurveyType ()=0
 Set survey type for the seismic (TWO_D or THREE_D)
 
virtual double getSRD (const std::string &lengthUnits="")=0
 Get Seismic Reference Datum.
 
virtual bool hasPKeySort (const std::string &pKeyName)=0
 Check if PKey sort is prepared.
 
virtual bool removePKeySort (const std::string &pKeyName)=0
 Remove PKey sorting.
 
virtual bool addPKeySort (const std::string &pKeyName)=0
 Prepare sorting.
 
virtual bool updateTraceHeaderSampRate ()=0
 Set trace header samp rate from binary header.
 
virtual bool updateTraceHeaderNSamp ()=0
 Set trace header number of samples from binary header.
 
virtual H5SeisContaineropenSeisContainer ()=0
 Open H5SeisContainer where current seismic resides.
 
virtual std::optional< h5gt::DataSet > getTextHeaderD ()=0
 Get text header DataSet.
 
virtual std::optional< h5gt::DataSet > getBinHeaderD ()=0
 Get binary header DataSet.
 
virtual std::optional< h5gt::DataSet > getTraceHeaderD ()=0
 Get trace header DataSet.
 
virtual std::optional< h5gt::DataSet > getTraceD ()=0
 Get trace DataSet.
 
virtual std::optional< h5gt::Group > getSortG ()=0
 Get sorting Group.
 
virtual std::optional< h5gt::Group > getUValG ()=0
 Get sorting unique values Group.
 
virtual std::optional< h5gt::Group > getIndexesG ()=0
 Get sorting indexes Group.
 
virtual std::optional< h5gt::Group > getSEGYG ()=0
 Get SEGY Group (for mapped H5Seis only)
 
virtual std::optional< h5gt::DataSet > getSEGYTextHeaderD ()=0
 Get SEGY text header DataSet (for mapped H5Seis only)
 
virtual std::optional< h5gt::DataSet > getSEGYBinHeader2BytesD ()=0
 Get SEGY 2-bytes binary header DataSet (for mapped H5Seis only)
 
virtual std::optional< h5gt::DataSet > getSEGYBinHeader4BytesD ()=0
 Get SEGY 4-bytes binary header DataSet (for mapped H5Seis only)
 
virtual std::optional< h5gt::DataSet > getSEGYTraceHeader2BytesD ()=0
 Get SEGY 2-bytes trace header DataSet (for mapped H5Seis only)
 
virtual std::optional< h5gt::DataSet > getSEGYTraceHeader4BytesD ()=0
 Get SEGY 4-bytes trace header DataSet (for mapped H5Seis only)
 
virtual std::optional< h5gt::DataSet > getSEGYTraceFloatD ()=0
 Get SEGY float trace DataSet (for mapped H5Seis only)
 
virtual bool updateTraceHeaderLimits (size_t nTrcBuffer=1e7)=0
 Calculate and write min/max trace headers.
 
virtual bool updatePKeySort (const std::string &pKeyName)=0
 Update sorting for prepared PKey
 
virtual Eigen::MatrixXd calcBoundary (const std::string &lengthUnits="", bool doCoordTransform=false)=0
 Calculate XY boundary around the survey.
 
virtual bool exportToVol (H5Vol *vol, const std::string &xHeader="CDP_X", const std::string &yHeader="CDP_Y", const std::string &ilHeader="INLINE", const std::string &xlHeader="XLINE", double ilMin=std::numeric_limits< double >::min(), double ilMax=std::numeric_limits< double >::max(), double xlMin=std::numeric_limits< double >::min(), double xlMax=std::numeric_limits< double >::max(), size_t fromSampInd=0, size_t nSamp=std::numeric_limits< size_t >::max(), std::function< void(double)> progressCallback=nullptr)=0
 Export seismic to H5Vol.
 
virtual bool exportToSEGY (const std::string &segyFile, size_t trcBuffer=10000, h5geo::Endian endian=h5geo::Endian::Big, std::function< void(double)> progressCallback=nullptr)=0
 
- Public Member Functions inherited from H5BaseObject
virtual H5BasePointsopenPoints (const std::string &name)=0
 Open H5BasePoints derived points.
 
virtual H5BasePointsopenPoints (h5gt::Group group)=0
 Open H5BasePoints derived points.
 
virtual H5HorizonopenHorizon (const std::string &name)=0
 
virtual H5HorizonopenHorizon (h5gt::Group group)=0
 
virtual H5Points1createPoints1 (std::string &name, H5PointsParam &p, h5geo::CreationType createFlag)=0
 
virtual H5Points1createPoints1 (h5gt::Group group, H5PointsParam &p, h5geo::CreationType createFlag)=0
 
virtual H5Points2createPoints2 (std::string &name, H5PointsParam &p, h5geo::CreationType createFlag)=0
 
virtual H5Points2createPoints2 (h5gt::Group group, H5PointsParam &p, h5geo::CreationType createFlag)=0
 
virtual H5Points3createPoints3 (std::string &name, H5PointsParam &p, h5geo::CreationType createFlag)=0
 
virtual H5Points3createPoints3 (h5gt::Group group, H5PointsParam &p, h5geo::CreationType createFlag)=0
 
virtual H5Points4createPoints4 (std::string &name, H5PointsParam &p, h5geo::CreationType createFlag)=0
 
virtual H5Points4createPoints4 (h5gt::Group group, H5PointsParam &p, h5geo::CreationType createFlag)=0
 
virtual H5HorizoncreateHorizon (std::string &name, H5HorizonParam &p, h5geo::CreationType createFlag)=0
 
virtual H5HorizoncreateHorizon (h5gt::Group group, H5HorizonParam &p, h5geo::CreationType createFlag)=0
 
virtual bool setSpatialReference (const std::string &str)=0
 Set spatial reference for current geo-object using authName:code form.
 
virtual bool setSpatialReference (const std::string &authName, const std::string &code)=0
 Set spatial reference for current geo-object.
 
virtual bool setLengthUnits (const std::string &str)=0
 Set length units for the current geo-object.
 
virtual bool setTemporalUnits (const std::string &str)=0
 Set temporal units for the current geo-object.
 
virtual bool setAngularUnits (const std::string &str)=0
 Set angular units for the current geo-object.
 
virtual bool setDataUnits (const std::string &str)=0
 Set data units for the current geo-object.
 
virtual bool setNullValue (double val)=0
 Set NULL value for the current geo-object.
 
virtual std::string getSpatialReference ()=0
 Get spatial reference for current geo-object.
 
virtual std::string getLengthUnits ()=0
 Get length units for the current geo-object.
 
virtual std::string getTemporalUnits ()=0
 Get temporal units for the current geo-object.
 
virtual std::string getAngularUnits ()=0
 Get angular units for the current geo-object.
 
virtual std::string getDataUnits ()=0
 Get data units for the current geo-object.
 
virtual double getNullValue ()=0
 Get NULL value for the current geo-object.
 
virtual h5gt::File getH5File () const =0
 Get HDF5 file.
 
virtual h5gt::Group getObjG () const =0
 Get HDF5 Group.
 
virtual std::string getName () const =0
 Get geo-object's name without path.
 
virtual std::string getFullName () const =0
 Get geo-object's name with full path to that object.
 
virtual std::vector< h5gt::Group > getObjGroupList (const h5geo::ObjectType &objType, bool recursive)=0
 Find all geo-objects of specified type within current geo-object and return them as vector of Groups.
 
virtual std::vector< std::string > getObjNameList (const h5geo::ObjectType &objType, bool recursive)=0
 Find all geo-objects of specified type within current geo-object and return them as vector of names.
 
virtual size_t getObjCount (const h5geo::ObjectType &objType, bool recursive)=0
 Get number of geo-objects of specified type within current geo-object.
 
virtual bool isEqual (H5BaseObject *other) const =0
 Check if geo-objects are the same.
 
virtual bool operator== (const H5BaseObject &other) const =0
 Check if geo-objects are the same (compares HDF5 Groups)
 
virtual bool operator!= (const H5BaseObject &other) const =0
 Check if geo-objects are not the same (compares HDF5 Groups)
 
- Public Member Functions inherited from H5Base
virtual H5Baseclone ()=0
 
virtual void Delete ()=0
 

Detailed Description

Provides API to work with seismic.

Seis is complex object that is represented by:

  • text header dataset
  • binary header dataset
  • trace header dataset
  • trace data dataset
  • optional datasets for sorting and other purposes

Dont forget to call:

It is designed to be fast and convenient.
SEGY reader is included.

Member Function Documentation

◆ addPKeySort()

virtual bool H5Seis::addPKeySort ( const std::string & pKeyName)
pure virtual

Prepare sorting.

Sorting is needed for any sorting related operations (like H5Seis::getSortedData()).
If you plan to get CDP-DSREG data, you have to call H5Seis::addPKeySort("CDP") first.

Implemented in H5SeisImpl.

◆ calcBoundary()

virtual Eigen::MatrixXd H5Seis::calcBoundary ( const std::string & lengthUnits = "",
bool doCoordTransform = false )
pure virtual

Calculate XY boundary around the survey.

Return two cols Eigen matrix. Use it to write as H5Horizon.

Implemented in H5SeisImpl.

◆ checkSampleLimits()

virtual bool H5Seis::checkSampleLimits ( const size_t & fromSampInd,
size_t & nSamp )
pure virtual

Check fromSampInd and nSamp (passed by reference) and diminish nSamp to fit in data limits (if fromSampInd is inside limit)

Parameters
fromSampIndfirst sample index
nSampnumber of samples (to read for example)
Returns

Implemented in H5SeisImpl.

◆ checkTraceHeaderLimits()

virtual bool H5Seis::checkTraceHeaderLimits ( const size_t & fromHdr,
size_t & nHdr )
pure virtual

Check fromHdr and nHdr (passed by reference) and diminish nHdr to fit in data limits (if fromTrc is inside limit)

Parameters
fromHdrfirst header (usually there are 78 headers so the value should be less then this value)
nHdrnumber of headers (to read for example)
Returns

Implemented in H5SeisImpl.

◆ checkTraceLimits()

virtual bool H5Seis::checkTraceLimits ( const size_t & fromTrc,
size_t & nTrc )
pure virtual

Check fromTrc, nTrc (passed by reference) and diminish nTrc to fit in data limits (if fromTrc is inside limit)

Parameters
fromTrcfirst trace (to read for example)
nTrcnumber of trace (to read for example)
Returns

Implemented in H5SeisImpl.

◆ exportToVol()

virtual bool H5Seis::exportToVol ( H5Vol * vol,
const std::string & xHeader = "CDP_X",
const std::string & yHeader = "CDP_Y",
const std::string & ilHeader = "INLINE",
const std::string & xlHeader = "XLINE",
double ilMin = std::numeric_limits< double >::min(),
double ilMax = std::numeric_limits< double >::max(),
double xlMin = std::numeric_limits< double >::min(),
double xlMax = std::numeric_limits< double >::max(),
size_t fromSampInd = 0,
size_t nSamp = std::numeric_limits< size_t >::max(),
std::function< void(double)> progressCallback = nullptr )
pure virtual

Export seismic to H5Vol.

Note
Selected traces must shape a rectangle.

Implemented in H5SeisImpl.

◆ generatePRESTKGeometry()

virtual bool H5Seis::generatePRESTKGeometry ( double src_x0,
double src_dx,
size_t src_nx,
double src_y0,
double src_dy,
size_t src_ny,
double src_z,
double rec_x0,
double rec_dx,
size_t rec_nx,
double rec_y0,
double rec_dy,
size_t rec_ny,
double rec_z,
double orientation,
bool moveRec,
const std::string & lengthUnits = "",
bool doCoordTransform = false )
pure virtual

Convenient function to prepare geometry for PRESTACK data.

Parameters
orientationcounter clock angle (in radians)

Implemented in H5SeisImpl.

◆ generateSTKGeometry()

virtual bool H5Seis::generateSTKGeometry ( double x0,
double dx,
size_t nx,
double y0,
double dy,
size_t ny,
double z,
double orientation,
const std::string & lengthUnits = "",
bool doCoordTransform = false )
pure virtual

Convenient function to prepare geometry for STACK data.

Parameters
orientationcounter clock angle (in radians)

Implemented in H5SeisImpl.

◆ getBinHeader() [1/2]

virtual std::map< std::string, double > H5Seis::getBinHeader ( )
pure virtual

Get binary header.

Implemented in H5SeisImpl.

◆ getBinHeader() [2/2]

virtual double H5Seis::getBinHeader ( const std::string & hdrName,
const std::string & unitsFrom = "",
const std::string & unitsTo = "" )
pure virtual

Get text header by name.

Implemented in H5SeisImpl.

◆ getBinHeaderD()

virtual std::optional< h5gt::DataSet > H5Seis::getBinHeaderD ( )
pure virtual

Get binary header DataSet.

Implemented in H5SeisImpl.

◆ getBinHeaderIndex()

virtual ptrdiff_t H5Seis::getBinHeaderIndex ( const std::string & hdrName)
pure virtual

Get index (position within 1D dataset) for a given binary header.

Implemented in H5SeisImpl.

◆ getDataType()

virtual h5geo::SeisDataType H5Seis::getDataType ( )
pure virtual

Get datatype for the seismic (STACK or PRESTACK)

Implemented in H5SeisImpl.

◆ getDomain()

virtual h5geo::Domain H5Seis::getDomain ( )
pure virtual

Get domain (TVD, TVDSS, TWT, OWT)

Implemented in H5SeisImpl.

◆ getFirstSample()

virtual double H5Seis::getFirstSample ( const size_t & trcInd,
const std::string & units = "" )
pure virtual

Get first sample.

Implemented in H5SeisImpl.

◆ getIndexesG()

virtual std::optional< h5gt::Group > H5Seis::getIndexesG ( )
pure virtual

Get sorting indexes Group.

Implemented in H5SeisImpl.

◆ getLastSample()

virtual double H5Seis::getLastSample ( const size_t & trcInd,
const std::string & units = "" )
pure virtual

Get last sample.

Implemented in H5SeisImpl.

◆ getNBinHdr()

virtual size_t H5Seis::getNBinHdr ( )
pure virtual

Get number of binary headers.

Implemented in H5SeisImpl.

◆ getNSamp()

virtual size_t H5Seis::getNSamp ( )
pure virtual

Get number of samples.

Implemented in H5SeisImpl.

◆ getNTextHdrRows()

virtual size_t H5Seis::getNTextHdrRows ( )
pure virtual

Get number of text header rows.

Implemented in H5SeisImpl.

◆ getNTrc()

virtual size_t H5Seis::getNTrc ( )
pure virtual

Get number of traces.

Implemented in H5SeisImpl.

◆ getNTrcHdr()

virtual size_t H5Seis::getNTrcHdr ( )
pure virtual

Get number trace headers (usually 78)

Implemented in H5SeisImpl.

◆ getParam()

virtual H5SeisParam H5Seis::getParam ( )
pure virtual

Get parameters that were used to create current seis.

Implemented in H5SeisImpl.

◆ getPKeyIndexes()

virtual Eigen::VectorX< size_t > H5Seis::getPKeyIndexes ( const std::string & pKey,
double pMin,
double pMax,
size_t pStep = 1 )
pure virtual

Get trace indexes for given PKey

Before using it one should prepare primary sort keys with H5Seis::addPKeySort() method.

Parameters
pKeyprimary key name
pMinprimary key minimal value
pMaxprimary key maximal value
pStepprimary key step (0 and 1 means every pKey; else every n-th selected pKey)

Implemented in H5SeisImpl.

◆ getPKeyNames()

virtual std::vector< std::string > H5Seis::getPKeyNames ( )
pure virtual

Get names of prepared PKeys (names of prepared sortings PKeys)

Implemented in H5SeisImpl.

◆ getPKeySize() [1/2]

virtual size_t H5Seis::getPKeySize ( const std::string & pKey)
pure virtual

Get number of unique values for a given PKey

Before using it one should prepare primary sort keys with H5Seis::addPKeySort() method.

Implemented in H5SeisImpl.

◆ getPKeySize() [2/2]

virtual size_t H5Seis::getPKeySize ( const std::string & pKey,
double pMin,
double pMax,
size_t pStep = 1 )
pure virtual

Get number of unique values for a given PKey

Before using it one should prepare primary sort keys with H5Seis::addPKeySort() method.

Implemented in H5SeisImpl.

◆ getPKeyTraceSize()

virtual size_t H5Seis::getPKeyTraceSize ( const std::string & pKey,
double pMin,
double pMax,
size_t pStep = 1 )
pure virtual

Get number of traces to be selected for a given PKey

Before using it one should prepare primary sort keys with H5Seis::addPKeySort() method.

Parameters
pKeyprimary key name
pMinprimary key minimal value
pMaxprimary key maximal value
pStepprimary key step (0 and 1 means every pKey; else every n-th selected pKey)

Implemented in H5SeisImpl.

◆ getPKeyValues()

virtual Eigen::VectorXd H5Seis::getPKeyValues ( const std::string & pKey,
const std::string & unitsFrom = "",
const std::string & unitsTo = "" )
pure virtual

Get PKey unique values.

Before using it one should prepare primary sort keys with H5Seis::addPKeySort() method.

Implemented in H5SeisImpl.

◆ getSamples()

virtual Eigen::VectorXd H5Seis::getSamples ( const size_t & trcInd,
const std::string & units = "" )
pure virtual

Get vector of equally spaced samples in specified units.

Implemented in H5SeisImpl.

◆ getSampRate()

virtual double H5Seis::getSampRate ( const std::string & units = "")
pure virtual

Get sampling rate.

Implemented in H5SeisImpl.

◆ getSEGYBinHeader2BytesD()

virtual std::optional< h5gt::DataSet > H5Seis::getSEGYBinHeader2BytesD ( )
pure virtual

Get SEGY 2-bytes binary header DataSet (for mapped H5Seis only)

Implemented in H5SeisImpl.

◆ getSEGYBinHeader4BytesD()

virtual std::optional< h5gt::DataSet > H5Seis::getSEGYBinHeader4BytesD ( )
pure virtual

Get SEGY 4-bytes binary header DataSet (for mapped H5Seis only)

Implemented in H5SeisImpl.

◆ getSEGYG()

virtual std::optional< h5gt::Group > H5Seis::getSEGYG ( )
pure virtual

Get SEGY Group (for mapped H5Seis only)

Implemented in H5SeisImpl.

◆ getSEGYTextHeaderD()

virtual std::optional< h5gt::DataSet > H5Seis::getSEGYTextHeaderD ( )
pure virtual

Get SEGY text header DataSet (for mapped H5Seis only)

Implemented in H5SeisImpl.

◆ getSEGYTraceFloatD()

virtual std::optional< h5gt::DataSet > H5Seis::getSEGYTraceFloatD ( )
pure virtual

Get SEGY float trace DataSet (for mapped H5Seis only)

Implemented in H5SeisImpl.

◆ getSEGYTraceHeader2BytesD()

virtual std::optional< h5gt::DataSet > H5Seis::getSEGYTraceHeader2BytesD ( )
pure virtual

Get SEGY 2-bytes trace header DataSet (for mapped H5Seis only)

Implemented in H5SeisImpl.

◆ getSEGYTraceHeader4BytesD()

virtual std::optional< h5gt::DataSet > H5Seis::getSEGYTraceHeader4BytesD ( )
pure virtual

Get SEGY 4-bytes trace header DataSet (for mapped H5Seis only)

Implemented in H5SeisImpl.

◆ getSortedData()

virtual Eigen::VectorX< size_t > H5Seis::getSortedData ( Eigen::MatrixXf & TRACE,
Eigen::MatrixXd & HDR,
const std::vector< std::string > & keyList,
const std::vector< double > & minList,
const std::vector< double > & maxList,
size_t pStep = 1,
size_t fromSampInd = 0,
size_t nSamp = std::numeric_limits< size_t >::max(),
const std::string & dataUnits = "",
const std::string & lengthUnits = "",
bool doCoordTransform = false )
pure virtual

Get sorted data based on precalculated primary sort keys.

Before using it one should prepare primary sort keys with H5Seis::addPKeySort() method.

Parameters
TRACEnot Eigen::Ref<> because Eigen::Ref<> doesn't allow to resize matrices
HDRnot Eigen::Ref<> because Eigen::Ref<> doesn't allow to resize matrices
keyListtrace header names to be worked with (first is treated as PKey)
minListminimal value for each key
maxListmaximal value for each key
pStepretrieve only n-th pkeys (for ex: SP=[1,5,6,8,9], pStep=2 -> SP_out=[1,6,9])
fromSampIndfirst sample index to read (in range [0, H5Seis::getNSamp()])
nSampnumber of samples to be read (if 0 then TRACE will be empty). By default all samples
dataUnitstansform data to these units
lengthUnitsworks only in pair with 'doCoordTransform'
doCoordTransformworks only if two header names are passed (one for X and another Y)
Returns
vector of trace indexes read

Implemented in H5SeisImpl.

◆ getSortG()

virtual std::optional< h5gt::Group > H5Seis::getSortG ( )
pure virtual

Get sorting Group.

Implemented in H5SeisImpl.

◆ getSRD()

virtual double H5Seis::getSRD ( const std::string & lengthUnits = "")
pure virtual

Get Seismic Reference Datum.

Implemented in H5SeisImpl.

◆ getSurveyType()

virtual h5geo::SurveyType H5Seis::getSurveyType ( )
pure virtual

Set survey type for the seismic (TWO_D or THREE_D)

Implemented in H5SeisImpl.

◆ getTextHeader()

virtual std::vector< std::string > H5Seis::getTextHeader ( )
pure virtual

Get text header.

Implemented in H5SeisImpl.

◆ getTextHeaderD()

virtual std::optional< h5gt::DataSet > H5Seis::getTextHeaderD ( )
pure virtual

Get text header DataSet.

Implemented in H5SeisImpl.

◆ getTrace() [1/2]

virtual Eigen::MatrixXf H5Seis::getTrace ( const Eigen::Ref< const Eigen::VectorX< size_t > > & trcInd,
const size_t & fromSampInd = 0,
size_t nSamp = std::numeric_limits< size_t >::max(),
const std::string & dataUnits = "" )
pure virtual

Get traces by indexes.

Return empty matrix if max trcInd exceeds nTrc.
That is done intensionally to keep trcInd size and returned number of traces equal.

Implemented in H5SeisImpl.

◆ getTrace() [2/2]

virtual Eigen::MatrixXf H5Seis::getTrace ( const size_t & fromTrc,
size_t nTrc = 1,
const size_t & fromSampInd = 0,
size_t nSamp = std::numeric_limits< size_t >::max(),
const std::string & dataUnits = "" )
pure virtual

Get block of traces.

If nTrc or nSamp exceed max values then these values are changed to max allowed (that is why they are not const)

Implemented in H5SeisImpl.

◆ getTraceD()

virtual std::optional< h5gt::DataSet > H5Seis::getTraceD ( )
pure virtual

Get trace DataSet.

Implemented in H5SeisImpl.

◆ getTraceHeader() [1/6]

virtual Eigen::MatrixXd H5Seis::getTraceHeader ( const Eigen::Ref< const Eigen::VectorX< size_t > > & trcInd,
const Eigen::Ref< const Eigen::VectorX< size_t > > & trcHdrInd,
const std::vector< std::string > & unitsFrom = std::vector< std::string >(),
const std::vector< std::string > & unitsTo = std::vector< std::string >() )
pure virtual

Get trace headers by indexes.

Return empty matrix if max trcInd exceeds nTrc.
That is done intensionally to keep trcInd size and returned number of traces in trace headers equal.

Implemented in H5SeisImpl.

◆ getTraceHeader() [2/6]

virtual Eigen::MatrixXd H5Seis::getTraceHeader ( const size_t & fromTrc,
size_t nTrc = 1,
const size_t & fromHdr = 0,
size_t nHdr = std::numeric_limits< size_t >::max(),
const std::vector< std::string > & unitsFrom = std::vector< std::string >(),
const std::vector< std::string > & unitsTo = std::vector< std::string >() )
pure virtual

Get block of trace headers.

If nTrc or nHdr exceed max values then these values are changed to max allowed (that is why they are not const)

Implemented in H5SeisImpl.

◆ getTraceHeader() [3/6]

virtual Eigen::VectorXd H5Seis::getTraceHeader ( const std::string & hdrName,
const size_t & fromTrc = 0,
size_t nTrc = 1,
const std::string & unitsFrom = "",
const std::string & unitsTo = "" )
pure virtual

Get block of trace header by name.

If nTrc exceeds max value then this value is changed to max allowed (that is why it is not const)

Implemented in H5SeisImpl.

◆ getTraceHeader() [4/6]

virtual Eigen::MatrixXd H5Seis::getTraceHeader ( const std::vector< size_t > & trcInd,
const std::vector< size_t > & trcHdrInd,
const std::vector< std::string > & unitsFrom = std::vector< std::string >(),
const std::vector< std::string > & unitsTo = std::vector< std::string >() )
pure virtual

Get trace headers by indexes.

Return empty matrix if max trcInd exceeds nTrc.
That is done intensionally to keep trcInd size and returned number of traces in trace headers equal.

Implemented in H5SeisImpl.

◆ getTraceHeader() [5/6]

virtual Eigen::MatrixXd H5Seis::getTraceHeader ( const std::vector< std::string > & hdrNames,
const Eigen::Ref< const Eigen::VectorX< size_t > > & trcInd,
const std::vector< std::string > & unitsFrom = std::vector< std::string >(),
const std::vector< std::string > & unitsTo = std::vector< std::string >() )
pure virtual

Get trace headers by indexes and names.

Return empty matrix if max trcInd exceeds nTrc.
That is done intensionally to keep trcInd size and returned number of traces in trace headers equal.

Implemented in H5SeisImpl.

◆ getTraceHeader() [6/6]

virtual Eigen::MatrixXd H5Seis::getTraceHeader ( const std::vector< std::string > & hdrNames,
const std::vector< size_t > & trcInd,
const std::vector< std::string > & unitsFrom = std::vector< std::string >(),
const std::vector< std::string > & unitsTo = std::vector< std::string >() )
pure virtual

Get trace header by indexes and name.

Return empty matrix if max trcInd exceeds nTrc.
That is done intensionally to keep trcInd size and returned number of traces in trace header equal.

Implemented in H5SeisImpl.

◆ getTraceHeaderD()

virtual std::optional< h5gt::DataSet > H5Seis::getTraceHeaderD ( )
pure virtual

Get trace header DataSet.

Implemented in H5SeisImpl.

◆ getTraceHeaderIndex()

virtual ptrdiff_t H5Seis::getTraceHeaderIndex ( const std::string & hdrName)
pure virtual

Get index (row/col within 2D dataset) for a given trace header.

Implemented in H5SeisImpl.

◆ getTraceHeaderMax() [1/2]

virtual std::map< std::string, double > H5Seis::getTraceHeaderMax ( )
pure virtual

Get trace header maximal values.

Implemented in H5SeisImpl.

◆ getTraceHeaderMax() [2/2]

virtual double H5Seis::getTraceHeaderMax ( const std::string & hdrName,
const std::string & unitsFrom = "",
const std::string & unitsTo = "" )
pure virtual

Get trace header maximal value for a given trace header.

Implemented in H5SeisImpl.

◆ getTraceHeaderMin() [1/2]

virtual std::map< std::string, double > H5Seis::getTraceHeaderMin ( )
pure virtual

Get trace header minimal values.

Implemented in H5SeisImpl.

◆ getTraceHeaderMin() [2/2]

virtual double H5Seis::getTraceHeaderMin ( const std::string & hdrName,
const std::string & unitsFrom = "",
const std::string & unitsTo = "" )
pure virtual

Get trace header minimal value for a given trace header.

Implemented in H5SeisImpl.

◆ getUValG()

virtual std::optional< h5gt::Group > H5Seis::getUValG ( )
pure virtual

Get sorting unique values Group.

Implemented in H5SeisImpl.

◆ getXYTraceHeaders() [1/2]

virtual Eigen::MatrixXd H5Seis::getXYTraceHeaders ( const std::vector< std::string > & xyHdrNames,
const Eigen::Ref< const Eigen::VectorX< size_t > > & trcInd,
const std::string & lengthUnits = "",
bool doCoordTransform = false )
pure virtual

Get XY trace headers (two columns in Eigen column-major matrix)

Same as H5Seis::getTraceHeader() but also able to do a coordinate transformation.
Return empty matrix if max trcInd exceeds nTrc.
That is done intensionally to keep trcInd size and returned number of traces in trace headers equal.

Implemented in H5SeisImpl.

◆ getXYTraceHeaders() [2/2]

virtual Eigen::MatrixXd H5Seis::getXYTraceHeaders ( const std::vector< std::string > & xyHdrNames,
const size_t & fromTrc = 0,
size_t nTrc = std::numeric_limits< size_t >::max(),
const std::string & lengthUnits = "",
bool doCoordTransform = false )
pure virtual

Get XY trace headers (two columns in Eigen column-major matrix)

Same as H5Seis::getTraceHeader() but also able to do a coordinate transformation.

Implemented in H5SeisImpl.

◆ hasPKeySort()

virtual bool H5Seis::hasPKeySort ( const std::string & pKeyName)
pure virtual

Check if PKey sort is prepared.

Implemented in H5SeisImpl.

◆ openSeisContainer()

virtual H5SeisContainer * H5Seis::openSeisContainer ( )
pure virtual

Open H5SeisContainer where current seismic resides.

Implemented in H5SeisImpl.

◆ readSEGYBinHeader()

virtual bool H5Seis::readSEGYBinHeader ( const std::string & segy,
h5geo::Endian endian = static_cast< h5geo::Endian >(0) )
pure virtual

Read binary header from SEGY file.

Implemented in H5SeisImpl.

◆ readSEGYTextHeader()

virtual bool H5Seis::readSEGYTextHeader ( const std::string & segy,
h5geo::TextEncoding encoding = static_cast< h5geo::TextEncoding >(0) )
pure virtual

Read text header from SEGY file.

Implemented in H5SeisImpl.

◆ readSEGYTraces()

virtual bool H5Seis::readSEGYTraces ( const std::vector< std::string > & segyFiles,
std::vector< h5geo::SegyFormat > formats = std::vector< h5geo::SegyFormat >(),
std::vector< h5geo::Endian > endians = std::vector< h5geo::Endian >(),
std::vector< std::vector< std::string > > trcHdrNamesArr = std::vector< std::vector< std::string > >(),
size_t trcBuffer = 10000,
std::function< void(double)> progressCallback = nullptr )
pure virtual

Read trace headers and trace data from SEGY file.

Parameters
segyFilessegy files to read traces and trace headers
formatssegy format for each segy (empty to autodefine)
endiansPC endian for each segy (empty to autodefine)
trcHdrNamesArruse only those defined in 'getTraceHeaderShortNames', but you can change their order thus fix probably messed up trace header bytes (empty to use defined in 'getTraceHeaderShortNames' func)
trcBuffernumber of traces per thread to read before writing them at once
progressCallbackcallback function of form void foo(double progress)

Implemented in H5SeisImpl.

◆ readSEGYTracesMMap()

virtual bool H5Seis::readSEGYTracesMMap ( const std::vector< std::string > & segyFiles,
std::vector< h5geo::SegyFormat > formats = std::vector< h5geo::SegyFormat >(),
std::vector< h5geo::Endian > endians = std::vector< h5geo::Endian >(),
std::vector< std::vector< std::string > > trcHdrNamesArr = std::vector< std::vector< std::string > >(),
size_t trcBuffer = 10000,
int nThreads = -1,
std::function< void(double)> progressCallback = nullptr )
pure virtual

Read trace headers and trace data from SEGY file using Memory Mapping.

Parameters
segyFilessegy files to read traces and trace headers
formatssegy format for each segy (empty to autodefine)
endiansPC endian for each segy (empty to autodefine)
trcHdrNamesArruse only those defined in 'getTraceHeaderShortNames', but you can change their order thus fix probably messed up trace header bytes (empty to use defined in 'getTraceHeaderShortNames' func)
trcBuffernumber of traces per thread to read before writing them at once
nThreadsnumber of threads (to use all threads pass any number <1)
progressCallbackcallback function of form void foo(double progress)
Note
Memory Mappings works only if the SEGY file resides on the internal hardware

Implemented in H5SeisImpl.

◆ removePKeySort()

virtual bool H5Seis::removePKeySort ( const std::string & pKeyName)
pure virtual

Remove PKey sorting.

Implemented in H5SeisImpl.

◆ setDataType()

virtual bool H5Seis::setDataType ( const h5geo::SeisDataType & seisType)
pure virtual

Set datatype for the seismic (STACK or PRESTACK)

Implemented in H5SeisImpl.

◆ setDomain()

virtual bool H5Seis::setDomain ( const h5geo::Domain & domain)
pure virtual

Set domain for the seismic (TVD, TVDSS, TWT, OWT)

Implemented in H5SeisImpl.

◆ setFirstSample()

virtual bool H5Seis::setFirstSample ( double val,
const std::string & units = "" )
pure virtual

Set first sample.

Implemented in H5SeisImpl.

◆ setNSamp()

virtual bool H5Seis::setNSamp ( size_t nSamp)
pure virtual

Resize trace DataSet.

Implemented in H5SeisImpl.

◆ setNTrc()

virtual bool H5Seis::setNTrc ( size_t nTrc)
pure virtual

Resize trace and trace header DataSets.

Implemented in H5SeisImpl.

◆ setSampRate()

virtual bool H5Seis::setSampRate ( double val,
const std::string & units = "" )
pure virtual

Set sampling rate.

Implemented in H5SeisImpl.

◆ setSRD()

virtual bool H5Seis::setSRD ( double val,
const std::string & lengthUnits = "" )
pure virtual

Set Seismic Reference Datum.

Implemented in H5SeisImpl.

◆ setSurveyType()

virtual bool H5Seis::setSurveyType ( const h5geo::SurveyType & surveyType)
pure virtual

Set survey type for the seismic (TWO_D or THREE_D)

Implemented in H5SeisImpl.

◆ updatePKeySort()

virtual bool H5Seis::updatePKeySort ( const std::string & pKeyName)
pure virtual

Update sorting for prepared PKey

Implemented in H5SeisImpl.

◆ updateTraceHeaderLimits()

virtual bool H5Seis::updateTraceHeaderLimits ( size_t nTrcBuffer = 1e7)
pure virtual

Calculate and write min/max trace headers.

Implemented in H5SeisImpl.

◆ updateTraceHeaderNSamp()

virtual bool H5Seis::updateTraceHeaderNSamp ( )
pure virtual

Set trace header number of samples from binary header.

Implemented in H5SeisImpl.

◆ updateTraceHeaderSampRate()

virtual bool H5Seis::updateTraceHeaderSampRate ( )
pure virtual

Set trace header samp rate from binary header.

Implemented in H5SeisImpl.

◆ writeBinHeader() [1/4]

virtual bool H5Seis::writeBinHeader ( const double(&) binHdr[30])
pure virtual

Write binary header.

Implemented in H5SeisImpl.

◆ writeBinHeader() [2/4]

virtual bool H5Seis::writeBinHeader ( const Eigen::Ref< const Eigen::VectorXd > & binHdrVec)
pure virtual

Write text header.

Vector length should be equal to H5Seis::getNBinHdr()

Implemented in H5SeisImpl.

◆ writeBinHeader() [3/4]

virtual bool H5Seis::writeBinHeader ( const std::string & hdrName,
const double & value,
const std::string & unitsFrom = "",
const std::string & unitsTo = "" )
pure virtual

Write text header.

Implemented in H5SeisImpl.

◆ writeBinHeader() [4/4]

virtual bool H5Seis::writeBinHeader ( const std::vector< double > & binHdrVec)
pure virtual

Write text header.

Vector length should be equal to H5Seis::getNBinHdr()

Implemented in H5SeisImpl.

◆ writeTextHeader() [1/2]

virtual bool H5Seis::writeTextHeader ( const char(&) txtHdr[40][80])
pure virtual

Write text header.

Implemented in H5SeisImpl.

◆ writeTextHeader() [2/2]

virtual bool H5Seis::writeTextHeader ( const std::vector< std::string > & txtHdr)
pure virtual

Write text header.

Maximum 40x80 chars are possible (vector of size 40 with string less or equal to 80)

Implemented in H5SeisImpl.

◆ writeTrace() [1/2]

virtual bool H5Seis::writeTrace ( Eigen::Ref< Eigen::MatrixXf > TRACE,
const Eigen::Ref< const Eigen::VectorX< size_t > > & trcInd,
const size_t & fromSampInd = 0,
const std::string & dataUnits = "" )
pure virtual

Write traces using indexes.

Return true even if max trcInd exceeds nTrc.

Implemented in H5SeisImpl.

◆ writeTrace() [2/2]

virtual bool H5Seis::writeTrace ( Eigen::Ref< Eigen::MatrixXf > TRACE,
const size_t & fromTrc = 0,
const size_t & fromSampInd = 0,
const std::string & dataUnits = "" )
pure virtual

Write block of traces starting from trace fromTrc and from sample fromSampInd

Implemented in H5SeisImpl.

◆ writeTraceHeader() [1/3]

virtual bool H5Seis::writeTraceHeader ( const Eigen::Ref< const Eigen::MatrixXd > & HDR,
const size_t & fromTrc = 0,
const size_t & fromHdrInd = 0 )
pure virtual

Write block of trace headers starting from trace fromTrc and from header index fromHdrInd

Implemented in H5SeisImpl.

◆ writeTraceHeader() [2/3]

virtual bool H5Seis::writeTraceHeader ( const std::string & hdrName,
Eigen::Ref< Eigen::MatrixXd > hdr,
const Eigen::Ref< const Eigen::VectorX< size_t > > & trcInd,
const std::string & unitsFrom = "",
const std::string & unitsTo = "" )
pure virtual

Write trace header by name and trace indexes.

Return true even if max trcInd exceeds nTrc.

Implemented in H5SeisImpl.

◆ writeTraceHeader() [3/3]

virtual bool H5Seis::writeTraceHeader ( const std::string & hdrName,
Eigen::Ref< Eigen::MatrixXd > hdr,
const size_t & fromTrc = 0,
const std::string & unitsFrom = "",
const std::string & unitsTo = "" )
pure virtual

Write trace header by name.

Implemented in H5SeisImpl.

◆ writeXYTraceHeaders() [1/2]

virtual bool H5Seis::writeXYTraceHeaders ( const std::vector< std::string > & xyHdrNames,
Eigen::Ref< Eigen::MatrixX2d > & xy,
const Eigen::Ref< const Eigen::VectorX< size_t > > & trcInd,
const std::string & lengthUnits = "",
bool doCoordTransform = false )
pure virtual

Write XY trace headers (two columns in Eigen column-major matrix)

Same as H5Seis::writeTraceHeaders() but also able to do a coordinate transformation.
Return true even if max trcInd exceeds nTrc.

Implemented in H5SeisImpl.

◆ writeXYTraceHeaders() [2/2]

virtual bool H5Seis::writeXYTraceHeaders ( const std::vector< std::string > & xyHdrNames,
Eigen::Ref< Eigen::MatrixX2d > & xy,
const size_t & fromTrc = 0,
const std::string & lengthUnits = "",
bool doCoordTransform = false )
pure virtual

Write XY trace headers (two columns in Eigen column-major matrix)

Same as H5Seis::writeTraceHeaders() but also able to do a coordinate transformation.

Implemented in H5SeisImpl.


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