kspaceFirstOrder3D-CUDA  1.1
The CUDA/C++ implementation of the k-wave toolbox for the time-domain simulation of acoustic wave fields in 3D
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
BaseOutputHDF5Stream.h
Go to the documentation of this file.
1 /**
2  * @file BaseOutputHDF5Stream.h
3  *
4  * @author Jiri Jaros \n
5  * Faculty of Information Technology \n
6  * Brno University of Technology \n
7  * jarosjir@fit.vutbr.cz
8  *
9  * @brief The header file of the class saving RealMatrix data into the output HDF5 file.
10  *
11  * @version kspaceFirstOrder3D 3.4
12  *
13  * @date 11 July 2012, 10:30 (created) \n
14  * 26 July 2016, 12:35 (revised)
15  *
16  * @section License
17  * This file is part of the C++ extension of the k-Wave Toolbox
18  * (http://www.k-wave.org).\n Copyright (C) 2016 Jiri Jaros and Bradley Treeby.
19  *
20  * This file is part of the k-Wave. k-Wave is free software: you can redistribute it and/or modify
21  * it under the terms of the GNU Lesser General Public License as published by the Free Software
22  * Foundation, either version 3 of the License, or (at your option) any later version.
23  *
24  * k-Wave is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even
25  * the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser
26  * General Public License for more details.
27  *
28  * You should have received a copy of the GNU Lesser General Public License along with k-Wave.
29  * If not, see http://www.gnu.org/licenses/.
30  */
31 
32 #ifndef BASE_OUTPUT_HDF5_STREAM_H
33 #define BASE_OUTPUT_HDF5_STREAM_H
34 
35 
38 #include <HDF5/HDF5_File.h>
39 
40 /**
41  * @class TBaseOutputHDF5Stream
42  * @brief Abstract base class for output data streams (sampled data).
43  * @details Abstract base class for output data streams (sampled data).
44  *
45  */
47 {
48  public:
49  /**
50  * @enum TReduceOperator
51  * @brief How to aggregate data.
52  * @brief How to aggregate data.
53  */
55  {
56  /// store actual data (time series)
58  /// calculate root mean square
59  RMS,
60  /// store maximum
61  MAX,
62  /// store minimum
64  };
65 
66  /// Constructor.
71 
72  /// Destructor
73  virtual ~TBaseOutputHDF5Stream() {};
74 
75  /// Create a HDF5 stream and allocate data for it.
76  virtual void Create() = 0;
77 
78  /// Reopen the output stream after restart.
79  virtual void Reopen() = 0;
80 
81  /// Sample data into buffer, apply reduction - based on a sensor mask (no data flushed to disk).
82  virtual void Sample() = 0;
83 
84  /// Flush raw data to disk.
85  virtual void FlushRaw() = 0;
86 
87  /// Apply post-processing on the buffer and flush it to the file.
88  virtual void PostProcess();
89 
90  /// Checkpoint the stream.
91  virtual void Checkpoint() = 0;
92 
93  /// Close stream (apply post-processing if necessary, flush data and close).
94  virtual void Close() = 0;
95 
96  protected:
97  /// Default constructor not allowed.
99  /// Copy constructor not allowed.
101  /// Operator = not allowed (we don't want any data movements).
103 
104  /// A generic function to allocate memory - not used in the base class.
105  virtual void AllocateMemory();
106  /// A generic function to free memory - not used in the base class.
107  virtual void FreeMemory();
108 
109  /// Copy data hostBuffer-> deviceBuffer
110  virtual void CopyDataToDevice();
111  /// Copy data deviceBuffer -> hostBuffer
112  virtual void CopyDataFromDevice();
113 
114  /// HDF5 file handle.
116 
117  /// Dataset name.
118  std::string rootObjectName;
119 
120  /// Source matrix to be sampled.
122 
123  /// HDF5 dataset handle.
124  hid_t dataset;
125  /// Reduce operator.
127 
128  /// Position in the dataset
130 
131  /// Buffer size
132  size_t bufferSize;
133 
134  /// Temporary buffer for store on the GPU side
135  float* hostBuffer;
136  /// Temporary buffer on the GPU side - only for aggregated quantities
137  float* deviceBuffer;
138 
139  /// chunk size of 4MB in number of float elements
140  static const size_t CHUNK_SIZE_4MB = 1048576;
141  /// The minimum number of elements to start sampling in parallel (4MB)
142  static const size_t MIN_GRIDPOINTS_TO_SAMPLE_IN_PARALLEL = 1048576;
143 
144 };// end of TBaseOutputHDF5Stream
145 //--------------------------------------------------------------------------------------------------
146 #endif /* BASE_OUTPUT_HDF5_STREAM_H */
147 
THDF5_File & file
HDF5 file handle.
TBaseOutputHDF5Stream()
Default constructor not allowed.
TBaseOutputHDF5Stream & operator=(const TBaseOutputHDF5Stream &src)
Operator = not allowed (we don't want any data movements).
virtual void CopyDataFromDevice()
Copy data deviceBuffer -> hostBuffer.
The header file containing the class for real matrices.
calculate root mean square
std::string rootObjectName
Dataset name.
virtual void Close()=0
Close stream (apply post-processing if necessary, flush data and close).
const TReduceOperator reduceOp
Reduce operator.
The header file containing the HDF5 related classes.
virtual void Sample()=0
Sample data into buffer, apply reduction - based on a sensor mask (no data flushed to disk)...
static const size_t CHUNK_SIZE_4MB
chunk size of 4MB in number of float elements
const std::string TMatrixName
Datatype for matrix names.
Definition: MatrixNames.h:45
virtual void Reopen()=0
Reopen the output stream after restart.
virtual void PostProcess()
Apply post-processing on the buffer and flush it to the file.
hid_t dataset
HDF5 dataset handle.
virtual void AllocateMemory()
A generic function to allocate memory - not used in the base class.
const TRealMatrix & sourceMatrix
Source matrix to be sampled.
The header file containing the class for 64b integer matrices.
virtual void Create()=0
Create a HDF5 stream and allocate data for it.
TReduceOperator
How to aggregate data.
float * deviceBuffer
Temporary buffer on the GPU side - only for aggregated quantities.
store actual data (time series)
Abstract base class for output data streams (sampled data).
virtual void FreeMemory()
A generic function to free memory - not used in the base class.
The class for real matrices.
Definition: RealMatrix.h:45
float * hostBuffer
Temporary buffer for store on the GPU side.
virtual void CopyDataToDevice()
Copy data hostBuffer-> deviceBuffer.
TDimensionSizes position
Position in the dataset.
virtual void FlushRaw()=0
Flush raw data to disk.
size_t bufferSize
Buffer size.
virtual ~TBaseOutputHDF5Stream()
Destructor.
virtual void Checkpoint()=0
Checkpoint the stream.
static const size_t MIN_GRIDPOINTS_TO_SAMPLE_IN_PARALLEL
The minimum number of elements to start sampling in parallel (4MB)
Class wrapping the HDF5 routines.
Definition: HDF5_File.h:500
Structure with 4D dimension sizes (3 in space and 1 in time).