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
Tuxyz_sgxyzMatrix Class Reference

The velocity matrix. More...

#include <UXYZ_SGXYZMatrix.h>

Inheritance diagram for Tuxyz_sgxyzMatrix:
Collaboration diagram for Tuxyz_sgxyzMatrix:

Public Member Functions

 Tuxyz_sgxyzMatrix (struct TDimensionSizes DimensionSizes)
 Constructor. More...
 
void Compute_dt_rho_sg_mul_ifft_div_2 (const TRealMatrix &dt_rho_0_sgx, TFFTWComplexMatrix &FFT)
 Compute dt ./ rho0_sgx .* ifft (FFT). More...
 
void Compute_dt_rho_sg_mul_ifft_div_2 (const float dt_rho_0_sgx, TFFTWComplexMatrix &FFT)
 Compute dt ./ rho0_sgx .* ifft (FFT), if rho0_sgx is scalar, uniform grid. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
void AddTransducerSource (const TIndexMatrix &u_source_index, TIndexMatrix &delay_mask, const TRealMatrix &transducer_signal)
 Add transducer data source to X component. More...
 
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. More...
 
virtual ~Tuxyz_sgxyzMatrix ()
 Destructor.
 
- Public Member Functions inherited from TRealMatrix
 TRealMatrix (const TDimensionSizes &DimensionSizes)
 Constructor. More...
 
virtual ~TRealMatrix ()
 Destructor.
 
virtual void ReadDataFromHDF5File (THDF5_File &HDF5_File, const char *MatrixName)
 Read data from the HDF5 file - only from the root group. More...
 
virtual void WriteDataToHDF5File (THDF5_File &HDF5_File, const char *MatrixName, const size_t CompressionLevel)
 Write data into the HDF5 file. More...
 
float & operator[] (const size_t &index)
 operator []. More...
 
const float & operator[] (const size_t &index) const
 operator [], constant version. More...
 
float & GetElementFrom3D (const size_t X, const size_t Y, const size_t Z)
 Get element from 3D matrix. More...
 
- Public Member Functions inherited from TBaseFloatMatrix
 TBaseFloatMatrix ()
 Default constructor.
 
virtual TDimensionSizes GetDimensionSizes () const
 Get dimension sizes of the matrix.
 
virtual size_t GetTotalElementCount () const
 Get element count of the matrix.
 
virtual size_t GetTotalAllocatedElementCount () const
 Get total allocated element count (might differ from total element count used for the simulation because of padding).
 
virtual ~TBaseFloatMatrix ()
 Destructor.
 
virtual void CopyData (const TBaseFloatMatrix &src)
 Copy data from other matrix with the same size. More...
 
virtual void ZeroMatrix ()
 Zero all elements of the matrix (NUMA first touch). More...
 
virtual void ScalarDividedBy (const float scalar)
 Divide scalar/ matrix_element[i]. More...
 
virtual float * GetRawData ()
 Get raw data out of the class (for direct kernel access).
 
virtual const float * GetRawData () const
 Get raw data out of the class (for direct kernel access).
 
- Public Member Functions inherited from TBaseMatrix
 TBaseMatrix ()
 Default constructor.
 
virtual ~TBaseMatrix ()
 Destructor.
 

Protected Member Functions

 Tuxyz_sgxyzMatrix (const Tuxyz_sgxyzMatrix &src)
 
Tuxyz_sgxyzMatrixoperator= (const Tuxyz_sgxyzMatrix &src)
 operator = not allowed for public.
 
- Protected Member Functions inherited from TRealMatrix
 TRealMatrix ()
 Default constructor is not allowed for public.
 
 TRealMatrix (const TRealMatrix &src)
 Copy constructor not allowed for public.
 
TRealMatrixoperator= (const TRealMatrix &src)
 Operator = is not allowed for public.
 
virtual void InitDimensions (const TDimensionSizes &DimensionSizes)
 Init dimension. More...
 
- Protected Member Functions inherited from TBaseFloatMatrix
virtual void AllocateMemory ()
 Memory allocation. More...
 
virtual void FreeMemory ()
 Memory deallocation. More...
 
 TBaseFloatMatrix (const TBaseFloatMatrix &src)
 Copy constructor is not directly allowed.
 
TBaseFloatMatrixoperator= (const TBaseFloatMatrix &src)
 operator = is not directly allowed.
 

Additional Inherited Members

- Protected Attributes inherited from TBaseFloatMatrix
size_t pTotalElementCount
 Total number of elements.
 
size_t pTotalAllocatedElementCount
 Total number of allocated elements (in terms of floats).
 
struct TDimensionSizes pDimensionSizes
 Dimension sizes.
 
size_t pDataRowSize
 Size of a 1D row in X dimension.
 
size_t p2DDataSliceSize
 Size of a 2D slab (X,Y).
 
float * pMatrixData
 Raw matrix data.
 

Detailed Description

The velocity matrix. This class implements a couple of kernels that modify the particle velocity.

Definition at line 48 of file UXYZ_SGXYZMatrix.h.

Constructor & Destructor Documentation

Tuxyz_sgxyzMatrix::Tuxyz_sgxyzMatrix ( struct TDimensionSizes  DimensionSizes)
inline

Constructor allocating memory.

Parameters
[in]DimensionSizes

Definition at line 57 of file UXYZ_SGXYZMatrix.h.

Tuxyz_sgxyzMatrix::Tuxyz_sgxyzMatrix ( const Tuxyz_sgxyzMatrix src)
protected

Copy constructor not allowed for public.

Parameters
src

Member Function Documentation

void Tuxyz_sgxyzMatrix::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.

Parameters
[in]u_source_input- Source input to add
[in]u_source_index- 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 639 of file UXYZ_SGXYZMatrix.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void Tuxyz_sgxyzMatrix::AddTransducerSource ( const TIndexMatrix u_source_index,
TIndexMatrix delay_mask,
const TRealMatrix transducer_signal 
)

Add transducer data to X dimension

Parameters
[in]u_source_index- Index matrix
[in,out]delay_mask- Index matrix - all elements incremented by one
[in]transducer_signal- Transducer signal

Definition at line 613 of file UXYZ_SGXYZMatrix.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void Tuxyz_sgxyzMatrix::Compute_dt_rho_sg_mul_ifft_div_2 ( const TRealMatrix dt_rho0_sg,
TFFTWComplexMatrix FFT 
)

Compute dt ./ rho0_sgx .* the content of the matrix.

Parameters
[in]dt_rho0_sg- matrix with the component of dt .* rho0_sg{x,y,z}
[in]FFT- FFT matrix (data will be overwritten)

Definition at line 56 of file UXYZ_SGXYZMatrix.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void Tuxyz_sgxyzMatrix::Compute_dt_rho_sg_mul_ifft_div_2 ( const float  dt_rho_0_sgx,
TFFTWComplexMatrix FFT 
)

Compute dt./rho0_sgx .* ifft (FFT) if rho0_sgx is scalar, uniform case

Parameters
[in]dt_rho_0_sgx- scalar value
[in]FFT- FFT matrix (data will be overwritten)

Definition at line 79 of file UXYZ_SGXYZMatrix.cpp.

Here is the call graph for this function:

void Tuxyz_sgxyzMatrix::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, nonuniform non uniform grid, x component.

Parameters
[in]dt_rho_0_sgx- scalar value
[in]dxudxn_sgx- non-uniform mapping
[in]FFT- FFT matrix (data will be overwritten)

Definition at line 104 of file UXYZ_SGXYZMatrix.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void Tuxyz_sgxyzMatrix::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_sgy .* ifft (FFT).

  • if rho0_sgx is scalar, nonuniform non uniform grid, y component.
    Parameters
    [in]dt_rho_0_sgy- scalar value
    [in]dyudyn_sgy- non-uniform mapping
    [in]FFT- FFT matrix (data will be overwritten)

Definition at line 138 of file UXYZ_SGXYZMatrix.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void Tuxyz_sgxyzMatrix::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_sgz .* ifft (FFT). if rho0_sgx is scalar, non uniform grid, z component.

Parameters
[in]dt_rho_0_sgz- scalar value
[in]dzudzn_sgz- non-uniform mapping
[in]FFT- FFT matrix (data will be overwritten)

Definition at line 172 of file UXYZ_SGXYZMatrix.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void Tuxyz_sgxyzMatrix::Compute_ux_sgx_normalize ( const TRealMatrix FFT_p,
const TRealMatrix dt_rho0,
const TRealMatrix pml 
)

Compute a new value for ux_sgx. Default case (heterogenous).

Parameters
[in]FFT_p- fft of pressure
[in]dt_rho0- dt_rho0_sgx
[in]pml- pml_x

Definition at line 209 of file UXYZ_SGXYZMatrix.cpp.

Here is the caller graph for this function:

void Tuxyz_sgxyzMatrix::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. This is the case for rho0 being a scalar and a non-uniform grid.

Parameters
[in]FFT_p- matrix
[in]dt_rho0- scalar
[in]dxudxn_sgx- scalar
[in]pml- matrix

Definition at line 295 of file UXYZ_SGXYZMatrix.cpp.

Here is the caller graph for this function:

void Tuxyz_sgxyzMatrix::Compute_ux_sgx_normalize_scalar_uniform ( const TRealMatrix FFT_p,
const float  dt_rho0,
const TRealMatrix pml 
)

Compute a new value of ux_sgx. This is the case for rho0 being a scalar and a uniform grid.

Parameters
[in]FFT_p- matrix
[in]dt_rho0- scalar
[in]pml- matrix

Definition at line 251 of file UXYZ_SGXYZMatrix.cpp.

Here is the caller graph for this function:

void Tuxyz_sgxyzMatrix::Compute_uy_sgy_normalize ( const TRealMatrix FFT_p,
const TRealMatrix dt_rho0,
const TRealMatrix pml 
)

Compute new value of uy_sgy. Default case (heterogenous).

Parameters
[in]FFT_p- fft of pressure
[in]dt_rho0- dt_rh0_sgy
[in]pml- pml_y

Definition at line 340 of file UXYZ_SGXYZMatrix.cpp.

Here is the caller graph for this function:

void Tuxyz_sgxyzMatrix::Compute_uy_sgy_normalize_scalar_nonuniform ( const TRealMatrix FFT_p,
const float  dt_rho0,
const TRealMatrix dyudyn_sgy,
const TRealMatrix pml 
)

Compute new value of uy_sgy, This is the case for rho0 being a scalar and a non-uniform grid.

Parameters
[in]FFT_p- matrix
[in]dt_rho0- scalar
[in]dyudyn_sgy- scalar
[in]pml- matrix

Definition at line 429 of file UXYZ_SGXYZMatrix.cpp.

Here is the caller graph for this function:

void Tuxyz_sgxyzMatrix::Compute_uy_sgy_normalize_scalar_uniform ( const TRealMatrix FFT_p,
const float  dt_rho0,
const TRealMatrix pml 
)

Compute new value of uy_sgy. This is the case for rho0 being a scalar and a uniform grid.

Parameters
[in]FFT_p- matrix
[in]dt_rho0- scalar
[in]pml- matrix

Definition at line 384 of file UXYZ_SGXYZMatrix.cpp.

Here is the caller graph for this function:

void Tuxyz_sgxyzMatrix::Compute_uz_sgz_normalize ( const TRealMatrix FFT_p,
const TRealMatrix dt_rho0,
const TRealMatrix pml 
)

Compute new value of uz_sgz. Default case (heterogenous).

Parameters
[in]FFT_p- fft of pressure
[in]dt_rho0- dt_rh0_sgz
[in]pml- pml_z

Definition at line 476 of file UXYZ_SGXYZMatrix.cpp.

Here is the caller graph for this function:

void Tuxyz_sgxyzMatrix::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. This is the case for rho0 being a scalar and a non-uniform grid.

Parameters
[in]FFT_p- matrix
[in]dt_rho0- scalar
[in]dzudzn_sgz- scalar
[in]pml- matrix

Definition at line 566 of file UXYZ_SGXYZMatrix.cpp.

Here is the caller graph for this function:

void Tuxyz_sgxyzMatrix::Compute_uz_sgz_normalize_scalar_uniform ( const TRealMatrix FFT_p,
const float  dt_rho0,
const TRealMatrix pml 
)

Compute a new value for uz_sgz. This is the case for rho0 being a scalar and a uniform grid.

Parameters
[in]FFT_p- matrix
[in]dt_rho0- scalar
[in]pml- matrix

Definition at line 520 of file UXYZ_SGXYZMatrix.cpp.

Here is the caller graph for this function:


The documentation for this class was generated from the following files: