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
h5map.h
1#ifndef H5MAP_H
2#define H5MAP_H
3
4#include "h5baseobject.h"
5
6#include <Eigen/Dense>
7
9
18class H5Map : public H5BaseObject
19{
20protected:
21 virtual ~H5Map() = default;
22
23public:
24
25 #ifdef H5GEO_USE_GDAL
30 virtual bool readRasterCoordinates(
31 const std::string& file,
32 const std::string& lengthUnits = "") = 0;
36 virtual bool readRasterSpatialReference(const std::string& file) = 0;
40 virtual bool readRasterLengthUnits(const std::string& file) = 0;
45 virtual bool readRasterData(
46 const std::string& file,
47 const std::string& dataUnits = "") = 0;
48 #endif // H5GEO_USE_GDAL
49
51 virtual bool writeData(
52 Eigen::Ref<Eigen::MatrixXd> M,
53 const std::string& dataUnits = "") = 0;
54
56 virtual Eigen::MatrixXd getData(const std::string& dataUnits = "") = 0;
57
59 virtual bool setDomain(const h5geo::Domain& domain) = 0;
61 virtual bool setOrigin(
62 Eigen::Ref<Eigen::Vector2d> v,
63 const std::string& lengthUnits = "",
64 bool doCoordTransform = false) = 0;
66 virtual bool setPoint1(
67 Eigen::Ref<Eigen::Vector2d> v,
68 const std::string& lengthUnits = "",
69 bool doCoordTransform = false) = 0;
71 virtual bool setPoint2(
72 Eigen::Ref<Eigen::Vector2d> v,
73 const std::string& lengthUnits = "",
74 bool doCoordTransform = false) = 0;
75
85 virtual std::optional<h5gt::Group> addAttributeMap(H5Map* map, std::string name = "") = 0;
93 virtual std::optional<h5gt::Group> addInternalAttributeMap(H5Map* map, std::string name = "") = 0;
101 virtual std::optional<h5gt::Group> addExternalAttributeMap(H5Map* map, std::string name = "") = 0;
103 virtual bool removeAttributeMap(const std::string& name) = 0;
105 virtual H5Map* openAttributeMap(const std::string& name) = 0;
107 virtual std::vector<h5gt::Group> getAttributeMapGroupList() = 0;
109 virtual std::vector<std::string> getAttributeMapNameList() = 0;
111 virtual size_t getAttributeMapCount() = 0;
112
114 virtual h5geo::Domain getDomain() = 0;
116 virtual Eigen::VectorXd getOrigin(
117 const std::string& lengthUnits = "",
118 bool doCoordTransform = false) = 0;
120 virtual Eigen::VectorXd getPoint1(
121 const std::string& lengthUnits = "",
122 bool doCoordTransform = false) = 0;
124 virtual Eigen::VectorXd getPoint2(
125 const std::string& lengthUnits = "",
126 bool doCoordTransform = false) = 0;
128 virtual size_t getNX() = 0;
130 virtual size_t getNY() = 0;
131
133 virtual H5MapParam getParam() = 0;
134
136 virtual H5MapContainer* openMapContainer() const = 0;
137
139 virtual std::optional<h5gt::DataSet> getMapD() const = 0;
140};
141
142using H5Map_ptr = std::unique_ptr<H5Map, h5geo::ObjectDeleter>;
143
144#endif // H5MAP_H
Base class for geo-objects.
Definition h5baseobject.h:13
A container built around HDF5 file and used for storing and manipulating H5Map objects.
Definition h5mapcontainer.h:15
Provides API to work with maps.
Definition h5map.h:19
virtual bool setDomain(const h5geo::Domain &domain)=0
Set domain for the map (TVD, TVDSS, TWT, OWT)
virtual std::optional< h5gt::Group > addExternalAttributeMap(H5Map *map, std::string name="")=0
Add internal attribute map.
virtual bool setPoint1(Eigen::Ref< Eigen::Vector2d > v, const std::string &lengthUnits="", bool doCoordTransform=false)=0
Set coordinates of upper-right matrix corner.
virtual bool setPoint2(Eigen::Ref< Eigen::Vector2d > v, const std::string &lengthUnits="", bool doCoordTransform=false)=0
Set coordinates of bottom-left matrix corner.
virtual std::vector< h5gt::Group > getAttributeMapGroupList()=0
Find and return vector of attribute maps Groups within current map.
virtual std::optional< h5gt::Group > addInternalAttributeMap(H5Map *map, std::string name="")=0
Add internal attribute map.
virtual bool setOrigin(Eigen::Ref< Eigen::Vector2d > v, const std::string &lengthUnits="", bool doCoordTransform=false)=0
Set coordinates of upper-left matrix corner.
virtual Eigen::VectorXd getPoint1(const std::string &lengthUnits="", bool doCoordTransform=false)=0
Get coordinates of upper-right matrix corner.
virtual Eigen::MatrixXd getData(const std::string &dataUnits="")=0
Read data from DataSet.
virtual std::optional< h5gt::DataSet > getMapD() const =0
Get current map's DataSet.
virtual bool removeAttributeMap(const std::string &name)=0
Remove attribute map.
virtual H5MapParam getParam()=0
Get parameters that were used to create current map.
virtual size_t getNX()=0
Get number of Eigen column-major matrix columns.
virtual std::optional< h5gt::Group > addAttributeMap(H5Map *map, std::string name="")=0
Add attribute map.
virtual Eigen::VectorXd getOrigin(const std::string &lengthUnits="", bool doCoordTransform=false)=0
Get coordinates of upper-left matrix corner.
virtual Eigen::VectorXd getPoint2(const std::string &lengthUnits="", bool doCoordTransform=false)=0
Get coordinates of bottom-left matrix corner.
virtual size_t getAttributeMapCount()=0
Get number of atribute maps within current map.
virtual size_t getNY()=0
Get number of Eigen column-major matrix rows.
virtual H5MapContainer * openMapContainer() const =0
Open H5MapContainer where current map resides.
virtual H5Map * openAttributeMap(const std::string &name)=0
Open attribute map.
virtual std::vector< std::string > getAttributeMapNameList()=0
Find and return vector of attribute maps names within current map.
virtual bool writeData(Eigen::Ref< Eigen::MatrixXd > M, const std::string &dataUnits="")=0
Write data to DataSet.
virtual h5geo::Domain getDomain()=0
Get domain (TVD, TVDSS, TWT, OWT)
Class for creating H5Map.
Definition h5base.h:98