kspaceFirstOrder3D-OMP  1.1
The C++ implementation of the k-wave toolbox for the time-domain simulation of acoustic wave fields in 3D
 All Classes Files Functions Variables Typedefs Enumerations Friends Pages
UXYZ_SGXYZMatrix.h
Go to the documentation of this file.
1 /**
2  * @file UXYZ_SGXYZMatrix.h
3  * @author Jiri Jaros
4  * Faculty of Information Technology\n
5  * Brno University of Technology \n
6  * jarosjir@fit.vutbr.cz
7  *
8  * @brief The header file containing the particle velocity matrix.
9  *
10  * @version kspaceFirstOrder3D 2.15
11  *
12  * @date 28 July 2011, 11:37 (created) \n
13  * 26 September 2014, 14:15 (revised)
14  *
15  * @section License
16  * This file is part of the C++ extension of the k-Wave Toolbox (http://www.k-wave.org).\n
17  * Copyright (C) 2014 Jiri Jaros and Bradley Treeby
18  *
19  * This file is part of k-Wave. k-Wave is free software: you can redistribute it
20  * and/or modify it under the terms of the GNU Lesser General Public License as
21  * published by the Free Software Foundation, either version 3 of the License,
22  * or (at your option) any later version.
23  *
24  * k-Wave is distributed in the hope that it will be useful, but
25  * WITHOUT ANY WARRANTY; without even the implied warranty of
26  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
27  * See the GNU Lesser General Public License for more details.
28  *
29  * You should have received a copy of the GNU Lesser General Public License
30  * along with k-Wave. If not, see <http://www.gnu.org/licenses/>.
31  */
32 
33 
34 #ifndef UXYZ_SGXYZREALMATRIX_H
35 #define UXYZ_SGXYZREALMATRIX_H
36 
37 
41 
42 /**
43  * @class Tuxyz_sgxyzMatrix.
44  * @brief The velocity matrix
45  * @details The velocity matrix. This class implements a couple of kernels that
46  * modify the particle velocity.
47  */
49 {
50  public:
51 
52  /**
53  * @brief Constructor.
54  * @details Constructor allocating memory.
55  * @param [in] DimensionSizes
56  */
57  Tuxyz_sgxyzMatrix(struct TDimensionSizes DimensionSizes) :
58  TRealMatrix(DimensionSizes)
59  {};
60 
61 
62  /// Compute dt ./ rho0_sgx .* ifft (FFT).
63  void Compute_dt_rho_sg_mul_ifft_div_2(const TRealMatrix & dt_rho_0_sgx,
64  TFFTWComplexMatrix& FFT);
65  /// Compute dt ./ rho0_sgx .* ifft (FFT), if rho0_sgx is scalar, uniform grid.
66  void Compute_dt_rho_sg_mul_ifft_div_2(const float dt_rho_0_sgx,
67  TFFTWComplexMatrix& FFT);
68 
69  /// Compute dt ./ rho0_sgx .* ifft (FFT), if rho0_sgx is scalar, non uniform grid, x component.
70  void Compute_dt_rho_sg_mul_ifft_div_2_scalar_nonuniform_x(const float dt_rho_0_sgx,
71  const TRealMatrix & dxudxn_sgx,
72  TFFTWComplexMatrix& FFT);
73  /// Compute dt./rho0_sgx .* ifft (FFT), if rho0_sgx is scalar, non uniform grid, y component.
74  void Compute_dt_rho_sg_mul_ifft_div_2_scalar_nonuniform_y(const float dt_rho_0_sgy,
75  const TRealMatrix & dyudyn_sgy,
76  TFFTWComplexMatrix& FFT);
77  /// Compute dt./rho0_sgx .* ifft (FFT), if rho0_sgx is scalar, non uniform grid, z component.
78  void Compute_dt_rho_sg_mul_ifft_div_2_scalar_nonuniform_z(const float dt_rho_0_sgz,
79  const TRealMatrix & dzudzn_sgz,
80  TFFTWComplexMatrix& FFT);
81 
82  //------------------------- X dimension -----------------------------------//
83  /// Compute a new value of ux_sgx, default case
84  void Compute_ux_sgx_normalize(const TRealMatrix& FFT_p,
85  const TRealMatrix& dt_rho0,
86  const TRealMatrix& pml);
87  /// Compute a new value of ux_sgx, scalar, uniform case
89  const float dt_rho0,
90  const TRealMatrix& pml);
91  /// Compute a new value of ux_sgx, scalar, non-uniform case
93  const float dt_rho0,
94  const TRealMatrix& dxudxn_sgx,
95  const TRealMatrix& pml);
96 
97  //------------------------- Y dimension ------------------------------------//
98  /// Compute a new value of uy_sgy, default case
99  void Compute_uy_sgy_normalize(const TRealMatrix& FFT_p,
100  const TRealMatrix& dt_rho0,
101  const TRealMatrix& pml);
102  /// Compute a new value of uy_sgy, scalar, uniform case
104  const float dt_rho0,
105  const TRealMatrix& pml);
106  /// Compute a new value of uy_sgy, scalar, non-uniform case
108  const float dt_rho0,
109  const TRealMatrix& dyudyn_sgy,
110  const TRealMatrix& pml);
111 
112  //------------------------- Z dimension -----------------------------------//
113  /// Compute a new value for uz_sgz, default case.
114  void Compute_uz_sgz_normalize(const TRealMatrix& FFT_p,
115  const TRealMatrix& dt_rho0,
116  const TRealMatrix& pml);
117  /// Compute a new value for uz_sgz, scalar, uniform case.
119  const float dt_rho0,
120  const TRealMatrix& pml);
121  /// Compute a new value for uz_sgz, scalar, non-uniform case.
123  const float dt_rho0,
124  const TRealMatrix& dzudzn_sgz,
125  const TRealMatrix& pml);
126 
127 
128  /// Add transducer data source to X component.
129  void AddTransducerSource(const TIndexMatrix& u_source_index,
130  TIndexMatrix& delay_mask,
131  const TRealMatrix & transducer_signal);
132 
133  /// Add in velocity source terms.
134  void Add_u_source(const TRealMatrix & u_source_input,
135  const TIndexMatrix& u_source_index,
136  const size_t t_index,
137  const size_t u_source_mode,
138  const size_t u_source_many);
139 
140  /// Destructor
141  virtual ~Tuxyz_sgxyzMatrix() {};
142 
143 protected:
144  // Default constructor not allowed for public.
146 
147  /**
148  * Copy constructor not allowed for public.
149  * @param src
150  */
152 
153  /// operator = not allowed for public.
155 private:
156 
157 }; // end of Tuxyz_sgxyzMatrix
158 
159 #endif /* UX_SGREALMATRIX_H */
160 
void Compute_dt_rho_sg_mul_ifft_div_2(const TRealMatrix &dt_rho_0_sgx, TFFTWComplexMatrix &FFT)
Compute dt ./ rho0_sgx .* ifft (FFT).
void Compute_uz_sgz_normalize(const TRealMatrix &FFT_p, const TRealMatrix &dt_rho0, const TRealMatrix &pml)
Compute a new value for uz_sgz, default case.
Tuxyz_sgxyzMatrix & operator=(const Tuxyz_sgxyzMatrix &src)
operator = not allowed for public.
void Compute_ux_sgx_normalize_scalar_nonuniform(const TRealMatrix &FFT_p, const float dt_rho0, const TRealMatrix &dxudxn_sgx, const TRealMatrix &pml)
Compute a new value of ux_sgx, scalar, non-uniform case.
void Compute_uy_sgy_normalize(const TRealMatrix &FFT_p, const TRealMatrix &dt_rho0, const TRealMatrix &pml)
Compute a new value of uy_sgy, default case.
The header file containing the class for real matrices.
void Add_u_source(const TRealMatrix &u_source_input, const TIndexMatrix &u_source_index, const size_t t_index, const size_t u_source_mode, const size_t u_source_many)
Add in velocity source terms.
void Compute_uz_sgz_normalize_scalar_uniform(const TRealMatrix &FFT_p, const float dt_rho0, const TRealMatrix &pml)
Compute a new value for uz_sgz, scalar, uniform case.
virtual ~Tuxyz_sgxyzMatrix()
Destructor.
void AddTransducerSource(const TIndexMatrix &u_source_index, TIndexMatrix &delay_mask, const TRealMatrix &transducer_signal)
Add transducer data source to X component.
Class implementing 3D Real-To-Complex and Complex-To-Real transforms using FFTW interface.
void Compute_dt_rho_sg_mul_ifft_div_2_scalar_nonuniform_x(const float dt_rho_0_sgx, const TRealMatrix &dxudxn_sgx, TFFTWComplexMatrix &FFT)
Compute dt ./ rho0_sgx .* ifft (FFT), if rho0_sgx is scalar, non uniform grid, x component.
Tuxyz_sgxyzMatrix(struct TDimensionSizes DimensionSizes)
Constructor.
The header file containing the class for 64b integer matrices.
void Compute_uy_sgy_normalize_scalar_nonuniform(const TRealMatrix &FFT_p, const float dt_rho0, const TRealMatrix &dyudyn_sgy, const TRealMatrix &pml)
Compute a new value of uy_sgy, scalar, non-uniform case.
void Compute_uz_sgz_normalize_scalar_nonuniform(const TRealMatrix &FFT_p, const float dt_rho0, const TRealMatrix &dzudzn_sgz, const TRealMatrix &pml)
Compute a new value for uz_sgz, scalar, non-uniform case.
The class for real matrices.
Definition: RealMatrix.h:48
The class for 64b unsigned integers (indices). It is used for sensor_mask_index or sensor_corners_mas...
Definition: IndexMatrix.h:50
void Compute_uy_sgy_normalize_scalar_uniform(const TRealMatrix &FFT_p, const float dt_rho0, const TRealMatrix &pml)
Compute a new value of uy_sgy, scalar, uniform case.
The velocity matrix.
void Compute_dt_rho_sg_mul_ifft_div_2_scalar_nonuniform_z(const float dt_rho_0_sgz, const TRealMatrix &dzudzn_sgz, TFFTWComplexMatrix &FFT)
Compute dt./rho0_sgx .* ifft (FFT), if rho0_sgx is scalar, non uniform grid, z component.
Structure with 4D dimension sizes (3 in space and 1 in time).
void Compute_ux_sgx_normalize_scalar_uniform(const TRealMatrix &FFT_p, const float dt_rho0, const TRealMatrix &pml)
Compute a new value of ux_sgx, scalar, uniform case.
void Compute_dt_rho_sg_mul_ifft_div_2_scalar_nonuniform_y(const float dt_rho_0_sgy, const TRealMatrix &dyudyn_sgy, TFFTWComplexMatrix &FFT)
Compute dt./rho0_sgx .* ifft (FFT), if rho0_sgx is scalar, non uniform grid, y component.
void Compute_ux_sgx_normalize(const TRealMatrix &FFT_p, const TRealMatrix &dt_rho0, const TRealMatrix &pml)
Compute a new value of ux_sgx, default case.
The header file containing the class that implements 3D FFT using the FFTW interface.