![]() |
kspaceFirstOrder3D-OMP 1.0
The C++ implementation of the k-wave toolbox for the time-domain simulation of acoustic wave fields in 3D
|
The velocity matrix.
#include <UXYZ_SGXYZMatrix.h>
Public Member Functions | |
Tuxyz_sgxyzMatrix (struct TDimensionSizes DimensionSizes) | |
Constructor. | |
void | Compute_dt_rho_sg_mul_ifft_div_2 (TRealMatrix &dt_rho_0_sgx, TFFTWComplexMatrix &FFT) |
compute this formula dt./rho0_sgx .* ifft (FFT) | |
void | Compute_dt_rho_sg_mul_ifft_div_2 (float dt_rho_0_sgx, TFFTWComplexMatrix &FFT) |
compute this formula dt./rho0_sgx .* ifft (FFT), if rho0_sgx is scalar, uniform grid | |
void | Compute_dt_rho_sg_mul_ifft_div_2_scalar_nonuniform_x (float dt_rho_0_sgx, TRealMatrix &dxudxn_sgx, TFFTWComplexMatrix &FFT) |
compute this formula dt./rho0_sgx .* ifft (FFT), if rho0_sgx is scalar, non uniform grid, x component | |
void | Compute_dt_rho_sg_mul_ifft_div_2_scalar_nonuniform_y (float dt_rho_0_sgy, TRealMatrix &dyudyn_sgy, TFFTWComplexMatrix &FFT) |
compute this formula dt./rho0_sgx .* ifft (FFT), if rho0_sgx is scalar, non uniform grid, y component | |
void | Compute_dt_rho_sg_mul_ifft_div_2_scalar_nonuniform_z (float dt_rho_0_sgz, TRealMatrix &dzudzn_sgz, TFFTWComplexMatrix &FFT) |
compute this formula dt./rho0_sgx .* ifft (FFT), if rho0_sgx is scalar, non uniform grid, z component | |
void | Compute_ux_sgx_normalize (TRealMatrix &FFT_p, TRealMatrix &dt_rho0, TRealMatrix &pml) |
Compute new value of ux_sgx, default case. | |
void | Compute_ux_sgx_normalize_scalar_uniform (TRealMatrix &FFT_p, float dt_rho0, TRealMatrix &pml) |
Compute new value of ux_sgx, scalar, uniform case. | |
void | Compute_ux_sgx_normalize_scalar_nonuniform (TRealMatrix &FFT_p, float dt_rho0, TRealMatrix &dxudxn_sgx, TRealMatrix &pml) |
Compute new value of ux_sgx, scalar, non-uniform case. | |
void | Compute_uy_sgy_normalize (TRealMatrix &FFT_p, TRealMatrix &dt_rho0, TRealMatrix &pml) |
Compute new value of uy_sgy, default case. | |
void | Compute_uy_sgy_normalize_scalar_uniform (TRealMatrix &FFT_p, float dt_rho0, TRealMatrix &pml) |
Compute new value of uy_sgy, scalar, uniform case. | |
void | Compute_uy_sgy_normalize_scalar_nonuniform (TRealMatrix &FFT_p, float dt_rho0, TRealMatrix &dyudyn_sgy, TRealMatrix &pml) |
Compute new value of uy_sgy, scalar, non-uniform case. | |
void | Compute_uz_sgz_normalize (TRealMatrix &FFT_p, TRealMatrix &dt_rho0, TRealMatrix &pml) |
Compute new value for uz_sgz, default case. | |
void | Compute_uz_sgz_normalize_scalar_uniform (TRealMatrix &FFT_p, float &dt_rho0, TRealMatrix &pml) |
Compute new value for uz_sgz, scalar, uniform case. | |
void | Compute_uz_sgz_normalize_scalar_nonuniform (TRealMatrix &FFT_p, float &dt_rho0, TRealMatrix &dzudzn_sgz, TRealMatrix &pml) |
Compute new value for uz_sgz, scalar, non-uniform case. | |
void | AddTransducerSource (TLongMatrix &u_source_index, TLongMatrix &delay_mask, TRealMatrix &transducer_signal) |
Add transducer data source to X component. | |
void | Add_u_source (TRealMatrix &u_source_input, TLongMatrix &u_source_index, int t_index, long u_source_mode, long u_source_many) |
Add in velocity source terms. | |
virtual | ~Tuxyz_sgxyzMatrix () |
Destructor. | |
Protected Member Functions | |
Tuxyz_sgxyzMatrix (const Tuxyz_sgxyzMatrix &src) | |
Tuxyz_sgxyzMatrix & | operator= (const Tuxyz_sgxyzMatrix &src) |
operator = not allowed for public |
Definition at line 43 of file UXYZ_SGXYZMatrix.h.
Tuxyz_sgxyzMatrix::Tuxyz_sgxyzMatrix | ( | struct TDimensionSizes | DimensionSizes | ) | [inline] |
[in] | DimensionSizes |
Definition at line 50 of file UXYZ_SGXYZMatrix.h.
Tuxyz_sgxyzMatrix::Tuxyz_sgxyzMatrix | ( | const Tuxyz_sgxyzMatrix & | src | ) | [protected] |
Copy constructor not allowed for public
src |
void Tuxyz_sgxyzMatrix::Add_u_source | ( | TRealMatrix & | u_source_input, |
TLongMatrix & | u_source_index, | ||
int | t_index, | ||
long | u_source_mode, | ||
long | u_source_many | ||
) |
Add in velocity source terms
[in] | u_source_input | - Source input to add |
[in] | u_source_index | - long index matrix |
[in] | t_index | - actual time step |
[in] | u_source_mode | - Mode 0 = dirichlet boundary, 1 = add in |
[in] | u_source_many | - 0 = One series, 1 = multiple series |
Definition at line 656 of file UXYZ_SGXYZMatrix.cpp.
void Tuxyz_sgxyzMatrix::AddTransducerSource | ( | TLongMatrix & | u_source_index, |
TLongMatrix & | delay_mask, | ||
TRealMatrix & | transducer_signal | ||
) |
Add transducer data to X dimension
[in] | u_source_index | - long index matrix |
[in,out] | delay_mask | - long index matrix - modified inside (+1) |
[in] | transducer_signal | - transducer signal |
Definition at line 629 of file UXYZ_SGXYZMatrix.cpp.
void Tuxyz_sgxyzMatrix::Compute_dt_rho_sg_mul_ifft_div_2 | ( | float | dt_rho_0_sgx, |
TFFTWComplexMatrix & | FFT | ||
) |
Compute dt./rho0_sgx .* ifft (FFT), if rho0_sgx is scalar, uniform case
[in] | dt_rho_0_sgx | - scalar value |
[in] | FFT | - FFT matrix |
Definition at line 81 of file UXYZ_SGXYZMatrix.cpp.
void Tuxyz_sgxyzMatrix::Compute_dt_rho_sg_mul_ifft_div_2 | ( | TRealMatrix & | dt_rho0_sg, |
TFFTWComplexMatrix & | FFT | ||
) |
Compute dt./rho0_sgx .* the content of the matrix
[in] | dt_rho0_sg | - matrix with the component of dt .* rho0_sg{x,y,z} |
[in] | FFT | - FFT matrix |
Definition at line 55 of file UXYZ_SGXYZMatrix.cpp.
void Tuxyz_sgxyzMatrix::Compute_dt_rho_sg_mul_ifft_div_2_scalar_nonuniform_x | ( | float | dt_rho_0_sgx, |
TRealMatrix & | dxudxn_sgx, | ||
TFFTWComplexMatrix & | FFT | ||
) |
Compute dt./rho0_sgx .* ifft (FFT), when rho0_sgx is scalar, nonuniform
[in] | dt_rho_0_sgx | - scalar value |
[in] | dxudxn_sgx | - non-uniform mapping |
[in] | FFT | - FFT matrix |
Definition at line 107 of file UXYZ_SGXYZMatrix.cpp.
void Tuxyz_sgxyzMatrix::Compute_dt_rho_sg_mul_ifft_div_2_scalar_nonuniform_y | ( | float | dt_rho_0_sgy, |
TRealMatrix & | dyudyn_sgy, | ||
TFFTWComplexMatrix & | FFT | ||
) |
Compute dt./rho0_sgy .* ifft (FFT), if rho0_sgx is scalar, nonuniform
[in] | dt_rho_0_sgy | - scalar value |
[in] | dyudyn_sgy | - non-uniform mapping |
[in] | FFT | - FFT matrix |
Definition at line 145 of file UXYZ_SGXYZMatrix.cpp.
void Tuxyz_sgxyzMatrix::Compute_dt_rho_sg_mul_ifft_div_2_scalar_nonuniform_z | ( | float | dt_rho_0_sgz, |
TRealMatrix & | dzudzn_sgz, | ||
TFFTWComplexMatrix & | FFT | ||
) |
Compute dt./rho0_sgz .* ifft (FFT), if rho0_sgx is scalar, uniform
[in] | dt_rho_0_sgz | - scalar value |
[in] | dzudzn_sgz | |
[in] | FFT | - FFT matrix |
Definition at line 183 of file UXYZ_SGXYZMatrix.cpp.
void Tuxyz_sgxyzMatrix::Compute_ux_sgx_normalize | ( | TRealMatrix & | FFT_p, |
TRealMatrix & | dt_rho0, | ||
TRealMatrix & | pml | ||
) |
Compute new value for ux_sgx.
[in] | FFT_p | - fft of pressure |
[in] | dt_rho0 | - dt_rho0_sgx |
[in] | pml | - pml_x |
Definition at line 224 of file UXYZ_SGXYZMatrix.cpp.
void Tuxyz_sgxyzMatrix::Compute_ux_sgx_normalize_scalar_nonuniform | ( | TRealMatrix & | FFT_p, |
float | dt_rho0, | ||
TRealMatrix & | dxudxn_sgx, | ||
TRealMatrix & | pml | ||
) |
Compute_ux_sgx_normalize if rho0 is a scalar, non uniform.
[in] | FFT_p | - matrix |
[in] | dt_rho0 | - scalar |
[in] | dxudxn_sgx | - scalar |
[in] | pml | - matrix |
Definition at line 309 of file UXYZ_SGXYZMatrix.cpp.
void Tuxyz_sgxyzMatrix::Compute_ux_sgx_normalize_scalar_uniform | ( | TRealMatrix & | FFT_p, |
float | dt_rho0, | ||
TRealMatrix & | pml | ||
) |
Compute_ux_sgx_normalize if rho0 is a scalar.
[in] | FFT_p | - matrix |
[in] | dt_rho0 | - scalar |
[in] | pml | - matrix |
Definition at line 265 of file UXYZ_SGXYZMatrix.cpp.
void Tuxyz_sgxyzMatrix::Compute_uy_sgy_normalize | ( | TRealMatrix & | FFT_p, |
TRealMatrix & | dt_rho0, | ||
TRealMatrix & | pml | ||
) |
Compute new value for uy_sgy
[in] | FFT_p | - fft of pressure |
[in] | dt_rho0 | - dt_rh0_sgy |
[in] | pml | - pml_y |
Definition at line 353 of file UXYZ_SGXYZMatrix.cpp.
void Tuxyz_sgxyzMatrix::Compute_uy_sgy_normalize_scalar_nonuniform | ( | TRealMatrix & | FFT_p, |
float | dt_rho0, | ||
TRealMatrix & | dyudyn_sgy, | ||
TRealMatrix & | pml | ||
) |
Compute_uy_sgy_normalize if rho0 is a scalar, non uniform
[in] | FFT_p | - matrix |
[in] | dt_rho0 | - scalar |
[in] | dyudyn_sgy | - scalar |
[in] | pml | - matrix |
Definition at line 444 of file UXYZ_SGXYZMatrix.cpp.
void Tuxyz_sgxyzMatrix::Compute_uy_sgy_normalize_scalar_uniform | ( | TRealMatrix & | FFT_p, |
float | dt_rho0, | ||
TRealMatrix & | pml | ||
) |
Compute_uy_sgy_normalize if rho0 is a scalar
[in] | FFT_p | - matrix |
[in] | dt_rho0 | - scalar |
[in] | pml | - matrix |
Definition at line 398 of file UXYZ_SGXYZMatrix.cpp.
void Tuxyz_sgxyzMatrix::Compute_uz_sgz_normalize | ( | TRealMatrix & | FFT_p, |
TRealMatrix & | dt_rho0, | ||
TRealMatrix & | pml | ||
) |
Compute new value of uz_sgz
[in] | FFT_p | - fft of pressure |
[in] | dt_rho0 | - dt_rh0_sgz |
[in] | pml | - pml_z |
Definition at line 490 of file UXYZ_SGXYZMatrix.cpp.
void Tuxyz_sgxyzMatrix::Compute_uz_sgz_normalize_scalar_nonuniform | ( | TRealMatrix & | FFT_p, |
float & | dt_rho0, | ||
TRealMatrix & | dzudzn_sgz, | ||
TRealMatrix & | pml | ||
) |
Compute_uz_sgz_normalize if rho0 is a scalar, non uniform
[in] | FFT_p | - matrix |
[in] | dt_rho0 | - scalar |
[in] | dzudzn_sgz | - scalar |
[in] | pml | - matrix |
Definition at line 582 of file UXYZ_SGXYZMatrix.cpp.
void Tuxyz_sgxyzMatrix::Compute_uz_sgz_normalize_scalar_uniform | ( | TRealMatrix & | FFT_p, |
float & | dt_rho0, | ||
TRealMatrix & | pml | ||
) |
Compute_uz_sgz_normalize if rho0 is a scalar
[in] | FFT_p | - matrix |
[in] | dt_rho0 | - scalar |
[in] | pml | - matrix |
Definition at line 534 of file UXYZ_SGXYZMatrix.cpp.