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
h5vol.h
1#ifndef H5VOL_H
2#define H5VOL_H
3
4#include "h5baseobject.h"
5
6#include <Eigen/Dense>
7
9
14class H5Vol : public H5BaseObject
15{
16protected:
17 virtual ~H5Vol() = default;
18
19public:
20
23 virtual bool writeData(
24 Eigen::Ref<Eigen::MatrixXf> data,
25 const size_t& iX0,
26 const size_t& iY0,
27 const size_t& iZ0,
28 const size_t& nX,
29 const size_t& nY,
30 const size_t& nZ,
31 const std::string& dataUnits = "") = 0;
32
51 virtual bool readSEGYSTACK(
52 const std::string& segy,
53 const size_t& ilHdrOffset,
54 const size_t& ilHdrSize,
55 const size_t& xlHdrOffset,
56 const size_t& xlHdrSize,
57 const size_t& xHdrOffset,
58 const size_t& xHdrSize,
59 const size_t& yHdrOffset,
60 const size_t& yHdrSize,
61 double sampRate,
62 size_t nSamp = 0,
63 size_t nTrc = 0,
64 h5geo::SegyFormat format = static_cast<h5geo::SegyFormat>(0),
65 h5geo::Endian endian = static_cast<h5geo::Endian>(0),
66 std::function<void(double)> progressCallback = nullptr) = 0;
67
69 virtual bool setDomain(const h5geo::Domain& domain) = 0;
71 virtual bool setOrigin(
72 Eigen::Ref<Eigen::Vector3d> v,
73 const std::string& lengthUnits = "",
74 const std::string& temporalUnits = "",
75 bool doCoordTransform = false) = 0;
77 virtual bool setSpacings(
78 Eigen::Ref<Eigen::Vector3d> v,
79 const std::string& lengthUnits = "",
80 const std::string& temporalUnits = "") = 0;
82 virtual bool setOrientation(
83 double val,
84 const std::string& angularUnits = "") = 0;
85
87 virtual bool resize(
88 size_t nx, size_t ny, size_t nz) = 0;
89
92 virtual Eigen::MatrixXf getData(
93 const size_t& iX0,
94 const size_t& iY0,
95 const size_t& iZ0,
96 const size_t& nX,
97 const size_t& nY,
98 const size_t& nZ,
99 const std::string& dataUnits = "") = 0;
100
102 virtual h5geo::Domain getDomain() = 0;
104 virtual Eigen::VectorXd getOrigin(
105 const std::string& lengthUnits = "",
106 const std::string& temporalUnits = "",
107 bool doCoordTransform = false) = 0;
109 virtual Eigen::VectorXd getSpacings(
110 const std::string& lengthUnits = "",
111 const std::string& temporalUnits = "") = 0;
113 virtual double getOrientation(
114 const std::string& angularUnits = "") = 0;
116 virtual size_t getNX() = 0;
118 virtual size_t getNY() = 0;
120 virtual size_t getNZ() = 0;
121
123 virtual H5VolParam getParam() = 0;
124
126 virtual H5VolContainer* openVolContainer() const = 0;
127
129 virtual std::optional<h5gt::DataSet> getVolD() const = 0;
130
131 virtual bool exportToSEGY(
132 const std::string& segyFile,
133 h5geo::Endian endian = h5geo::Endian::Big,
134 std::function<void(double)> progressCallback = nullptr) = 0;
135
137 virtual bool recreateVolD(
138 size_t nX, size_t nY, size_t nZ,
139 size_t xChunk, size_t yChunk, size_t zChunk,
140 unsigned compressionLevel) = 0;
141};
142
143using H5Vol_ptr = std::unique_ptr<H5Vol, h5geo::ObjectDeleter>;
144
145#endif // H5VOL_H
Base class for geo-objects.
Definition h5baseobject.h:13
A container built around HDF5 file and used for storing and manipulating H5Vol objects.
Definition h5volcontainer.h:15
Provides API to work with volumes.
Definition h5vol.h:15
virtual size_t getNZ()=0
Get number of Z values.
virtual Eigen::VectorXd getOrigin(const std::string &lengthUnits="", const std::string &temporalUnits="", bool doCoordTransform=false)=0
Get coordinates of origin.
virtual size_t getNX()=0
Get number of X values.
virtual double getOrientation(const std::string &angularUnits="")=0
Get XY plane orientation.
virtual Eigen::MatrixXf getData(const size_t &iX0, const size_t &iY0, const size_t &iZ0, const size_t &nX, const size_t &nY, const size_t &nZ, const std::string &dataUnits="")=0
Read subvolume starting from iX0, iY0, iZ0 indices. data matrix is of size: nRows=nX,...
virtual Eigen::VectorXd getSpacings(const std::string &lengthUnits="", const std::string &temporalUnits="")=0
Get X,Y,Z unrotated spacings.
virtual h5geo::Domain getDomain()=0
Get domain (TVD, TVDSS, TWT, OWT)
virtual bool setDomain(const h5geo::Domain &domain)=0
Set domain for the map (TVD, TVDSS, TWT, OWT)
virtual H5VolContainer * openVolContainer() const =0
Open H5VolContainer where current vol resides.
virtual std::optional< h5gt::DataSet > getVolD() const =0
Get current vol's DataSet.
virtual bool resize(size_t nx, size_t ny, size_t nz)=0
Resize volume.
virtual size_t getNY()=0
Get number of Y values.
virtual bool readSEGYSTACK(const std::string &segy, const size_t &ilHdrOffset, const size_t &ilHdrSize, const size_t &xlHdrOffset, const size_t &xlHdrSize, const size_t &xHdrOffset, const size_t &xHdrSize, const size_t &yHdrOffset, const size_t &yHdrSize, double sampRate, size_t nSamp=0, size_t nTrc=0, h5geo::SegyFormat format=static_cast< h5geo::SegyFormat >(0), h5geo::Endian endian=static_cast< h5geo::Endian >(0), std::function< void(double)> progressCallback=nullptr)=0
Read SEGY STACK data, i.e. nTrc should be equal to nil*nxl. After reading origin, spacings,...
virtual bool recreateVolD(size_t nX, size_t nY, size_t nZ, size_t xChunk, size_t yChunk, size_t zChunk, unsigned compressionLevel)=0
Unlink and create new dataset without copying data.
virtual H5VolParam getParam()=0
Get parameters that were used to create current map.
virtual bool writeData(Eigen::Ref< Eigen::MatrixXf > data, const size_t &iX0, const size_t &iY0, const size_t &iZ0, const size_t &nX, const size_t &nY, const size_t &nZ, const std::string &dataUnits="")=0
Write subvolume starting from iX0, iY0, iZ0 indices. data matrix is of size: nRows=nX,...
virtual bool setOrientation(double val, const std::string &angularUnits="")=0
Set XY plane orientation.
virtual bool setOrigin(Eigen::Ref< Eigen::Vector3d > v, const std::string &lengthUnits="", const std::string &temporalUnits="", bool doCoordTransform=false)=0
Set coordinates of origin.
virtual bool setSpacings(Eigen::Ref< Eigen::Vector3d > v, const std::string &lengthUnits="", const std::string &temporalUnits="")=0
Set X,Y,Z unrotated spacings.
Class for creating H5Vol.
Definition h5base.h:116