5#include "h5geo/h5core.h"
17template <
class To,
class From>
18typename std::enable_if<(
sizeof(To) ==
sizeof(From)) &&
19std::is_trivially_copyable<From>::value &&
20std::is_trivial<To>::value,
25inline bit_cast(
const From &src)
noexcept {
27 std::memcpy(&dst, &src,
sizeof(To));
31inline float ibm2ieee(
const int &from) {
32 std::bitset<32> bits_from(from);
34 const int constnum1 = 2130706432;
35 const int constnum2 = 16777215;
36 const float constnum3 = 16777216;
38 int sgn = (1 - 2 * bits_from[31]);
39 int expon = (from & constnum1) >> 24;
40 int mantiss = (from & constnum2);
42 to = float(sgn) * float(pow(16, expon - 64)) * float(mantiss) / constnum3;
47inline unsigned char ebc_to_ascii_table(
unsigned char ascii) {
50 int Tableebc2Asc[256] = {
51 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32,
52 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32,
53 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32,
54 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32,
55 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 32, 91,
56 46, 60, 40, 43, 33, 38, 32, 32, 32, 32, 32, 32, 32, 32, 32,
57 93, 36, 42, 41, 59, 94, 45, 47, 32, 32, 32, 32, 32, 32, 32,
58 32, 124, 44, 37, 95, 62, 63, 32, 32, 32, 32, 32, 32, 238, 160,
59 161, 96, 58, 35, 64, 39, 61, 34, 230, 97, 98, 99, 100, 101, 102,
60 103, 104, 105, 164, 165, 228, 163, 229, 168, 169, 106, 107, 108, 109, 110,
61 111, 112, 113, 114, 170, 171, 172, 173, 174, 175, 239, 126, 115, 116, 117,
62 118, 119, 120, 121, 122, 224, 225, 226, 227, 166, 162, 236, 235, 167, 232,
63 237, 233, 231, 234, 158, 128, 129, 150, 132, 133, 148, 131, 123, 65, 66,
64 67, 68, 69, 70, 71, 72, 73, 149, 136, 137, 138, 139, 140, 125, 74,
65 75, 76, 77, 78, 79, 80, 81, 82, 141, 142, 143, 159, 144, 145, 92,
66 32, 83, 84, 85, 86, 87, 88, 89, 90, 146, 147, 134, 130, 156, 155,
67 48, 49, 50, 51, 52, 53, 54, 55, 56, 57, 135, 152, 157, 153, 151,
70 asc = Tableebc2Asc[ascii];
76inline T bswap(T val) {
78 char *pVal = (
char*) &val;
79 char *pRetVal = (
char*)&retVal;
81 for(
int i=0; i<size; i++) {
82 pRetVal[size-1-i] = pVal[i];
92template <
typename Scalar,
typename II,
typename OI>
93inline void bswap(II begin, II end, OI dest) {
94 std::transform(begin, end, dest, bswap<Scalar>);
98inline T to_big_endian(T val, h5geo::Endian fromEndian) {
99 if (fromEndian == h5geo::Endian::Little){
105template <
typename Scalar,
typename II,
typename OI>
106inline void to_big_endian(II begin, II end, OI dest, h5geo::Endian fromEndian) {
107 if (fromEndian == h5geo::Endian::Little){
108 bswap<Scalar>(begin, end, dest);
113inline T to_little_endian(T val, h5geo::Endian fromEndian) {
114 if (fromEndian == h5geo::Endian::Big){
120template <
typename Scalar,
typename II,
typename OI>
121inline void to_little_endian(II begin, II end, OI dest, h5geo::Endian fromEndian) {
122 if (fromEndian == h5geo::Endian::Big){
123 bswap<Scalar>(begin, end, dest);
128inline T to_native_endian(T val, h5geo::Endian fromEndian) {
129 if (O32_HOST_ORDER == O32_LITTLE_ENDIAN){
130 if (fromEndian == h5geo::Endian::Big){
133 }
else if (O32_HOST_ORDER == O32_BIG_ENDIAN){
134 if (fromEndian == h5geo::Endian::Little){
141template <
typename Scalar,
typename II,
typename OI>
142inline void to_native_endian(II begin, II end, OI dest, h5geo::Endian fromEndian) {
143 if (O32_HOST_ORDER == O32_LITTLE_ENDIAN){
144 if (fromEndian == h5geo::Endian::Big){
145 bswap<Scalar>(begin, end, dest);
147 }
else if (O32_HOST_ORDER == O32_BIG_ENDIAN){
148 if (fromEndian == h5geo::Endian::Little){
149 bswap<Scalar>(begin, end, dest);
154H5GEO_EXPORT
bool isSEGY(
const std::string& segy);
156H5GEO_EXPORT TextEncoding getSEGYTextEncoding(
const std::string& segy);
158H5GEO_EXPORT Endian getSEGYEndian(
const std::string& segy);
160H5GEO_EXPORT SegyFormat getSEGYFormat(
161 const std::string& segy, h5geo::Endian endian =
static_cast<h5geo::Endian
>(0));
163H5GEO_EXPORT
bool readSEGYTextHeader(
164 const std::string& segy,
165 char txtHdr[40][80], h5geo::TextEncoding encoding =
static_cast<h5geo::TextEncoding
>(0));
173 const std::string& segy,
174 char txtHdr[40][80],
bool truncate);
176H5GEO_EXPORT
bool readSEGYBinHeader(
177 const std::string& segy,
178 double binHdr[30], h5geo::Endian endian =
static_cast<h5geo::Endian
>(0));
186 const std::string& segy,
187 double binHdr[30],
bool truncate,
188 h5geo::Endian endian = h5geo::Endian::Big);
190H5GEO_EXPORT
double getSEGYSampRate(
191 const std::string& segy, h5geo::Endian endian =
static_cast<h5geo::Endian
>(0));
193H5GEO_EXPORT
size_t getSEGYNSamp(
194 const std::string& segy, h5geo::Endian endian =
static_cast<h5geo::Endian
>(0));
196H5GEO_EXPORT
size_t getSEGYNTrc(
197 const std::string& segy,
size_t nSamp = 0, h5geo::Endian endian =
static_cast<h5geo::Endian
>(0));
211 const std::string& segy,
212 const size_t& hdrOffset,
213 const size_t& hdrSize,
215 size_t toTrc = std::numeric_limits<size_t>::max(),
218 h5geo::Endian endian =
static_cast<h5geo::Endian
>(0),
219 std::function<
void(
double)> progressCallback =
nullptr);
232 h5geo::SegyFormat format,
233 h5geo::Endian endian,
234 Eigen::Ref<Eigen::VectorXf> trace);
242 const std::string& segy,
243 Eigen::Ref<Eigen::MatrixXd> HDR,
244 Eigen::Ref<Eigen::MatrixXf> TRACE,
245 h5geo::Endian endian = h5geo::Endian::Big);
260 const std::string& segy,
262 size_t toSamp = std::numeric_limits<size_t>::max(),
264 size_t toTrc = std::numeric_limits<size_t>::max(),
267 h5geo::SegyFormat format =
static_cast<h5geo::SegyFormat
>(0),
268 h5geo::Endian endian =
static_cast<h5geo::Endian
>(0),
269 std::function<
void(
double)> progressCallback =
nullptr);
291 const std::string& segy,
292 bool appendTraces =
false,
295 h5geo::SegyFormat format =
static_cast<h5geo::SegyFormat
>(0),
296 h5geo::Endian endian =
static_cast<h5geo::Endian
>(0),
297 std::vector<std::string> trcHdrNames = std::vector<std::string>(),
298 size_t trcBuffer = 10000,
300 std::function<
void(
double)> progressCallback =
nullptr);
320 const std::string& segy,
321 bool appendTraces =
false,
324 h5geo::SegyFormat format =
static_cast<h5geo::SegyFormat
>(0),
325 h5geo::Endian endian =
static_cast<h5geo::Endian
>(0),
326 std::vector<std::string> trcHdrNames = std::vector<std::string>(),
327 size_t trcBuffer = 10000,
328 std::function<
void(
double)> progressCallback =
nullptr);
351 const std::string& segy,
352 const size_t& ilHdrOffset,
353 const size_t& ilHdrSize,
354 const size_t& xlHdrOffset,
355 const size_t& xlHdrSize,
356 const size_t& xHdrOffset,
357 const size_t& xHdrSize,
358 const size_t& yHdrOffset,
359 const size_t& yHdrSize,
363 h5geo::SegyFormat format =
static_cast<h5geo::SegyFormat
>(0),
364 h5geo::Endian endian =
static_cast<h5geo::Endian
>(0),
365 std::function<
void(
double)> progressCallback =
nullptr);
Provides API to work with seismic.
Definition h5seis.h:33
Provides API to work with volumes.
Definition h5vol.h:15
Basic namespace.
Definition h5base.h:29
H5GEO_EXPORT Eigen::MatrixXf readSEGYTraces(const std::string &segy, size_t fromSamp=0, size_t toSamp=std::numeric_limits< size_t >::max(), size_t fromTrc=0, size_t toTrc=std::numeric_limits< size_t >::max(), 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)
readSEGYTraces read traces
Definition h5core_segy.cpp:569
H5GEO_EXPORT Eigen::VectorX< ptrdiff_t > readSEGYTraceHeader(const std::string &segy, const size_t &hdrOffset, const size_t &hdrSize, size_t fromTrc=0, size_t toTrc=std::numeric_limits< size_t >::max(), size_t nSamp=0, size_t nTrc=0, h5geo::Endian endian=static_cast< h5geo::Endian >(0), std::function< void(double)> progressCallback=nullptr)
readSEGYTraceHeader read selected header from all the traces (4 byte SEGY supported only)
Definition h5core_segy.cpp:359
H5GEO_EXPORT bool writeSEGYTextHeader(const std::string &segy, char txtHdr[40][80], bool truncate)
create new file or open existing SEGY
Definition h5core_segy.cpp:175
H5GEO_EXPORT bool writeSEGYTraces(const std::string &segy, Eigen::Ref< Eigen::MatrixXd > HDR, Eigen::Ref< Eigen::MatrixXf > TRACE, h5geo::Endian endian=h5geo::Endian::Big)
writeSEGYTraces write traces and their headers to the end of SEGY file
Definition h5core_segy.cpp:479
H5GEO_EXPORT void readSEGYTrace(std::ifstream &file, size_t trcInd, h5geo::SegyFormat format, h5geo::Endian endian, Eigen::Ref< Eigen::VectorXf > trace)
low level api. No any checks are done. User is responsible for that.
Definition h5core_segy.cpp:451
H5GEO_EXPORT bool readSEGYTracesMMap(H5Seis *seis, const std::string &segy, bool appendTraces=false, 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::vector< std::string > trcHdrNames=std::vector< std::string >(), size_t trcBuffer=10000, int nThreads=-1, std::function< void(double)> progressCallback=nullptr)
readSEGYTracesMMap read and write SEGY traces and trace headers to H5Seis object using Memory Mapping...
Definition h5core_segy.cpp:683
H5GEO_EXPORT bool readSEGYSTACK(H5Vol *vol, 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)
Read SEGY STACK data, i.e. nTrc should be equal to nil*nxl. After reading origin, spacings,...
Definition h5core_segy.cpp:1066
H5GEO_EXPORT bool writeSEGYBinHeader(const std::string &segy, double binHdr[30], bool truncate, h5geo::Endian endian=h5geo::Endian::Big)
create new file or open existing SEGY
Definition h5core_segy.cpp:258