21 #ifndef GEOS_PHYSICSSOLVERS_WAVEPROPAGATION_WAVESOLVERBASE_HPP_
22 #define GEOS_PHYSICSSOLVERS_WAVEPROPAGATION_WAVESOLVERBASE_HPP_
27 #include "common/LifoStorage.hpp"
29 #if !defined( GEOS_USE_HIP )
34 #if !defined( GEOS_USE_HIP )
35 #define SEM_FE_TYPES \
36 finiteElement::Q1_Hexahedron_Lagrange_GaussLobatto, \
37 finiteElement::Q2_Hexahedron_Lagrange_GaussLobatto, \
38 finiteElement::Q3_Hexahedron_Lagrange_GaussLobatto, \
39 finiteElement::Q4_Hexahedron_Lagrange_GaussLobatto, \
40 finiteElement::Q5_Hexahedron_Lagrange_GaussLobatto
45 #define SELECTED_FE_TYPES SEM_FE_TYPES
54 static constexpr
real64 epsilonLoc = WaveSolverUtils::epsilonLoc;
55 using EXEC_POLICY = WaveSolverUtils::EXEC_POLICY;
56 using wsCoordType = WaveSolverUtils::wsCoordType;
59 Group *
const parent );
85 static constexpr
char const * sourceCoordinatesString() {
return "sourceCoordinates"; }
87 static constexpr
char const * timeSourceFrequencyString() {
return "timeSourceFrequency"; }
88 static constexpr
char const * timeSourceDelayString() {
return "timeSourceDelay"; }
89 static constexpr
char const * rickerOrderString() {
return "rickerOrder"; }
91 static constexpr
char const * receiverCoordinatesString() {
return "receiverCoordinates"; }
93 static constexpr
char const * sourceNodeIdsString() {
return "sourceNodeIds"; }
94 static constexpr
char const * sourceConstantsString() {
return "sourceConstants"; }
95 static constexpr
char const * sourceIsAccessibleString() {
return "sourceIsAccessible"; }
97 static constexpr
char const * receiverNodeIdsString() {
return "receiverNodeIds"; }
98 static constexpr
char const * receiverConstantsString() {
return "receiverConstants"; }
99 static constexpr
char const * receiverIsLocalString() {
return "receiverIsLocal"; }
101 static constexpr
char const * outputSeismoTraceString() {
return "outputSeismoTrace"; }
102 static constexpr
char const * dtSeismoTraceString() {
return "dtSeismoTrace"; }
103 static constexpr
char const * indexSeismoTraceString() {
return "indexSeismoTrace"; }
104 static constexpr
char const * forwardString() {
return "forward"; }
105 static constexpr
char const * saveFieldsString() {
return "saveFields"; }
106 static constexpr
char const * shotIndexString() {
return "shotIndex"; }
107 static constexpr
char const * enableLifoString() {
return "enableLifo"; }
108 static constexpr
char const * lifoSizeString() {
return "lifoSize"; }
109 static constexpr
char const * lifoOnDeviceString() {
return "lifoOnDevice"; }
110 static constexpr
char const * lifoOnHostString() {
return "lifoOnHost"; }
112 static constexpr
char const * useDASString() {
return "useDAS"; }
113 static constexpr
char const * linearDASSamplesString() {
return "linearDASSamples"; }
114 static constexpr
char const * linearDASGeometryString() {
return "linearDASGeometry"; }
115 static constexpr
char const * linearDASVectorXString() {
return "linearDASVectorX"; }
116 static constexpr
char const * linearDASVectorYString() {
return "linearDASVectorY"; }
117 static constexpr
char const * linearDASVectorZString() {
return "linearDASVectorZ"; }
119 static constexpr
char const * usePMLString() {
return "usePML"; }
120 static constexpr
char const * parametersPMLString() {
return "parametersPML"; }
122 static constexpr
char const * receiverElemString() {
return "receiverElem"; }
123 static constexpr
char const * receiverRegionString() {
return "receiverRegion"; }
124 static constexpr
char const * freeSurfaceString() {
return "FreeSurface"; }
126 static constexpr
char const * timestepStabilityLimitString() {
return "timestepStabilityLimit"; }
127 static constexpr
char const * timeStepString() {
return "timeStep"; }
129 static constexpr
char const * attenuationTypeString() {
return "attenuationType"; }
130 static constexpr
char const * slsReferenceAngularFrequenciesString() {
return "slsReferenceAngularFrequencies"; }
131 static constexpr
char const * slsAnelasticityCoefficientsString() {
return "slsAnelasticityCoefficients"; }
133 static constexpr
char const * sourceWaveletTableNames() {
return "sourceWaveletTableNames"; }
145 localIndex const numQuadraturePointsPerElem );
174 virtual void incrementIndexSeismoTrace(
real64 const time_n );
236 bool const computeGradient ) = 0;
250 bool const computeGradient ) = 0;
368 std::unique_ptr< LifoStorage< real32, localIndex > >
m_lifo;
408 "referencePosition32",
413 "Copy of the referencePosition from NodeManager in 32 bits integer" );
#define DECLARE_FIELD(NAME, KEY, TYPE, DEFAULT, PLOTLEVEL, RESTARTFLAG, DESCRIPTION)
Generates a traits struct.
Partition of the decomposed physical domain. It also manages the connexion information to its neighbo...
Class facilitating the representation of a multi-level discretization of a MeshBody.
Base class for all physics solvers.
localIndex m_nsamplesSeismoTrace
Amount of seismoTrace that will be recorded for each receiver.
array1d< real32 > m_linearDASVectorX
X component of the linear DAS direction vector.
real32 m_timeSourceFrequency
Central frequency for the Ricker time source.
virtual void precomputeSourceAndReceiverTerm(MeshLevel &baseMesh, MeshLevel &mesh, arrayView1d< string const > const ®ionNames)=0
Locate sources and receivers positions in the mesh elements, evaluate the basis functions at each poi...
integer m_usePML
Flag to apply PML.
virtual real64 explicitStepForward(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain, bool const computeGradient)=0
Perform forward explicit step.
array1d< localIndex > m_receiverIsLocal
Flag that indicates whether the receiver is local or not to the MPI rank.
std::unique_ptr< LifoStorage< real32, localIndex > > m_lifo
LIFO to store p_dt2.
virtual void initializePML()=0
Initialize Perfectly Matched Layer (PML) information.
array1d< localIndex > m_sourceIsAccessible
Flag that indicates whether the source is local or not to the MPI rank.
virtual void applyPML(real64 const time, DomainPartition &domain)=0
Apply Perfectly Matched Layer (PML) to the regions defined in the geometry box from the xml.
virtual void computeAllSeismoTraces(real64 const time_n, real64 const dt, arrayView1d< real32 const > const var_np1, arrayView1d< real32 const > const var_n, arrayView2d< real32 > varAtReceivers, arrayView1d< real32 > coeffs={}, bool add=false)
Computes the traces on all receivers (see @computeSeismoTraces) up to time_n+dt.
localIndex m_lifoSize
lifo size (should be the total number of buffer to save in LIFO)
localIndex m_lifoOnDevice
Number of buffers to store on device by LIFO (if negative, opposite of percentage of remaining memory...
virtual real64 explicitStep(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain) override
Entry function for an explicit time integration step.
array2d< real64 > m_linearDASGeometry
Geometry parameters for a linear DAS fiber (dip, azimuth, gauge length)
SortedArray< localIndex > m_solverTargetNodesSet
A set of target nodes IDs that will be handled by the current solver.
integer m_timestepStabilityLimit
virtual void registerDataOnMesh(Group &meshBodies) override
Register wrappers that contain data on the mesh objects.
array2d< real64 > m_sourceCoordinates
Coordinates of the sources in the mesh.
localIndex m_saveFields
Indicate if we want to save fields to restore them during backward.
localIndex m_enableLifo
Flag to enable LIFO.
localIndex m_indexSeismoTrace
Cycle number for output SeismoTrace.
virtual void initializePreSubGroups() override
Called by Initialize() prior to initializing sub-Groups.
array1d< localIndex > m_receiverElem
Array containing the elements which contain a receiver.
array1d< real32 > m_slsReferenceAngularFrequencies
Vector containing the reference frequencies for the standard-linear-solid (SLS) anelasticity model.
void reinit() override final
Re-initialize source and receivers positions in the mesh, and resize the pressureNp1_at_receivers arr...
integer m_linearDASSamples
Number of points used for strain integration for dipole DAS.
array1d< string > m_sourceWaveletTableNames
Names of table functions for source wavelet (time dependency)
array1d< real32 > m_slsAnelasticityCoefficients
Vector containing the anelasticity coefficients for the standard-linear-solid (SLS) anelasticity mode...
localIndex m_rickerOrder
Flag that indicates the order of the Ricker to be used, order 2 by default.
virtual void postInputInitialization() override
WaveSolverUtils::AttenuationType m_attenuationType
Flag to indicate which attenuation type will be modeled.
real64 m_dtSeismoTrace
Time step for seismoTrace output.
localIndex m_lifoOnHost
Number of buffers to store on host by LIFO (if negative, opposite of percentage of remaining memory)
integer m_outputSeismoTrace
Flag that indicates if we write the seismo trace in a file .txt, 0 no output, 1 otherwise.
array1d< real32 > m_linearDASVectorZ
Z component of the linear DAS direction vector.
array2d< localIndex > m_receiverNodeIds
Indices of the element nodes (in the right order) for each receiver point.
array2d< real64 > m_sourceConstants
Constant part of the source for the nodes listed in m_sourceNodeIds.
array2d< real64 > m_receiverCoordinates
Coordinates of the receivers in the mesh.
WaveSolverUtils::DASType m_useDAS
Flag to indicate which DAS type will be modeled.
virtual void applyFreeSurfaceBC(real64 const time, DomainPartition &domain)=0
Apply free surface condition to the face defined in the geometry box of the xml.
virtual real64 solverStep(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain) override
entry function to perform a solver step
bool directoryExists(std::string const &directoryName)
Utility function to check if a directory exists.
array2d< localIndex > m_sourceNodeIds
Indices of the nodes (in the right order) for each source point.
localIndex m_forward
Indicate if we want to compute forward ou backward.
real32 m_timeSourceDelay
Source time delay (1 / f0 by default)
virtual void compute2dVariableAllSeismoTraces(localIndex const regionIndex, real64 const time_n, real64 const dt, arrayView2d< real32 const > const var_np1, arrayView2d< real32 const > const var_n, arrayView2d< real32 > varAtReceivers)
Computes the traces on all receivers (see @computeSeismoTraces) up to time_n+dt.
bool m_useSourceWaveletTables
Flag to indicate if source wavelet table functions are used.
array2d< real64 > m_receiverConstants
Basis function evaluated at the receiver for the nodes listed in m_receiverNodeIds.
array1d< localIndex > m_receiverRegion
Array containing the elements which contain the region which the receiver belongs.
array1d< real32 > m_linearDASVectorY
Y component of the linear DAS direction vector.
virtual real64 explicitStepBackward(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain, bool const computeGradient)=0
Perform backward explicit step.
array1d< TableFunction::KernelWrapper > m_sourceWaveletTableWrappers
Wrappers of table functions for source wavelet (time dependency)
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
float real32
32-bit floating point type.
Array< T, 2, PERMUTATION > array2d
Alias for 2D array.
std::string string
String type.
double real64
64-bit floating point type.
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
std::int32_t integer
Signed integer type.
LvArray::SortedArray< T, localIndex, LvArray::ChaiBuffer > SortedArray
A sorted array of local indices.
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Array< T, 1 > array1d
Alias for 1D array.
Structure to hold scoped key names.
R1Tensor32 xMaxPML
Maximum (x,y,z) coordinates of inner PML boundaries.
R1Tensor32 xMinPML
Mininum (x,y,z) coordinates of inner PML boundaries.
R1Tensor32 waveSpeedMinXYZPML
Wave speed in the PML region, used to compute the damping profile.
R1Tensor32 thicknessMinXYZPML
Thickness of the PML region, used to compute the damping profile.
real32 reflectivityPML
Desired reflectivity of the PML region, used to compute the damping profile.