GEOS
WaveSolverBase.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2016-2024 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2024 TotalEnergies
7  * Copyright (c) 2018-2024 The Board of Trustees of the Leland Stanford Junior University
8  * Copyright (c) 2023-2024 Chevron
9  * Copyright (c) 2019- GEOS/GEOSX Contributors
10  * All rights reserved
11  *
12  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
13  * ------------------------------------------------------------------------------------------------------------
14  */
15 
16 
21 #ifndef GEOS_PHYSICSSOLVERS_WAVEPROPAGATION_WAVESOLVERBASE_HPP_
22 #define GEOS_PHYSICSSOLVERS_WAVEPROPAGATION_WAVESOLVERBASE_HPP_
23 
24 
25 #include "mesh/MeshFields.hpp"
27 #include "common/LifoStorage.hpp"
29 #if !defined( GEOS_USE_HIP )
31 #endif
32 #include "WaveSolverUtils.hpp"
33 
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
41 #else
42 #define SEM_FE_TYPES
43 #endif
44 
45 #define SELECTED_FE_TYPES SEM_FE_TYPES
46 
47 namespace geos
48 {
49 
51 {
52 public:
53 
54  static constexpr real64 epsilonLoc = WaveSolverUtils::epsilonLoc;
55  using EXEC_POLICY = WaveSolverUtils::EXEC_POLICY;
56  using wsCoordType = WaveSolverUtils::wsCoordType;
57 
58  WaveSolverBase( const std::string & name,
59  Group * const parent );
60 
61  virtual ~WaveSolverBase() override;
62 
63  WaveSolverBase() = delete;
64  WaveSolverBase( WaveSolverBase const & ) = delete;
65  WaveSolverBase( WaveSolverBase && ) = default;
66 
67  WaveSolverBase & operator=( WaveSolverBase const & ) = delete;
68  WaveSolverBase & operator=( WaveSolverBase && ) = delete;
69 
70  virtual void initializePreSubGroups() override;
71 
72  virtual real64 solverStep( real64 const & time_n,
73  real64 const & dt,
74  integer const cycleNumber,
75  DomainPartition & domain ) override;
76 
77 
78  virtual real64 explicitStep( real64 const & time_n,
79  real64 const & dt,
80  integer const cycleNumber,
81  DomainPartition & domain ) override;
82 
84  {
85  static constexpr char const * sourceCoordinatesString() { return "sourceCoordinates"; }
86  static constexpr char const * sourceValueString() { return "sourceValue"; }
87 
88  static constexpr char const * timeSourceFrequencyString() { return "timeSourceFrequency"; }
89  static constexpr char const * timeSourceDelayString() { return "timeSourceDelay"; }
90  static constexpr char const * rickerOrderString() { return "rickerOrder"; }
91 
92  static constexpr char const * receiverCoordinatesString() { return "receiverCoordinates"; }
93 
94  static constexpr char const * sourceNodeIdsString() { return "sourceNodeIds"; }
95  static constexpr char const * sourceConstantsString() { return "sourceConstants"; }
96  static constexpr char const * sourceIsAccessibleString() { return "sourceIsAccessible"; }
97 
98  static constexpr char const * receiverNodeIdsString() { return "receiverNodeIds"; }
99  static constexpr char const * receiverConstantsString() {return "receiverConstants"; }
100  static constexpr char const * receiverIsLocalString() { return "receiverIsLocal"; }
101 
102  static constexpr char const * outputSeismoTraceString() { return "outputSeismoTrace"; }
103  static constexpr char const * dtSeismoTraceString() { return "dtSeismoTrace"; }
104  static constexpr char const * indexSeismoTraceString() { return "indexSeismoTrace"; }
105  static constexpr char const * forwardString() { return "forward"; }
106  static constexpr char const * saveFieldsString() { return "saveFields"; }
107  static constexpr char const * shotIndexString() { return "shotIndex"; }
108  static constexpr char const * enableLifoString() { return "enableLifo"; }
109  static constexpr char const * lifoSizeString() { return "lifoSize"; }
110  static constexpr char const * lifoOnDeviceString() { return "lifoOnDevice"; }
111  static constexpr char const * lifoOnHostString() { return "lifoOnHost"; }
112 
113  static constexpr char const * useDASString() { return "useDAS"; }
114  static constexpr char const * linearDASSamplesString() { return "linearDASSamples"; }
115  static constexpr char const * linearDASGeometryString() { return "linearDASGeometry"; }
116  static constexpr char const * linearDASVectorXString() { return "linearDASVectorX"; }
117  static constexpr char const * linearDASVectorYString() { return "linearDASVectorY"; }
118  static constexpr char const * linearDASVectorZString() { return "linearDASVectorZ"; }
119 
120  static constexpr char const * usePMLString() { return "usePML"; }
121  static constexpr char const * parametersPMLString() { return "parametersPML"; }
122 
123  static constexpr char const * receiverElemString() { return "receiverElem"; }
124  static constexpr char const * receiverRegionString() { return "receiverRegion"; }
125  static constexpr char const * freeSurfaceString() { return "FreeSurface"; }
126 
127  static constexpr char const * timestepStabilityLimitString() { return "timestepStabilityLimit"; }
128  static constexpr char const * timeStepString() { return "timeStep"; }
129 
130  static constexpr char const * attenuationTypeString() { return "attenuationType"; }
131  static constexpr char const * slsReferenceAngularFrequenciesString() { return "slsReferenceAngularFrequencies"; }
132  static constexpr char const * slsAnelasticityCoefficientsString() { return "slsAnelasticityCoefficients"; }
133 
134  static constexpr char const * sourceWaveletTableNames() { return "sourceWaveletTableNames"; }
135  };
136 
140  void reinit() override final;
141 
142  SortedArray< localIndex > const & getSolverNodesSet() { return m_solverTargetNodesSet; }
143 
144  void computeTargetNodeSet( arrayView2d< localIndex const, cells::NODE_MAP_USD > const & elemsToNodes,
145  localIndex const subRegionSize,
146  localIndex const numQuadraturePointsPerElem );
147 
148 protected:
149 
150  virtual void postInputInitialization() override;
151 
157  bool directoryExists( std::string const & directoryName );
158 
164  virtual void applyFreeSurfaceBC( real64 const time, DomainPartition & domain ) = 0;
165 
168  virtual real64 computeTimeStep( real64 & dtOut ) = 0;
169 
173  virtual void initializePML() = 0;
174 
175  virtual void incrementIndexSeismoTrace( real64 const time_n );
176 
187  virtual void computeAllSeismoTraces( real64 const time_n,
188  real64 const dt,
189  arrayView1d< real32 const > const var_np1,
190  arrayView1d< real32 const > const var_n,
191  arrayView2d< real32 > varAtReceivers,
192  arrayView1d< real32 > coeffs = {},
193  bool add = false );
202  virtual void compute2dVariableAllSeismoTraces( localIndex const regionIndex,
203  real64 const time_n,
204  real64 const dt,
205  arrayView2d< real32 const > const var_np1,
206  arrayView2d< real32 const > const var_n,
207  arrayView2d< real32 > varAtReceivers );
208 
214  virtual void applyPML( real64 const time, DomainPartition & domain ) = 0;
215 
222  virtual void precomputeSourceAndReceiverTerm( MeshLevel & baseMesh, MeshLevel & mesh, string_array const & regionNames ) = 0;
223 
233  virtual real64 explicitStepForward( real64 const & time_n,
234  real64 const & dt,
235  integer const cycleNumber,
236  DomainPartition & domain,
237  bool const computeGradient ) = 0;
247  virtual real64 explicitStepBackward( real64 const & time_n,
248  real64 const & dt,
249  integer const cycleNumber,
250  DomainPartition & domain,
251  bool const computeGradient ) = 0;
252 
253 
254  virtual void registerDataOnMesh( Group & meshBodies ) override;
255 
256  localIndex getNumNodesPerElem();
257 
260 
263 
266 
269 
272 
275 
278 
281 
284 
287 
290 
293 
296 
299 
302 
305 
308 
311 
314 
317 
320 
321  // Indicate the current shot computed for naming saved temporary data
322  integer m_shotIndex;
323 
326 
331 
332  //Time step computed with power iteration
333  real64 m_timeStep;
334 
337 
340 
343 
346 
349 
352 
355 
358 
361 
364 
367 
370 
372  std::unique_ptr< LifoStorage< real32, localIndex > > m_lifo;
373 
376 
379 
382 
385 
387  {
390 
393 
396 
399  R1Tensor32 thicknessMaxXYZPML;
400 
403  R1Tensor32 waveSpeedMaxXYZPML;
404  };
405 
406 };
407 
408 namespace fields
409 {
411 DECLARE_FIELD( referencePosition32,
412  "referencePosition32",
413  reference32Type,
414  0,
415  NOPLOT,
416  WRITE_AND_READ,
417  "Copy of the referencePosition from NodeManager in 32 bits integer" );
418 }
419 } /* namespace geos */
420 
421 #endif /* GEOS_PHYSICSSOLVERS_WAVEPROPAGATION_WAVESOLVERBASE_HPP_ */
#define DECLARE_FIELD(NAME, KEY, TYPE, DEFAULT, PLOTLEVEL, RESTARTFLAG, DESCRIPTION)
Generates a traits struct.
Definition: MeshFields.hpp:39
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.
Definition: MeshLevel.hpp:42
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.
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.
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.
virtual void precomputeSourceAndReceiverTerm(MeshLevel &baseMesh, MeshLevel &mesh, string_array const &regionNames)=0
Locate sources and receivers positions in the mesh elements, evaluate the basis functions at each poi...
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< 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.
array2d< real32 > m_sourceValue
Precomputed value of the source terms.
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
string_array m_sourceWaveletTableNames
Names of table functions for source wavelet (time dependency)
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.
Definition: DataTypes.hpp:180
float real32
32-bit floating point type.
Definition: DataTypes.hpp:97
Array< T, 2, PERMUTATION > array2d
Alias for 2D array.
Definition: DataTypes.hpp:192
std::string string
String type.
Definition: DataTypes.hpp:91
double real64
64-bit floating point type.
Definition: DataTypes.hpp:99
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:85
std::int32_t integer
Signed integer type.
Definition: DataTypes.hpp:82
LvArray::SortedArray< T, localIndex, LvArray::ChaiBuffer > SortedArray
A sorted array of local indices.
Definition: DataTypes.hpp:267
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Definition: DataTypes.hpp:196
std::vector< string > string_array
A 1-dimensional array of geos::string types.
Definition: DataTypes.hpp:393
Array< T, 1 > array1d
Alias for 1D array.
Definition: DataTypes.hpp:176
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.