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
IndexOutputHDF5Stream.h
Go to the documentation of this file.
1 /**
2  * @file IndexOutputHDF5Stream.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 data based on the index senor mask into the
10  * output HDF5 file.
11  *
12  * @version kspaceFirstOrder3D 3.4
13  *
14  * @date 28 August 2014, 10:00 (created)
15  * 26 July 2016, 13:55 (revised)
16  *
17  * @section License
18  * This file is part of the C++ extension of the k-Wave Toolbox
19  * (http://www.k-wave.org).\n Copyright (C) 2016 Jiri Jaros and Bradley Treeby.
20  *
21  * This file is part of the k-Wave. k-Wave is free software: you can redistribute it and/or modify
22  * it under the terms of the GNU Lesser General Public License as published by the Free Software
23  * Foundation, either version 3 of the License, or (at your option) any later version.
24  *
25  * k-Wave is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even
26  * the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser
27  * General Public License for more details.
28  *
29  * You should have received a copy of the GNU Lesser General Public License along with k-Wave.
30  * If not, see http://www.gnu.org/licenses/.
31  */
32 #ifndef INDEX_OUTPUT_HDF5_STREAM_H
33 #define INDEX_OUTPUT_HDF5_STREAM_H
34 
35 #include <cuda_runtime.h>
36 
38 
39 /**
40  * @class TIndexOutputHDF5Stream.
41  * @brief Output stream for quantities sampled by an index sensor mask.
42  * @details Output stream for quantities sampled by an index sensor mask.
43  * This class writes data to a single dataset in a root group of the HDF5 file (time-series
44  * as well as aggregations).
45  *
46  */
48 {
49  public:
50 
51  /// Constructor.
53  TMatrixName& datasetName,
55  const TIndexMatrix& sensorMask,
57 
58 
59  /// Destructor.
60  virtual ~TIndexOutputHDF5Stream();
61 
62  /// Create a HDF5 stream and allocate data for it.
63  virtual void Create();
64 
65  /// Reopen the output stream after restart and reload data.
66  virtual void Reopen();
67 
68  /// Sample data into buffer, apply reduction or copy to the CPU side
69  virtual void Sample();
70 
71  /// Flush data to disk (from raw streams only).
72  virtual void FlushRaw();
73 
74  /// Apply post-processing on the buffer and flush it to the file.
75  virtual void PostProcess();
76 
77  /// Checkpoint the stream.
78  virtual void Checkpoint();
79 
80  /// Close stream (apply post-processing if necessary, flush data and close).
81  virtual void Close();
82 
83  protected:
84 
85  /// Flush the buffer to the file.
86  virtual void FlushBufferToFile();
87 
88  /// Sensor mask to sample data.
90 
91  /// Handle to a HDF5 dataset.
92  hid_t dataset;
93 
94  /// Time step to store (N/A for aggregated)
96 
97  /// Has the sampling finished?
98  cudaEvent_t eventSamplingFinished;
99 
100 } ; // end of TIndexOutputHDF5Stream
101 //------------------------------------------------------------------------------
102 
103 #endif /* INDEX_OUTPUT_HDF5_STREAM_H */
104 
THDF5_File & file
HDF5 file handle.
The header file of the class saving RealMatrix data into the output HDF5 file.
const TReduceOperator reduceOp
Reduce operator.
virtual void Reopen()
Reopen the output stream after restart and reload data.
virtual void Sample()
Sample data into buffer, apply reduction or copy to the CPU side.
const std::string TMatrixName
Datatype for matrix names.
Definition: MatrixNames.h:45
const TRealMatrix & sourceMatrix
Source matrix to be sampled.
TIndexOutputHDF5Stream(THDF5_File &file, TMatrixName &datasetName, const TRealMatrix &sourceMatrix, const TIndexMatrix &sensorMask, const TReduceOperator reduceOp)
Constructor.
virtual void Checkpoint()
Checkpoint the stream.
virtual void FlushBufferToFile()
Flush the buffer to the file.
Output stream for quantities sampled by an index sensor mask.
TReduceOperator
How to aggregate data.
cudaEvent_t eventSamplingFinished
Has the sampling finished?
Abstract base class for output data streams (sampled data).
virtual void Close()
Close stream (apply post-processing if necessary, flush data and close).
hid_t dataset
Handle to a HDF5 dataset.
The class for real matrices.
Definition: RealMatrix.h:45
The class for 64b unsigned integers (indices). It is used for sensor_mask_index or sensor_corners_mas...
Definition: IndexMatrix.h:46
virtual void FlushRaw()
Flush data to disk (from raw streams only).
size_t sampledTimeStep
Time step to store (N/A for aggregated)
virtual void PostProcess()
Apply post-processing on the buffer and flush it to the file.
virtual ~TIndexOutputHDF5Stream()
Destructor.
Class wrapping the HDF5 routines.
Definition: HDF5_File.h:500
virtual void Create()
Create a HDF5 stream and allocate data for it.
const TIndexMatrix & sensorMask
Sensor mask to sample data.