| 
    GEOS
    
   | 
 
Public Types | |
| enum class | DASType : integer { none , dipole , strainIntegration } | 
| enum class | AttenuationType : integer { none , sls } | 
| using | EXEC_POLICY = parallelDevicePolicy<> | 
| using | wsCoordType = real32 | 
Static Public Member Functions | |
| static GEOS_HOST_DEVICE real32 | evaluateRicker (real64 const time_n, real32 const f0, real32 const t0, localIndex const order) | 
| static void | initTrace (char const *prefix, string const &name, bool const outputSeismoTrace, localIndex const nReceivers, arrayView1d< localIndex const > const receiverIsLocal) | 
| Initialize (clear) the trace file.  More... | |
| static void | writeSeismoTraceVector (char const *prefix, string const &name, bool const outputSeismoTrace, localIndex const nReceivers, arrayView1d< localIndex const > const receiverIsLocal, localIndex const nsamplesSeismoTrace, arrayView2d< real32 const > const varAtReceiversx, arrayView2d< real32 const > const varAtReceiversy, arrayView2d< real32 const > const varAtReceiversz) | 
| Convenient helper for 3D vectors calling 3 times the scalar version with only the sampled variable argument changed.  More... | |
| static void | writeSeismoTrace (char const *prefix, string const &name, bool const outputSeismoTrace, localIndex const nReceivers, arrayView1d< localIndex const > const receiverIsLocal, localIndex const nsamplesSeismoTrace, arrayView2d< real32 const > const varAtReceivers) | 
| Write the seismo traces to a file.  More... | |
| static void | computeSeismoTrace (real64 const time_n, real64 const dt, real64 const timeSeismo, localIndex const iSeismo, arrayView2d< localIndex const > const receiverNodeIds, arrayView2d< real64 const > const receiverConstants, arrayView1d< localIndex const > const receiverIsLocal, arrayView1d< real32 const > const var_np1, arrayView1d< real32 const > const var_n, arrayView2d< real32 > varAtReceivers, arrayView1d< real32 > coeffs={}, bool add=false) | 
| Compute the seismo traces.  More... | |
| static void | compute2dVariableSeismoTrace (real64 const time_n, real64 const dt, localIndex const regionIndex, arrayView1d< localIndex const > const receiverRegion, real64 const timeSeismo, localIndex const iSeismo, arrayView1d< localIndex const > const receiverElem, arrayView2d< real64 const > const receiverConstants, arrayView1d< localIndex const > const receiverIsLocal, arrayView2d< real32 const > const var_np1, arrayView2d< real32 const > const var_n, arrayView2d< real32 > varAtReceivers) | 
| Compute the seismo traces for 2d arrays.  More... | |
| static GEOS_HOST_DEVICE bool | locateSourceElement (real64 const numFacesPerElem, real64 const (&elemCenter)[3], arrayView2d< real64 const > const faceNormal, arrayView2d< real64 const > const faceCenter, arraySlice1d< localIndex const > const elemsToFaces, real64 const (&coords)[3]) | 
| Check if the source point is inside an element or not.  More... | |
| template<typename FE_TYPE > | |
| static GEOS_HOST_DEVICE void | computeCoordinatesOnReferenceElement (real64 const (&coords)[3], arraySlice1d< localIndex const, cells::NODE_MAP_USD - 1 > const elemsToNodes, arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const nodeCoords, real64(&coordsOnRefElem)[3]) | 
| Convert a mesh element point coordinate into a coordinate on the reference element.  More... | |
| static void | dotProduct (localIndex const size, arrayView1d< real32 > const &vector1, arrayView1d< real32 > const &vector2, real64 &res) | 
| Compute dotProduct between two vectors.  More... | |
| static GEOS_HOST_DEVICE R1Tensor | computeDASVector (real64 const dip, real64 const azimuth) | 
| Converts the DAS direction from dip/azimuth to a 3D unit vector.  More... | |
| template<typename REAL > | |
| static GEOS_HOST_DEVICE REAL | computeReferenceLengthForPenalty (arraySlice1d< localIndex const, cells::NODE_MAP_USD - 1 > const elemsToNodes, arrayView2d< REAL const, nodes::REFERENCE_POSITION_USD > const nodeCoords) | 
| Compute the reference size used for the tetrahederal DG penalty coefficient. This implementation uses the radius of the sphere inscribed to the tetrahedron.  More... | |
Static Public Attributes | |
| static constexpr real64 | epsilonLoc = 1e-8 | 
Definition at line 30 of file WaveSolverUtils.hpp.
      
  | 
  strong | 
| Enumerator | |
|---|---|
| none | deactivate attenuation (default)  | 
| sls | istandard-linear-solid description [Fichtner 2014]  | 
Definition at line 44 of file WaveSolverUtils.hpp.
      
  | 
  strong | 
| Enumerator | |
|---|---|
| none | deactivate DAS computation  | 
| dipole | use dipole formulation for DAS  | 
| strainIntegration | use strain integration for DAS  | 
Definition at line 37 of file WaveSolverUtils.hpp.
      
  | 
  inlinestatic | 
Compute the seismo traces for 2d arrays.
| [in] | time_n | Current time iteration | 
| [in] | dt | time-step | 
| [in] | regionIndex | Index of the current region | 
| [in] | receiverRegion | Array containing the region in which the receiver is located | 
| [in] | timeSeismo | time when the seismo is computed | 
| [in] | iSeismo | i-th seismo trace | 
| [in] | receiverElem | Array containing the element on which the receiver is located | 
| [in] | receiverConstants | constant part of the receiver term | 
| [in] | receiverIsLocal | flag indicating whether the receiver is local or not | 
| [in] | var_np1 | Array containing the variable at time n+1 | 
| [in] | var_n | Array containing the variable at time n | 
| [out] | varAtReceivers | Array containing the the variable computed at the receivers | 
Definition at line 270 of file WaveSolverUtils.hpp.
      
  | 
  inlinestatic | 
Convert a mesh element point coordinate into a coordinate on the reference element.
| FE_TYPE | finite element type | 
| [in] | coords | coordinate of the point | 
| [in] | elemsToNodes | element to node map for the base mesh | 
| [in] | nodeCoords | array of base mesh nodes coordinates | 
| [out] | coordsOnRefElem | to contain the coordinate computed in the reference element | 
Definition at line 373 of file WaveSolverUtils.hpp.
      
  | 
  inlinestatic | 
Converts the DAS direction from dip/azimuth to a 3D unit vector.
| [in] | dip | the dip of the linear DAS | 
| [in] | azimuth | the azimuth of the linear DAS | 
| [out] | a | unit vector pointing in the DAS direction | 
Definition at line 429 of file WaveSolverUtils.hpp.
      
  | 
  inlinestatic | 
Compute the reference size used for the tetrahederal DG penalty coefficient. This implementation uses the radius of the sphere inscribed to the tetrahedron.
| [in] | elemsToNodes | element to node map for the base mesh | 
| [in] | nodeCoords | array of base mesh nodes coordinates | 
Definition at line 448 of file WaveSolverUtils.hpp.
      
  | 
  inlinestatic | 
Compute the seismo traces.
| [in] | time_n | Current time iteration | 
| [in] | dt | time-step | 
| [in] | timeSeismo | time when the seismo is computed | 
| [in] | iSeismo | i-th seismo trace | 
| [in] | receiverNodeIds | indices of the nodes of the element where the receiver is located | 
| [in] | receiverConstants | constant part of the receiver term | 
| [in] | receiverIsLocal | flag indicating whether the receiver is local or not | 
| [in] | var_np1 | Array containing the variable at time n+1 | 
| [in] | var_n | Array containing the variable at time n | 
| [out] | varAtReceivers | Array containing the the variable computed at the receivers | 
| [in] | coeffs | Coefficients array for receivers | 
| [in] | add | Boolean to say if you want to add the value of interpolation to the same receiver coefficient or not | 
Definition at line 206 of file WaveSolverUtils.hpp.
      
  | 
  inlinestatic | 
Compute dotProduct between two vectors.
| numFacesPerElem | number of face on an element | 
| elemCenter | array containing the center of the elements | 
| faceNormal | array containing the normal of all faces | 
| faceCenter | array containing the center of all faces | 
| elemsToFaces | map to get the global faces from element index and local face index | 
| coords | coordinate of the point | 
Definition at line 409 of file WaveSolverUtils.hpp.
      
  | 
  inlinestatic | 
Initialize (clear) the trace file.
| [in] | prefix | Prefix of the output file | 
| [in] | name | Name of the solver on which you write the seismo trace | 
| [in] | outputSeismoTrace | Boolean equals to 1 if you want to output the seismotrace on a txt file 0 either | 
| [in] | nReceivers | Number of receivers | 
| [in] | receiverIsLocal | Array to check if the receiver is local to the MPI partition | 
Definition at line 92 of file WaveSolverUtils.hpp.
      
  | 
  inlinestatic | 
Check if the source point is inside an element or not.
| numFacesPerElem | number of face on an element | 
| elemCenter | array containing the center of the elements | 
| faceNormal | array containing the normal of all faces | 
| faceCenter | array containing the center of all faces | 
| elemsToFaces | map to get the global faces from element index and local face index | 
| coords | coordinate of the point | 
Definition at line 324 of file WaveSolverUtils.hpp.
      
  | 
  inlinestatic | 
Write the seismo traces to a file.
| [in] | prefix | Prefix of the output file | 
| [in] | name | Name of the solver on which you write the seismo trace | 
| [in] | outputSeismoTrace | Boolean equals to 1 if you want to output the seismotrace on a txt file 0 either | 
| [in] | nReceivers | Number of receivers | 
| [in] | receiverIsLocal | Array to check if the receiver is local to the MPI partition | 
| [in] | nsamplesSeismoTrace | Number of samples per seismo trace | 
| [in] | varAtReceivers | Array containing the the variable computed at the receivers | 
Definition at line 155 of file WaveSolverUtils.hpp.
      
  | 
  inlinestatic | 
Convenient helper for 3D vectors calling 3 times the scalar version with only the sampled variable argument changed.
| [in] | prefix | Prefix of the output file | 
| [in] | name | Name of the solver on which you write the seismo trace | 
| [in] | outputSeismoTrace | Boolean equals to 1 if you want to output the seismotrace on a txt file 0 either | 
| [in] | nReceivers | Number of receivers | 
| [in] | receiverIsLocal | Array to check if the receiver is local to the MPI partition | 
| [in] | nsamplesSeismoTrace | Number of samples per seismo trace | 
| [out] | varAtReceiversx | Array containing the variable (x-direction) computed at the receivers | 
| [out] | varAtReceiversy | Array containing the variable (y-direction) computed at the receivers | 
| [out] | varAtReceiversz | Array containing the variable (z-direction) computed at the receivers | 
Definition at line 130 of file WaveSolverUtils.hpp.