kspaceFirstOrder3D-OMP  1.2
The C++ implementation of the k-wave toolbox for the time-domain simulation of acoustic wave fields in 3D
KSpaceFirstOrder3DSolver Class Reference

Class responsible for running the k-space first order 3D method. More...

#include <KSpaceFirstOrder3DSolver.h>

Collaboration diagram for KSpaceFirstOrder3DSolver:

Public Member Functions

 KSpaceFirstOrder3DSolver ()
 Constructor. More...
 
 KSpaceFirstOrder3DSolver (const KSpaceFirstOrder3DSolver &)=delete
 Copy constructor not allowed for public.
 
virtual ~KSpaceFirstOrder3DSolver ()
 Destructor. More...
 
KSpaceFirstOrder3DSolveroperator= (const KSpaceFirstOrder3DSolver &)=delete
 operator= not allowed for public.
 
virtual void allocateMemory ()
 Memory allocation. More...
 
virtual void freeMemory ()
 Memory deallocation. More...
 
virtual void loadInputData ()
 Load simulation data. More...
 
virtual void compute ()
 This method computes k-space First Order 3D simulation. More...
 
virtual size_t getMemoryUsage () const
 Get memory usage in MB on the host side. More...
 
std::string getCodeName () const
 Get code name - release code version. More...
 
void printFullCodeNameAndLicense () const
 Print the code name and license. More...
 
void setProcessorAffinity ()
 Set processor affinity. More...
 
double getTotalTime () const
 Get total simulation time. More...
 
double getPreProcessingTime () const
 Get pre-processing time. More...
 
double getDataLoadTime () const
 Get data load time. More...
 
double getSimulationTime () const
 Get simulation time (time loop). More...
 
double getPostProcessingTime () const
 Get post-processing time. More...
 
double getCumulatedTotalTime () const
 Get total simulation time accumulated over all legs. More...
 
double getCumulatedPreProcessingTime () const
 Get pre-processing time accumulated over all legs. More...
 
double getCumulatedDataLoadTime () const
 Get simulation time (time loop) accumulated over all legs. More...
 
double getCumulatedSimulationTime () const
 Get simulation time (time loop) accumulated over all legs. More...
 
double getCumulatedPostProcessingTime () const
 Get post-processing time accumulated over all legs. More...
 

Protected Member Functions

void InitializeFftwPlans ()
 Initialize FFTW plans. More...
 
void preProcessing ()
 Compute pre-processing phase. More...
 
void computeMainLoop ()
 Compute the main time loop of the kspaceFirstOrder3D. More...
 
void postProcessing ()
 Post processing, and closing the output streams. More...
 
void storeSensorData ()
 Store sensor data. More...
 
void writeOutputDataInfo ()
 Write statistics and header into the output file. More...
 
void saveCheckpointData ()
 Save checkpoint data and flush aggregated outputs into the output file. More...
 
void computeVelocity ()
 Compute new values of acoustic velocity. More...
 
void computeVelocityGradient ()
 Compute new values of acoustic velocity gradients. More...
 
void computeDensityNonliner ()
 Compute new values of acoustic density for nonlinear case. More...
 
void computeDensityLinear ()
 Compute new values of acoustic density for linear case. More...
 
void computePressureNonlinear ()
 Compute acoustic pressure for nonlinear case. More...
 
void computePressureLinear ()
 Compute acoustic pressure for linear case. More...
 
void addVelocitySource ()
 Add in all velocity sources. More...
 
void computeVelocitySourceTerm (RealMatrix &velocityMatrix, const RealMatrix &velocitySourceInput, const IndexMatrix &velocitySourceIndex)
 Add in velocity source terms. More...
 
void addPressureSource ()
 Add in pressure source. More...
 
void addInitialPressureSource ()
 Calculate initial pressure source. More...
 
void addTransducerSource ()
 Add transducer data source to velocity x component. More...
 
void generateKappa ()
 Generate kappa matrix for lossless media. More...
 
void generateKappaAndNablas ()
 Generate kappa matrix, absorbNabla1, absorbNabla2 for absorbing medium. More...
 
void generateTauAndEta ()
 Generate absorbTau, absorbEta for heterogenous medium. More...
 
void generateInitialDenisty ()
 Calculate dt ./ rho0 for nonuniform grids. More...
 
void computeC2 ()
 Calculate square of velocity. More...
 
void computeInitialVelocityHeterogeneous ()
 Compute velocity for the initial pressure problem, heterogeneous medium, uniform grid. More...
 
void computeInitialVelocityHomogeneousUniform ()
 Compute velocity for the initial pressure problem, homogeneous medium, uniform grid. More...
 
void computeInitialVelocityHomogeneousNonuniform ()
 Compute acoustic velocity for initial pressure problem, homogenous medium, nonuniform grid. More...
 
void computeVelocityHeterogeneous ()
 Compute acoustic velocity for heterogeneous medium and a uniform grid. More...
 
void computeVelocityHomogeneousUniform ()
 Compute acoustic velocity for homogeneous medium and a uniform grid. More...
 
void computeVelocityHomogeneousNonuniform ()
 Compute acoustic velocity for homogenous medium and nonuniform grid. More...
 
void computePressureGradient ()
 Compute part of the new velocity term - gradient of pressure. More...
 
template<bool bOnAScalarFlag, bool rho0ScalarFlag>
void computePressureTermsNonlinear (RealMatrix &densitySum, RealMatrix &nonlinearTerm, RealMatrix &velocityGradientSum)
 Calculate three temporary sums in the new pressure formula before taking the FFT, nonlinear absorbing case. More...
 
void computePressureTermsLinear (RealMatrix &densitySum, RealMatrix &velocityGradientSum)
 Calculate two temporary sums in the new pressure formula before taking the FFT, linear absorbing case. More...
 
void computeAbsorbtionTerm (FftwComplexMatrix &fftPart1, FftwComplexMatrix &fftPart2)
 Compute absorbing term with abosrbNabla1 and absorbNabla2. More...
 
template<bool c0ScalarFlag, bool areTauAndEtaScalars>
void sumPressureTermsNonlinear (const RealMatrix &absorbTauTerm, const RealMatrix &absorbEtaTerm, const RealMatrix &nonlinearTerm)
 Sum sub-terms to calculate new pressure, after FFTs, nonlinear case. More...
 
template<bool c0ScalarFlag, bool areTauAndEtaScalars>
void sumPressureTermsLinear (const RealMatrix &absorbTauTerm, const RealMatrix &absorbEtaTerm, const RealMatrix &densitySum)
 Sum sub-terms to calculate new pressure, after FFTs, linear case. More...
 
template<bool c0ScalarFlag, bool nonlinearFlag, bool rho0ScalarFlag>
void sumPressureTermsNonlinearLossless ()
 Sum sub-terms for new pressure, linear lossless case. More...
 
void sumPressureTermsLinearLossless ()
 Sum sub-terms for new pressure, linear lossless case. More...
 
void computeShiftedVelocity ()
 compute shifted velocity for –u_non_staggered flag. More...
 
void printStatistics ()
 Print progress statistics. More...
 
bool isTimeToCheckpoint ()
 Is time to checkpoint (save actual state on disk). More...
 
bool isCheckpointInterruption () const
 Was the loop interrupted to checkpoint? More...
 
void checkOutputFile ()
 Check the output file has the correct format and version. More...
 
void checkCheckpointFile ()
 Check the file type and the version of the checkpoint file. More...
 
void loadElapsedTimeFromOutputFile ()
 Reads the header of the output file and sets the cumulative elapsed time from the first log. More...
 
size_t get1DIndex (const size_t z, const size_t y, const size_t x, const DimensionSizes &dimensionSizes)
 Compute 1D index using 3 spatial coordinates and the size of the matrix. More...
 
RealMatrixgetKappa ()
 Get the kappa matrix from the container. More...
 
RealMatrixgetC2 ()
 Get the c^2 matrix from the container. More...
 
RealMatrixgetP ()
 Get pressure matrix. More...
 
RealMatrixgetUxSgx ()
 Get velocity matrix on staggered grid in x direction. More...
 
RealMatrixgetUySgy ()
 Get velocity matrix on staggered grid in y direction. More...
 
RealMatrixgetUzSgz ()
 Get velocity matrix on staggered grid in z direction. More...
 
RealMatrixgetUxShifted ()
 Get velocity shifted on normal grid in x direction. More...
 
RealMatrixgetUyShifted ()
 Get velocity shifted on normal grid in y direction. More...
 
RealMatrixgetUzShifted ()
 Get velocity shifted on normal grid in z direction. More...
 
RealMatrixgetDuxdx ()
 Get velocity gradient on in x direction. More...
 
RealMatrixgetDuydy ()
 Get velocity gradient on in y direction. More...
 
RealMatrixgetDuzdz ()
 Get velocity gradient on in z direction. More...
 
RealMatrixgetDtRho0Sgx ()
 Get dt * rho0Sgx matrix (time step size * ambient velocity on staggered grid in x direction). More...
 
RealMatrixgetDtRho0Sgy ()
 Get dt * rho0Sgy matrix (time step size * ambient velocity on staggered grid in y direction). More...
 
RealMatrixgetDtRho0Sgz ()
 Get dt * rho0Sgz matrix (time step size * ambient velocity on staggered grid in z direction). More...
 
RealMatrixgetRhoX ()
 Get density matrix in x direction. More...
 
RealMatrixgetRhoY ()
 Get density matrix in y direction. More...
 
RealMatrixgetRhoZ ()
 Get density matrix in z direction. More...
 
RealMatrixgetRho0 ()
 Get ambient density matrix. More...
 
ComplexMatrixgetDdxKShiftPos ()
 Get positive Fourier shift in x. More...
 
ComplexMatrixgetDdyKShiftPos ()
 Get positive Fourier shift in y. More...
 
ComplexMatrixgetDdzKShiftPos ()
 Get positive Fourier shift in z. More...
 
ComplexMatrixgetDdxKShiftNeg ()
 Get negative Fourier shift in x. More...
 
ComplexMatrixgetDdyKShiftNeg ()
 Get negative Fourier shift in y. More...
 
ComplexMatrixgetDdzKShiftNeg ()
 Get negative Fourier shift in z. More...
 
ComplexMatrixgetXShiftNegR ()
 Get negative shift for non-staggered velocity in x. More...
 
ComplexMatrixgetYShiftNegR ()
 Get negative shift for non-staggered velocity in y. More...
 
ComplexMatrixgetZShiftNegR ()
 Get negative shift for non-staggered velocity in z. More...
 
RealMatrixgetPmlXSgx ()
 Get PML on staggered grid x. More...
 
RealMatrixgetPmlYSgy ()
 Get PML on staggered grid y. More...
 
RealMatrixgetPmlZSgz ()
 Get PML on staggered grid z. More...
 
RealMatrixgetPmlX ()
 Get PML in x. More...
 
RealMatrixgetPmlY ()
 Get PML in y. More...
 
RealMatrixgetPmlZ ()
 Get PML in z. More...
 
RealMatrixgetDxudxn ()
 Non uniform grid acoustic velocity in x. More...
 
RealMatrixgetDyudyn ()
 Non uniform grid acoustic velocity in y. More...
 
RealMatrixgetDzudzn ()
 Non uniform grid acoustic velocity in z. More...
 
RealMatrixgetDxudxnSgx ()
 Non uniform grid acoustic velocity on staggered grid x. More...
 
RealMatrixgetDyudynSgy ()
 Non uniform grid acoustic velocity on staggered grid x. More...
 
RealMatrixgetDzudznSgz ()
 Non uniform grid acoustic velocity on staggered grid x. More...
 
RealMatrixgetBOnA ()
 Get B on A (nonlinear coefficient). More...
 
RealMatrixgetAbsorbTau ()
 Get absorbing coefficient Tau. More...
 
RealMatrixgetAbsorbEta ()
 Get absorbing coefficient Eta. More...
 
RealMatrixgetAbsorbNabla1 ()
 Get absorbing coefficient Nabla1. More...
 
RealMatrixgetAbsorbNabla2 ()
 Get absorbing coefficient Nabla2. More...
 
IndexMatrixgetSensorMaskIndex ()
 Get linear sensor mask (spatial geometry of the sensor). More...
 
IndexMatrixgetSensorMaskCorners ()
 Get cuboid corners sensor mask. (Spatial geometry of multiple sensors). More...
 
IndexMatrixgetVelocitySourceIndex ()
 Get velocity source geometry data. More...
 
IndexMatrixgetPressureSourceIndex ()
 Get pressure source geometry data. More...
 
IndexMatrixgetDelayMask ()
 Get delay mask for many types sources. More...
 
RealMatrixgetTransducerSourceInput ()
 Get transducer source input data (signal). More...
 
RealMatrixgetPressureSourceInput ()
 Get pressure source input data (signal). More...
 
RealMatrixgetInitialPressureSourceInput ()
 Get initial pressure source input data (whole matrix). More...
 
RealMatrixGetVelocityXSourceInput ()
 Get Velocity source input data in x direction. More...
 
RealMatrixGetVelocityYSourceInput ()
 Get Velocity source input data in y direction. More...
 
RealMatrixgetVelocityZSourceInput ()
 Get Velocity source input data in z direction. More...
 
RealMatrixgetTemp1Real3D ()
 Get first real 3D temporary matrix. More...
 
RealMatrixgetTemp2Real3D ()
 Get second real 3D temporary matrix. More...
 
RealMatrixgetTemp3Real3D ()
 Get third real 3D temporary matrix. More...
 
FftwComplexMatrixgetTempFftwX ()
 Get temporary matrix for 1D fft in x. More...
 
FftwComplexMatrixgetTempFftwY ()
 Get temporary matrix for 1D fft in y. More...
 
FftwComplexMatrixgetTempFftwZ ()
 Get temporary matrix for 1D fft in z. More...
 
FftwComplexMatrixgetTempFftwShift ()
 Get temporary matrix for fft shift. More...
 

Private Attributes

MatrixContainer mMatrixContainer
 Matrix container with all the matrix classes.
 
OutputStreamContainer mOutputStreamContainer
 Output stream container.
 
ParametersmParameters
 Global parameters of the simulation.
 
size_t mActPercent
 Percentage of the simulation done.
 
TimeMeasure mTotalTime
 Total time of the simulation.
 
TimeMeasure mPreProcessingTime
 Pre-processing time of the simulation.
 
TimeMeasure mDataLoadTime
 Data load time of the simulation.
 
TimeMeasure mSimulationTime
 Simulation time of the simulation.
 
TimeMeasure mPostProcessingTime
 Post-processing time of the simulation.
 
TimeMeasure mIterationTime
 Iteration time of the simulation.
 

Detailed Description

Class responsible for running the k-space first order 3D method. This class maintain the whole k-wave (implements the time loop).

Definition at line 57 of file KSpaceFirstOrder3DSolver.h.

Constructor & Destructor Documentation

◆ KSpaceFirstOrder3DSolver()

KSpaceFirstOrder3DSolver::KSpaceFirstOrder3DSolver ( )

Constructor of the class.

Definition at line 74 of file KSpaceFirstOrder3DSolver.cpp.

◆ ~KSpaceFirstOrder3DSolver()

KSpaceFirstOrder3DSolver::~KSpaceFirstOrder3DSolver ( )
virtual

Destructor of the class.

Definition at line 92 of file KSpaceFirstOrder3DSolver.cpp.

Member Function Documentation

◆ addInitialPressureSource()

void KSpaceFirstOrder3DSolver::addInitialPressureSource ( )
protected

Calculate p0 source when necessary.

Matlab code:

   % add the initial pressure to rho as a mass source
   p = source.p0;
   rhox = source.p0 ./ (3 .* c.^2);
   rhoy = source.p0 ./ (3 .* c.^2);
   rhoz = source.p0 ./ (3 .* c.^2);

   % compute u(t = t1 + dt/2) based on the assumption u(dt/2) = -u(-dt/2)
   % which forces u(t = t1) = 0
   ux_sgx = dt .* rho0_sgx_inv .* real(ifftn( bsxfun(@times, ddx_k_shift_pos, kappa .* fftn(p)) )) / 2;
   uy_sgy = dt .* rho0_sgy_inv .* real(ifftn( bsxfun(@times, ddy_k_shift_pos, kappa .* fftn(p)) )) / 2;
   uz_sgz = dt .* rho0_sgz_inv .* real(ifftn( bsxfun(@times, ddz_k_shift_pos, kappa .* fftn(p)) )) / 2;

Definition at line 1465 of file KSpaceFirstOrder3DSolver.cpp.

◆ addPressureSource()

void KSpaceFirstOrder3DSolver::addPressureSource ( )
protected

Add in pressure source.

Definition at line 1399 of file KSpaceFirstOrder3DSolver.cpp.

◆ addTransducerSource()

void KSpaceFirstOrder3DSolver::addTransducerSource ( )
protected

Add transducer data source to velocity x component.

Definition at line 1518 of file KSpaceFirstOrder3DSolver.cpp.

◆ addVelocitySource()

void KSpaceFirstOrder3DSolver::addVelocitySource ( )
protected

Add u source to the particle velocity.

Definition at line 1338 of file KSpaceFirstOrder3DSolver.cpp.

◆ allocateMemory()

void KSpaceFirstOrder3DSolver::allocateMemory ( )
virtual

The method allocates the matrix container, creates all matrices and creates all output streams (however not allocating memory).

Definition at line 102 of file KSpaceFirstOrder3DSolver.cpp.

◆ checkCheckpointFile()

void KSpaceFirstOrder3DSolver::checkCheckpointFile ( )
protected
Exceptions
ios::failure- If an error happens

Check the file type and the version of the checkpoint file.

Definition at line 2630 of file KSpaceFirstOrder3DSolver.cpp.

◆ checkOutputFile()

void KSpaceFirstOrder3DSolver::checkOutputFile ( )
protected
Exceptions
ios::failure- If an error happens.

Check the output file has the correct format and version.

Definition at line 2574 of file KSpaceFirstOrder3DSolver.cpp.

◆ compute()

void KSpaceFirstOrder3DSolver::compute ( )
virtual

This method computes k-space First Order 3D simulation. It launches calculation on a given dataset going through FFT initialization, pre-processing, main loop and post-processing phases.

This method computes k-space First Order 3D simulation.

Definition at line 226 of file KSpaceFirstOrder3DSolver.cpp.

◆ computeAbsorbtionTerm()

void KSpaceFirstOrder3DSolver::computeAbsorbtionTerm ( FftwComplexMatrix fftPart1,
FftwComplexMatrix fftPart2 
)
protected
Parameters
[in,out]fftPart1- fftPart1 = absorbNabla1 .* fftPart1
[in,out]fftPart2- fftPart1 = absorbNabla1 .* fftPart2

Compute absorbing term with abosrbNabla1 and absorbNabla2.

Definition at line 2245 of file KSpaceFirstOrder3DSolver.cpp.

◆ computeC2()

void KSpaceFirstOrder3DSolver::computeC2 ( )
protected

Compute c^2.

Definition at line 1763 of file KSpaceFirstOrder3DSolver.cpp.

◆ computeDensityLinear()

void KSpaceFirstOrder3DSolver::computeDensityLinear ( )
protected

Calculate new values of acoustic density for linear case (rhoX, rhoy and rhoZ).

Matlab code:

    rhox = bsxfun(@times, pml_x, bsxfun(@times, pml_x, rhox) - dt .* rho0 .* duxdx);
    rhoy = bsxfun(@times, pml_y, bsxfun(@times, pml_y, rhoy) - dt .* rho0 .* duydy);
    rhoz = bsxfun(@times, pml_z, bsxfun(@times, pml_z, rhoz) - dt .* rho0 .* duzdz);

Definition at line 1060 of file KSpaceFirstOrder3DSolver.cpp.

◆ computeDensityNonliner()

void KSpaceFirstOrder3DSolver::computeDensityNonliner ( )
protected

Calculate new values of acoustic density for nonlinear case (rhoX, rhoy and rhoZ).

Matlab code:

   rho0_plus_rho = 2 .* (rhox + rhoy + rhoz) + rho0;
   rhox = bsxfun(@times, pml_x, bsxfun(@times, pml_x, rhox) - dt .* rho0_plus_rho .* duxdx);
   rhoy = bsxfun(@times, pml_y, bsxfun(@times, pml_y, rhoy) - dt .* rho0_plus_rho .* duydy);
   rhoz = bsxfun(@times, pml_z, bsxfun(@times, pml_z, rhoz) - dt .* rho0_plus_rho .* duzdz);

Definition at line 981 of file KSpaceFirstOrder3DSolver.cpp.

◆ computeInitialVelocityHeterogeneous()

void KSpaceFirstOrder3DSolver::computeInitialVelocityHeterogeneous ( )
protected

Matlab code:

   ux_sgx = dt ./ rho0_sgx .* ifft(ux_sgx).
   uy_sgy = dt ./ rho0_sgy .* ifft(uy_sgy).
   uz_sgz = dt ./ rho0_sgz .* ifft(uz_sgz).

Compute acoustic velocity for initial pressure problem.

Definition at line 1781 of file KSpaceFirstOrder3DSolver.cpp.

◆ computeInitialVelocityHomogeneousNonuniform()

void KSpaceFirstOrder3DSolver::computeInitialVelocityHomogeneousNonuniform ( )
protected

Matlab code:

    ux_sgx = dt ./ rho0_sgx .* dxudxn_sgx .* ifft(ux_sgx)
    uy_sgy = dt ./ rho0_sgy .* dyudxn_sgy .* ifft(uy_sgy)
    uz_sgz = dt ./ rho0_sgz .* dzudzn_sgz .* ifft(uz_sgz)

Compute acoustic velocity for initial pressure problem, homogenous medium, nonuniform grid.

Definition at line 1841 of file KSpaceFirstOrder3DSolver.cpp.

◆ computeInitialVelocityHomogeneousUniform()

void KSpaceFirstOrder3DSolver::computeInitialVelocityHomogeneousUniform ( )
protected

Matlab code:

    ux_sgx = dt ./ rho0_sgx .* ifft(ux_sgx).
    uy_sgy = dt ./ rho0_sgy .* ifft(uy_sgy).
    uz_sgz = dt ./ rho0_sgz .* ifft(uz_sgz).

Compute velocity for the initial pressure problem, homogeneous medium, uniform grid.

Definition at line 1813 of file KSpaceFirstOrder3DSolver.cpp.

◆ computeMainLoop()

void KSpaceFirstOrder3DSolver::computeMainLoop ( )
protected

Compute the main time loop of KSpaceFirstOrder3D.

Definition at line 610 of file KSpaceFirstOrder3DSolver.cpp.

◆ computePressureGradient()

void KSpaceFirstOrder3DSolver::computePressureGradient ( )
protected

Compute part of the new velocity term - gradient of pressure. Matlab code:

  bsxfun(\@times, ddx_k_shift_pos, kappa .* fftn(p))

Definition at line 2112 of file KSpaceFirstOrder3DSolver.cpp.

◆ computePressureLinear()

void KSpaceFirstOrder3DSolver::computePressureLinear ( )
protected

Compute new p for linear case.

Matlab code:

   case 'lossless'

       % calculate p using a linear adiabatic equation of state
       p = c.^2 .* (rhox + rhoy + rhoz);

   case 'absorbing'

       % calculate p using a linear absorbing equation of state
       p = c.^2 .* ( ...
           (rhox + rhoy + rhoz) ...
           + absorb_tau .* real(ifftn( absorb_nabla1 .* fftn(rho0 .* (duxdx + duydy + duzdz)) )) ...
           - absorb_eta .* real(ifftn( absorb_nabla2 .* fftn(rhox + rhoy + rhoz) )) ...
           );

Definition at line 1287 of file KSpaceFirstOrder3DSolver.cpp.

◆ computePressureNonlinear()

void KSpaceFirstOrder3DSolver::computePressureNonlinear ( )
protected

Compute acoustic pressure for non-linear case.

Matlab code:

   case 'lossless'
       % calculate p using a nonlinear adiabatic equation of state
       p = c.^2 .* (rhox + rhoy + rhoz + medium.BonA .* (rhox + rhoy + rhoz).^2 ./ (2 .* rho0));

   case 'absorbing'
       % calculate p using a nonlinear absorbing equation of state
       p = c.^2 .* (...
           (rhox + rhoy + rhoz) ...
           + absorb_tau .* real(ifftn( absorb_nabla1 .* fftn(rho0 .* (duxdx + duydy + duzdz)) ))...
           - absorb_eta .* real(ifftn( absorb_nabla2 .* fftn(rhox + rhoy + rhoz) ))...
           + medium.BonA .*(rhox + rhoy + rhoz).^2 ./ (2 .* rho0) ...
           );

Definition at line 1146 of file KSpaceFirstOrder3DSolver.cpp.

◆ computePressureTermsLinear()

void KSpaceFirstOrder3DSolver::computePressureTermsLinear ( RealMatrix densitySum,
RealMatrix velocityGradientSum 
)
protected
Parameters
[out]densitySum- rhoxSgx + rhoySgy + rhozSgz
[out]velocityGradientSum- rho0* (duxdx + duydy + duzdz)

Calculate two temporary sums in the new pressure formula, linear absorbing case.

Definition at line 2198 of file KSpaceFirstOrder3DSolver.cpp.

◆ computePressureTermsNonlinear()

template<bool bOnAScalarFlag, bool rho0ScalarFlag>
void KSpaceFirstOrder3DSolver::computePressureTermsNonlinear ( RealMatrix densitySum,
RealMatrix nonlinearTerm,
RealMatrix velocityGradientSum 
)
protected
Template Parameters
bOnAScalarFlag- is nonlinearity homogenous?
rho0ScalarFlag- is density homogeneous?
Parameters
[out]densitySum- rhoX + rhoY + rhoZ
[out]nonlinearTerm- BOnA + densitySum ^2 / 2 * rho0
[out]velocityGradientSum- rho0* (duxdx + duydy + duzdz)

Calculate three temporary sums in the new pressure formula non-linear absorbing case.

Definition at line 2154 of file KSpaceFirstOrder3DSolver.cpp.

◆ computeShiftedVelocity()

void KSpaceFirstOrder3DSolver::computeShiftedVelocity ( )
protected

Calculated shifted velocities.

Matlab code:

  ux_shifted = real(ifft(bsxfun(\@times, x_shift_neg, fft(ux_sgx, [], 1)), [], 1));
  uy_shifted = real(ifft(bsxfun(\@times, y_shift_neg, fft(uy_sgy, [], 2)), [], 2));
  uz_shifted = real(ifft(bsxfun(\@times, z_shift_neg, fft(uz_sgz, [], 3)), [], 3));

Definition at line 2428 of file KSpaceFirstOrder3DSolver.cpp.

◆ computeVelocity()

void KSpaceFirstOrder3DSolver::computeVelocity ( )
protected

Compute new values of acoustic velocity in all three dimensions (UxSgx, UySgy, UzSgz).

Matlab code:

  p_k = fftn(p);
  ux_sgx = bsxfun(@times, pml_x_sgx, ...
      bsxfun(@times, pml_x_sgx, ux_sgx) ...
      - dt .* rho0_sgx_inv .* real(ifftn( bsxfun(@times, ddx_k_shift_pos, kappa .* fftn(p)) )) ...
      );
  uy_sgy = bsxfun(@times, pml_y_sgy, ...
      bsxfun(@times, pml_y_sgy, uy_sgy) ...
      - dt .* rho0_sgy_inv .* real(ifftn( bsxfun(@times, ddy_k_shift_pos, kappa .* fftn(p)) )) ...
      );
  uz_sgz = bsxfun(@times, pml_z_sgz, ...
      bsxfun(@times, pml_z_sgz, uz_sgz) ...
      - dt .* rho0_sgz_inv .* real(ifftn( bsxfun(@times, ddz_k_shift_pos, kappa .* fftn(p)) )) ...
      );

Definition at line 864 of file KSpaceFirstOrder3DSolver.cpp.

◆ computeVelocityGradient()

void KSpaceFirstOrder3DSolver::computeVelocityGradient ( )
protected

Compute new values for duxdx, duydy, duzdz.

Definition at line 895 of file KSpaceFirstOrder3DSolver.cpp.

◆ computeVelocityHeterogeneous()

void KSpaceFirstOrder3DSolver::computeVelocityHeterogeneous ( )
protected

Matlab code:

  ux_sgx = bsxfun(@times, pml_x_sgx, bsxfun(@times, pml_x_sgx, ux_sgx) - dt .* rho0_sgx_inv .* real(ifftX)
  uy_sgy = bsxfun(@times, pml_y_sgy, bsxfun(@times, pml_y_sgy, uy_sgy) - dt .* rho0_sgy_inv .* real(ifftY)
  uz_sgz = bsxfun(@times, pml_z_sgz, bsxfun(@times, pml_z_sgz, uz_sgz) - dt .* rho0_sgz_inv .* real(ifftZ)

Compute acoustic velocity for heterogeneous medium and a uniform grid, x direction.

Definition at line 1884 of file KSpaceFirstOrder3DSolver.cpp.

◆ computeVelocityHomogeneousNonuniform()

void KSpaceFirstOrder3DSolver::computeVelocityHomogeneousNonuniform ( )
protected

Matlab code:

  ux_sgx = bsxfun(@times, pml_x_sgx, bsxfun(@times, pml_x_sgx, ux_sgx)  ...
                  - dt .* rho0_sgx_inv .* dxudxnSgx.* real(ifftX))
  uy_sgy = bsxfun(@times, pml_y_sgy, bsxfun(@times, pml_y_sgy, uy_sgy) ...
                  - dt .* rho0_sgy_inv .* dyudynSgy.* real(ifftY)
  uz_sgz = bsxfun(@times, pml_z_sgz, bsxfun(@times, pml_z_sgz, uz_sgz)
                  - dt .* rho0_sgz_inv .* dzudznSgz.* real(ifftZ)

Compute acoustic velocity for homogenous medium and nonuniform grid, x direction.

Definition at line 2031 of file KSpaceFirstOrder3DSolver.cpp.

◆ computeVelocityHomogeneousUniform()

void KSpaceFirstOrder3DSolver::computeVelocityHomogeneousUniform ( )
protected

Matlab code:

  ux_sgx = bsxfun(@times, pml_x_sgx, bsxfun(@times, pml_x_sgx, ux_sgx) - dt .* rho0_sgx_inv .* real(ifftX)
  uy_sgy = bsxfun(@times, pml_y_sgy, bsxfun(@times, pml_y_sgy, uy_sgy) - dt .* rho0_sgy_inv .* real(ifftY)
  uz_sgz = bsxfun(@times, pml_z_sgz, bsxfun(@times, pml_z_sgz, uz_sgz) - dt .* rho0_sgz_inv .* real(ifftZ)

Compute acoustic velocity for homogeneous medium and a uniform grid.

Definition at line 1959 of file KSpaceFirstOrder3DSolver.cpp.

◆ computeVelocitySourceTerm()

void KSpaceFirstOrder3DSolver::computeVelocitySourceTerm ( RealMatrix velocityMatrix,
const RealMatrix velocitySourceInput,
const IndexMatrix velocitySourceIndex 
)
protected
Parameters
[in]velocityMatrix- Velocity matrix to add the source in.
[in]velocitySourceInput- Source input to add.
[in]velocitySourceIndex- Source geometry index matrix.

Add in velocity source terms.

Definition at line 1362 of file KSpaceFirstOrder3DSolver.cpp.

◆ freeMemory()

void KSpaceFirstOrder3DSolver::freeMemory ( )
virtual

The method frees all memory allocated by the class.

Definition at line 122 of file KSpaceFirstOrder3DSolver.cpp.

◆ generateInitialDenisty()

void KSpaceFirstOrder3DSolver::generateInitialDenisty ( )
protected

Prepare dt./ rho0 for non-uniform grid.

Definition at line 1726 of file KSpaceFirstOrder3DSolver.cpp.

◆ generateKappa()

void KSpaceFirstOrder3DSolver::generateKappa ( )
protected

Generate kappa matrix for lossless medium.

Definition at line 1540 of file KSpaceFirstOrder3DSolver.cpp.

◆ generateKappaAndNablas()

void KSpaceFirstOrder3DSolver::generateKappaAndNablas ( )
protected

Generate kappa, absorb_nabla1, absorb_nabla2 for absorbing medium.

Definition at line 1591 of file KSpaceFirstOrder3DSolver.cpp.

◆ generateTauAndEta()

void KSpaceFirstOrder3DSolver::generateTauAndEta ( )
protected

Generate absorbTau and absorbEta in for heterogenous medium.

Definition at line 1658 of file KSpaceFirstOrder3DSolver.cpp.

◆ get1DIndex()

size_t KSpaceFirstOrder3DSolver::get1DIndex ( const size_t  z,
const size_t  y,
const size_t  x,
const DimensionSizes dimensionSizes 
)
inlineprotected
Parameters
[in]z- z coordinate
[in]y- y coordinate
[in]x- x coordinate
[in]dimensionSizes- Size of the matrix.
Returns

Definition at line 2704 of file KSpaceFirstOrder3DSolver.cpp.

◆ getAbsorbEta()

RealMatrix& KSpaceFirstOrder3DSolver::getAbsorbEta ( )
inlineprotected
Returns
Absorbing coefficient.

Definition at line 798 of file KSpaceFirstOrder3DSolver.h.

◆ getAbsorbNabla1()

RealMatrix& KSpaceFirstOrder3DSolver::getAbsorbNabla1 ( )
inlineprotected
Returns
Absorbing coefficient.

Definition at line 807 of file KSpaceFirstOrder3DSolver.h.

◆ getAbsorbNabla2()

RealMatrix& KSpaceFirstOrder3DSolver::getAbsorbNabla2 ( )
inlineprotected
Returns
Absorbing coefficient.

Definition at line 815 of file KSpaceFirstOrder3DSolver.h.

◆ getAbsorbTau()

RealMatrix& KSpaceFirstOrder3DSolver::getAbsorbTau ( )
inlineprotected
Returns
Absorbing coefficient.

Definition at line 790 of file KSpaceFirstOrder3DSolver.h.

◆ getBOnA()

RealMatrix& KSpaceFirstOrder3DSolver::getBOnA ( )
inlineprotected
Returns
Nonlinear coefficient.

Definition at line 782 of file KSpaceFirstOrder3DSolver.h.

◆ getC2()

RealMatrix& KSpaceFirstOrder3DSolver::getC2 ( )
inlineprotected
Returns
c^2 matrix.

Definition at line 451 of file KSpaceFirstOrder3DSolver.h.

◆ getCodeName()

std::string KSpaceFirstOrder3DSolver::getCodeName ( ) const
Returns
Release code version.

Get release code version.

Definition at line 370 of file KSpaceFirstOrder3DSolver.cpp.

◆ getCumulatedDataLoadTime()

double KSpaceFirstOrder3DSolver::getCumulatedDataLoadTime ( ) const
inline
Returns
Time to execute the simulation in seconds accumulated over all legs.

Definition at line 150 of file KSpaceFirstOrder3DSolver.h.

◆ getCumulatedPostProcessingTime()

double KSpaceFirstOrder3DSolver::getCumulatedPostProcessingTime ( ) const
inline
Returns
Time to post-processing simulation data in seconds accumulated over all legs.

Definition at line 160 of file KSpaceFirstOrder3DSolver.h.

◆ getCumulatedPreProcessingTime()

double KSpaceFirstOrder3DSolver::getCumulatedPreProcessingTime ( ) const
inline
Returns
Time to load data in seconds accumulated over all legs.

Definition at line 145 of file KSpaceFirstOrder3DSolver.h.

◆ getCumulatedSimulationTime()

double KSpaceFirstOrder3DSolver::getCumulatedSimulationTime ( ) const
inline
Returns
Time to execute the simulation in seconds accumulated over all legs.

Definition at line 155 of file KSpaceFirstOrder3DSolver.h.

◆ getCumulatedTotalTime()

double KSpaceFirstOrder3DSolver::getCumulatedTotalTime ( ) const
inline
Returns
Total execution time in seconds accumulated over all legs.

Definition at line 140 of file KSpaceFirstOrder3DSolver.h.

◆ getDataLoadTime()

double KSpaceFirstOrder3DSolver::getDataLoadTime ( ) const
inline
Returns
Time to load data in seconds.

Definition at line 124 of file KSpaceFirstOrder3DSolver.h.

◆ getDdxKShiftNeg()

ComplexMatrix& KSpaceFirstOrder3DSolver::getDdxKShiftNeg ( )
inlineprotected
Returns
Shift matrix.

Definition at line 630 of file KSpaceFirstOrder3DSolver.h.

◆ getDdxKShiftPos()

ComplexMatrix& KSpaceFirstOrder3DSolver::getDdxKShiftPos ( )
inlineprotected
Returns
Shift matrix.

Definition at line 606 of file KSpaceFirstOrder3DSolver.h.

◆ getDdyKShiftNeg()

ComplexMatrix& KSpaceFirstOrder3DSolver::getDdyKShiftNeg ( )
inlineprotected
Returns
Shift matrix.

Definition at line 638 of file KSpaceFirstOrder3DSolver.h.

◆ getDdyKShiftPos()

ComplexMatrix& KSpaceFirstOrder3DSolver::getDdyKShiftPos ( )
inlineprotected
Returns
Shift matrix.

Definition at line 614 of file KSpaceFirstOrder3DSolver.h.

◆ getDdzKShiftNeg()

ComplexMatrix& KSpaceFirstOrder3DSolver::getDdzKShiftNeg ( )
inlineprotected
Returns
shift matrix.

Definition at line 646 of file KSpaceFirstOrder3DSolver.h.

◆ getDdzKShiftPos()

ComplexMatrix& KSpaceFirstOrder3DSolver::getDdzKShiftPos ( )
inlineprotected
Returns
Shift matrix.

Definition at line 622 of file KSpaceFirstOrder3DSolver.h.

◆ getDelayMask()

IndexMatrix& KSpaceFirstOrder3DSolver::getDelayMask ( )
inlineprotected
Returns
delay mask.

Definition at line 857 of file KSpaceFirstOrder3DSolver.h.

◆ getDtRho0Sgx()

RealMatrix& KSpaceFirstOrder3DSolver::getDtRho0Sgx ( )
inlineprotected
Returns
Density matrix

Definition at line 547 of file KSpaceFirstOrder3DSolver.h.

◆ getDtRho0Sgy()

RealMatrix& KSpaceFirstOrder3DSolver::getDtRho0Sgy ( )
inlineprotected
Returns
Density matrix

Definition at line 555 of file KSpaceFirstOrder3DSolver.h.

◆ getDtRho0Sgz()

RealMatrix& KSpaceFirstOrder3DSolver::getDtRho0Sgz ( )
inlineprotected
Returns
Density matrix

Definition at line 563 of file KSpaceFirstOrder3DSolver.h.

◆ getDuxdx()

RealMatrix& KSpaceFirstOrder3DSolver::getDuxdx ( )
inlineprotected
Returns
Velocity gradient matrix.

Definition at line 521 of file KSpaceFirstOrder3DSolver.h.

◆ getDuydy()

RealMatrix& KSpaceFirstOrder3DSolver::getDuydy ( )
inlineprotected
Returns
Velocity gradient matrix.

Definition at line 529 of file KSpaceFirstOrder3DSolver.h.

◆ getDuzdz()

RealMatrix& KSpaceFirstOrder3DSolver::getDuzdz ( )
inlineprotected
Returns
Velocity gradient matrix.

Definition at line 537 of file KSpaceFirstOrder3DSolver.h.

◆ getDxudxn()

RealMatrix& KSpaceFirstOrder3DSolver::getDxudxn ( )
inlineprotected
Returns
Velocity matrix.

Definition at line 732 of file KSpaceFirstOrder3DSolver.h.

◆ getDxudxnSgx()

RealMatrix& KSpaceFirstOrder3DSolver::getDxudxnSgx ( )
inlineprotected
Returns
Velocity matrix.

Definition at line 756 of file KSpaceFirstOrder3DSolver.h.

◆ getDyudyn()

RealMatrix& KSpaceFirstOrder3DSolver::getDyudyn ( )
inlineprotected
Returns
Velocity matrix.

Definition at line 740 of file KSpaceFirstOrder3DSolver.h.

◆ getDyudynSgy()

RealMatrix& KSpaceFirstOrder3DSolver::getDyudynSgy ( )
inlineprotected
Returns
Velocity matrix.

Definition at line 764 of file KSpaceFirstOrder3DSolver.h.

◆ getDzudzn()

RealMatrix& KSpaceFirstOrder3DSolver::getDzudzn ( )
inlineprotected
Returns
Velocity matrix.

Definition at line 748 of file KSpaceFirstOrder3DSolver.h.

◆ getDzudznSgz()

RealMatrix& KSpaceFirstOrder3DSolver::getDzudznSgz ( )
inlineprotected
Returns
Velocity matrix.

Definition at line 772 of file KSpaceFirstOrder3DSolver.h.

◆ getInitialPressureSourceInput()

RealMatrix& KSpaceFirstOrder3DSolver::getInitialPressureSourceInput ( )
inlineprotected
Returns
Initial pressure source input data.

Definition at line 883 of file KSpaceFirstOrder3DSolver.h.

◆ getKappa()

RealMatrix& KSpaceFirstOrder3DSolver::getKappa ( )
inlineprotected
Returns
kappa matrix

Definition at line 442 of file KSpaceFirstOrder3DSolver.h.

◆ getMemoryUsage()

size_t KSpaceFirstOrder3DSolver::getMemoryUsage ( ) const
virtual
Returns
Memory consumed on the host side in MB.

Get peak memory usage.

Definition at line 340 of file KSpaceFirstOrder3DSolver.cpp.

◆ getP()

RealMatrix& KSpaceFirstOrder3DSolver::getP ( )
inlineprotected
Returns
Pressure matrix

Definition at line 460 of file KSpaceFirstOrder3DSolver.h.

◆ getPmlX()

RealMatrix& KSpaceFirstOrder3DSolver::getPmlX ( )
inlineprotected
Returns
PML matrix.

Definition at line 706 of file KSpaceFirstOrder3DSolver.h.

◆ getPmlXSgx()

RealMatrix& KSpaceFirstOrder3DSolver::getPmlXSgx ( )
inlineprotected
Returns
PML matrix.

Definition at line 681 of file KSpaceFirstOrder3DSolver.h.

◆ getPmlY()

RealMatrix& KSpaceFirstOrder3DSolver::getPmlY ( )
inlineprotected
Returns
PML matrix.

Definition at line 714 of file KSpaceFirstOrder3DSolver.h.

◆ getPmlYSgy()

RealMatrix& KSpaceFirstOrder3DSolver::getPmlYSgy ( )
inlineprotected
Returns
PML matrix.

Definition at line 689 of file KSpaceFirstOrder3DSolver.h.

◆ getPmlZ()

RealMatrix& KSpaceFirstOrder3DSolver::getPmlZ ( )
inlineprotected
Returns
PML matrix.

Definition at line 722 of file KSpaceFirstOrder3DSolver.h.

◆ getPmlZSgz()

RealMatrix& KSpaceFirstOrder3DSolver::getPmlZSgz ( )
inlineprotected
Returns
PML matrix.

Definition at line 697 of file KSpaceFirstOrder3DSolver.h.

◆ getPostProcessingTime()

double KSpaceFirstOrder3DSolver::getPostProcessingTime ( ) const
inline
Returns
Time to postprocess simulation data in seconds.

Definition at line 134 of file KSpaceFirstOrder3DSolver.h.

◆ getPreProcessingTime()

double KSpaceFirstOrder3DSolver::getPreProcessingTime ( ) const
inline
Returns
Pre-processing time in seconds.

Definition at line 119 of file KSpaceFirstOrder3DSolver.h.

◆ getPressureSourceIndex()

IndexMatrix& KSpaceFirstOrder3DSolver::getPressureSourceIndex ( )
inlineprotected
Returns
Source geometry indices

Definition at line 849 of file KSpaceFirstOrder3DSolver.h.

◆ getPressureSourceInput()

RealMatrix& KSpaceFirstOrder3DSolver::getPressureSourceInput ( )
inlineprotected
Returns
Pressure source input data.

Definition at line 875 of file KSpaceFirstOrder3DSolver.h.

◆ getRho0()

RealMatrix& KSpaceFirstOrder3DSolver::getRho0 ( )
inlineprotected
Returns
Density matrix.

Definition at line 596 of file KSpaceFirstOrder3DSolver.h.

◆ getRhoX()

RealMatrix& KSpaceFirstOrder3DSolver::getRhoX ( )
inlineprotected
Returns
Density matrix.

Definition at line 572 of file KSpaceFirstOrder3DSolver.h.

◆ getRhoY()

RealMatrix& KSpaceFirstOrder3DSolver::getRhoY ( )
inlineprotected
Returns
Density matrix.

Definition at line 580 of file KSpaceFirstOrder3DSolver.h.

◆ getRhoZ()

RealMatrix& KSpaceFirstOrder3DSolver::getRhoZ ( )
inlineprotected
Returns
Density matrix.

Definition at line 588 of file KSpaceFirstOrder3DSolver.h.

◆ getSensorMaskCorners()

IndexMatrix& KSpaceFirstOrder3DSolver::getSensorMaskCorners ( )
inlineprotected
Returns
Sensor mask data.

Definition at line 833 of file KSpaceFirstOrder3DSolver.h.

◆ getSensorMaskIndex()

IndexMatrix& KSpaceFirstOrder3DSolver::getSensorMaskIndex ( )
inlineprotected
Returns
Sensor mask data.

Definition at line 825 of file KSpaceFirstOrder3DSolver.h.

◆ getSimulationTime()

double KSpaceFirstOrder3DSolver::getSimulationTime ( ) const
inline
Returns
Time to execute the simulation in seconds.

Definition at line 129 of file KSpaceFirstOrder3DSolver.h.

◆ getTemp1Real3D()

RealMatrix& KSpaceFirstOrder3DSolver::getTemp1Real3D ( )
inlineprotected
Returns
Temporary real 3D matrix.

Definition at line 919 of file KSpaceFirstOrder3DSolver.h.

◆ getTemp2Real3D()

RealMatrix& KSpaceFirstOrder3DSolver::getTemp2Real3D ( )
inlineprotected
Returns
Temporary real 3D matrix.

Definition at line 927 of file KSpaceFirstOrder3DSolver.h.

◆ getTemp3Real3D()

RealMatrix& KSpaceFirstOrder3DSolver::getTemp3Real3D ( )
inlineprotected
Returns
Temporary real 3D matrix.

Definition at line 935 of file KSpaceFirstOrder3DSolver.h.

◆ getTempFftwShift()

FftwComplexMatrix& KSpaceFirstOrder3DSolver::getTempFftwShift ( )
inlineprotected
Returns
Temporary complex 3D matrix.

Definition at line 968 of file KSpaceFirstOrder3DSolver.h.

◆ getTempFftwX()

FftwComplexMatrix& KSpaceFirstOrder3DSolver::getTempFftwX ( )
inlineprotected
Returns
Temporary complex 3D matrix.

Definition at line 944 of file KSpaceFirstOrder3DSolver.h.

◆ getTempFftwY()

FftwComplexMatrix& KSpaceFirstOrder3DSolver::getTempFftwY ( )
inlineprotected
Returns
Temporary complex 3D matrix.

Definition at line 952 of file KSpaceFirstOrder3DSolver.h.

◆ getTempFftwZ()

FftwComplexMatrix& KSpaceFirstOrder3DSolver::getTempFftwZ ( )
inlineprotected
Returns
Temporary complex 3D matrix.

Definition at line 960 of file KSpaceFirstOrder3DSolver.h.

◆ getTotalTime()

double KSpaceFirstOrder3DSolver::getTotalTime ( ) const
inline
Returns
Total simulation time in seconds.

Definition at line 114 of file KSpaceFirstOrder3DSolver.h.

◆ getTransducerSourceInput()

RealMatrix& KSpaceFirstOrder3DSolver::getTransducerSourceInput ( )
inlineprotected
Returns
Transducer source input data.

Definition at line 867 of file KSpaceFirstOrder3DSolver.h.

◆ getUxSgx()

RealMatrix& KSpaceFirstOrder3DSolver::getUxSgx ( )
inlineprotected
Returns
Velocity matrix on staggered grid.

Definition at line 470 of file KSpaceFirstOrder3DSolver.h.

◆ getUxShifted()

RealMatrix& KSpaceFirstOrder3DSolver::getUxShifted ( )
inlineprotected
Returns
Unstaggeted velocity matrix.

Definition at line 495 of file KSpaceFirstOrder3DSolver.h.

◆ getUySgy()

RealMatrix& KSpaceFirstOrder3DSolver::getUySgy ( )
inlineprotected
Returns
Velocity matrix on staggered grid.

Definition at line 478 of file KSpaceFirstOrder3DSolver.h.

◆ getUyShifted()

RealMatrix& KSpaceFirstOrder3DSolver::getUyShifted ( )
inlineprotected
Returns
Unstaggered velocity matrix.

Definition at line 503 of file KSpaceFirstOrder3DSolver.h.

◆ getUzSgz()

RealMatrix& KSpaceFirstOrder3DSolver::getUzSgz ( )
inlineprotected
Returns
Velocity matrix on staggered grid.

Definition at line 486 of file KSpaceFirstOrder3DSolver.h.

◆ getUzShifted()

RealMatrix& KSpaceFirstOrder3DSolver::getUzShifted ( )
inlineprotected
Returns
Unstaggered velocity matrix.

Definition at line 511 of file KSpaceFirstOrder3DSolver.h.

◆ getVelocitySourceIndex()

IndexMatrix& KSpaceFirstOrder3DSolver::getVelocitySourceIndex ( )
inlineprotected
Returns
Source geometry indices

Definition at line 841 of file KSpaceFirstOrder3DSolver.h.

◆ GetVelocityXSourceInput()

RealMatrix& KSpaceFirstOrder3DSolver::GetVelocityXSourceInput ( )
inlineprotected
Returns
Velocity source input data.

Definition at line 892 of file KSpaceFirstOrder3DSolver.h.

◆ GetVelocityYSourceInput()

RealMatrix& KSpaceFirstOrder3DSolver::GetVelocityYSourceInput ( )
inlineprotected
Returns
Velocity source input data.

Definition at line 900 of file KSpaceFirstOrder3DSolver.h.

◆ getVelocityZSourceInput()

RealMatrix& KSpaceFirstOrder3DSolver::getVelocityZSourceInput ( )
inlineprotected
Returns
Velocity source input data.

Definition at line 908 of file KSpaceFirstOrder3DSolver.h.

◆ getXShiftNegR()

ComplexMatrix& KSpaceFirstOrder3DSolver::getXShiftNegR ( )
inlineprotected
Returns
Shift matrix.

Definition at line 655 of file KSpaceFirstOrder3DSolver.h.

◆ getYShiftNegR()

ComplexMatrix& KSpaceFirstOrder3DSolver::getYShiftNegR ( )
inlineprotected
Returns
Shift matrix.

Definition at line 663 of file KSpaceFirstOrder3DSolver.h.

◆ getZShiftNegR()

ComplexMatrix& KSpaceFirstOrder3DSolver::getZShiftNegR ( )
inlineprotected
Returns
Shift matrix.

Definition at line 671 of file KSpaceFirstOrder3DSolver.h.

◆ InitializeFftwPlans()

void KSpaceFirstOrder3DSolver::InitializeFftwPlans ( )
protected

Initialize FFTW plans.

Definition at line 463 of file KSpaceFirstOrder3DSolver.cpp.

◆ isCheckpointInterruption()

bool KSpaceFirstOrder3DSolver::isCheckpointInterruption ( ) const
protected
Returns
true if the simulation has been interrupted.

Was the loop interrupted to checkpoint?

Definition at line 2565 of file KSpaceFirstOrder3DSolver.cpp.

◆ isTimeToCheckpoint()

bool KSpaceFirstOrder3DSolver::isTimeToCheckpoint ( )
protected
Returns
true if it is time to interrupt the simulation and checkpoint.

Is time to checkpoint?

Definition at line 2551 of file KSpaceFirstOrder3DSolver.cpp.

◆ loadElapsedTimeFromOutputFile()

void KSpaceFirstOrder3DSolver::loadElapsedTimeFromOutputFile ( )
protected

Restore cumulated elapsed time from the output file.

Definition at line 2684 of file KSpaceFirstOrder3DSolver.cpp.

◆ loadInputData()

void KSpaceFirstOrder3DSolver::loadInputData ( )
virtual

If checkpointing is enabled, this may include reading data from checkpoint and output file.

Load data from the input file provided by the Parameter class and creates the output time series streams.

Definition at line 132 of file KSpaceFirstOrder3DSolver.cpp.

◆ postProcessing()

void KSpaceFirstOrder3DSolver::postProcessing ( )
protected

Post processing the quantities, closing the output streams and storing the sensor mask.

Definition at line 678 of file KSpaceFirstOrder3DSolver.cpp.

◆ preProcessing()

void KSpaceFirstOrder3DSolver::preProcessing ( )
protected

Initialize all indices, pre-compute constants such as c^2, rho0Sgx * dt and create kappa, absorbEta, absorbTau, absorbNabla1, absorbNabla2 matrices.

Compute pre-processing phase.

Definition at line 537 of file KSpaceFirstOrder3DSolver.cpp.

◆ printFullCodeNameAndLicense()

void KSpaceFirstOrder3DSolver::printFullCodeNameAndLicense ( ) const

Print full code name and the license.

Definition at line 380 of file KSpaceFirstOrder3DSolver.cpp.

◆ printStatistics()

void KSpaceFirstOrder3DSolver::printStatistics ( )
protected

Print progress statistics.

Definition at line 2515 of file KSpaceFirstOrder3DSolver.cpp.

◆ saveCheckpointData()

void KSpaceFirstOrder3DSolver::saveCheckpointData ( )
protected

Save checkpoint data into the checkpoint file, flush aggregated outputs into the output file.

Definition at line 773 of file KSpaceFirstOrder3DSolver.cpp.

◆ setProcessorAffinity()

void KSpaceFirstOrder3DSolver::setProcessorAffinity ( )
Warning
This may not work on some OS, it should be done by user before launching the code. Moreover, the user may want to change the thread placement, e.g. on NUMA systems. This the routine was disabled in ver 2.16

Set processor affinity.

Definition at line 434 of file KSpaceFirstOrder3DSolver.cpp.

◆ storeSensorData()

void KSpaceFirstOrder3DSolver::storeSensorData ( )
protected

Store sensor data.

Definition at line 719 of file KSpaceFirstOrder3DSolver.cpp.

◆ sumPressureTermsLinear()

template<bool c0ScalarFlag, bool areTauAndEtaScalars>
void KSpaceFirstOrder3DSolver::sumPressureTermsLinear ( const RealMatrix absorbTauTerm,
const RealMatrix absorbEtaTerm,
const RealMatrix densitySum 
)
protected
Template Parameters
c0ScalarFlag- is sound speed homogeneous?
areTauAndEtaScalars- are absorbTau and absorbEeta homogeneous?
Parameters
[in]absorbTauTerm- tau component
[in]absorbEtaTerm- eta component of the pressure term
[in]densitySum- Sum of three components of density (rhoXSgx + rhoYSgy + rhoZSgx)

Sum sub-terms to calculate new pressure, after FFTs, linear case.

Definition at line 2308 of file KSpaceFirstOrder3DSolver.cpp.

◆ sumPressureTermsLinearLossless()

void KSpaceFirstOrder3DSolver::sumPressureTermsLinearLossless ( )
protected

Sum sub-terms for new pressure, linear lossless case.

Definition at line 2383 of file KSpaceFirstOrder3DSolver.cpp.

◆ sumPressureTermsNonlinear()

template<bool c0ScalarFlag, bool areTauAndEtaScalars>
void KSpaceFirstOrder3DSolver::sumPressureTermsNonlinear ( const RealMatrix absorbTauTerm,
const RealMatrix absorbEtaTerm,
const RealMatrix nonlinearTerm 
)
protected

Sum sub-terms to calculate new pressure, after FFTs, non-linear case.

Template Parameters
c0ScalarFlag- is sound speed homogeneous?
areTauAndEtaScalars- are absorbTau and absorbEeta scalars
Parameters
[in]absorbTauTerm- tau component
[in]absorbEtaTerm- eta component of the pressure term
[in]nonlinearTerm- rho0 * (duxdx + duydy + duzdz)

Definition at line 2269 of file KSpaceFirstOrder3DSolver.cpp.

◆ sumPressureTermsNonlinearLossless()

template<bool c0ScalarFlag, bool nonlinearFlag, bool rho0ScalarFlag>
void KSpaceFirstOrder3DSolver::sumPressureTermsNonlinearLossless ( )
protected
Template Parameters
c0ScalarFlag- is sound speed homogeneous?
nonlinearFlag- is nonlinearity homogeneous?
rho0ScalarFlag- is density homogeneous?

Sum sub-terms for new p, nonlinear lossless case.

Definition at line 2347 of file KSpaceFirstOrder3DSolver.cpp.

◆ writeOutputDataInfo()

void KSpaceFirstOrder3DSolver::writeOutputDataInfo ( )
protected

Write statistics and the header into the output file.

Definition at line 736 of file KSpaceFirstOrder3DSolver.cpp.


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