GEOS
Public Types | Static Public Member Functions | Static Public Attributes | List of all members
geos::WaveSolverUtils Struct Reference

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...
 

Static Public Attributes

static constexpr real64 epsilonLoc = 1e-8
 

Detailed Description

Definition at line 31 of file WaveSolverUtils.hpp.

Member Enumeration Documentation

◆ AttenuationType

Enumerator
none 

deactivate attenuation (default)

sls 

istandard-linear-solid description [Fichtner 2014]

Definition at line 45 of file WaveSolverUtils.hpp.

◆ DASType

Enumerator
none 

deactivate DAS computation

dipole 

use dipole formulation for DAS

strainIntegration 

use strain integration for DAS

Definition at line 38 of file WaveSolverUtils.hpp.

Member Function Documentation

◆ compute2dVariableSeismoTrace()

static void geos::WaveSolverUtils::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 
)
inlinestatic

Compute the seismo traces for 2d arrays.

Parameters
[in]time_nCurrent time iteration
[in]dttime-step
[in]regionIndexIndex of the current region
[in]receiverRegionArray containing the region in which the receiver is located
[in]timeSeismotime when the seismo is computed
[in]iSeismoi-th seismo trace
[in]receiverElemArray containing the element on which the receiver is located
[in]receiverConstantsconstant part of the receiver term
[in]receiverIsLocalflag indicating whether the receiver is local or not
[in]var_np1Array containing the variable at time n+1
[in]var_nArray containing the variable at time n
[out]varAtReceiversArray containing the the variable computed at the receivers

Definition at line 270 of file WaveSolverUtils.hpp.

◆ computeCoordinatesOnReferenceElement()

template<typename FE_TYPE >
static GEOS_HOST_DEVICE void geos::WaveSolverUtils::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] 
)
inlinestatic

Convert a mesh element point coordinate into a coordinate on the reference element.

Template Parameters
FE_TYPEfinite element type
Parameters
[in]coordscoordinate of the point
[in]elemsToNodeselement to node map for the base mesh
[in]nodeCoordsarray of base mesh nodes coordinates
[out]coordsOnRefElemto contain the coordinate computed in the reference element

Definition at line 374 of file WaveSolverUtils.hpp.

◆ computeDASVector()

static GEOS_HOST_DEVICE R1Tensor geos::WaveSolverUtils::computeDASVector ( real64 const  dip,
real64 const  azimuth 
)
inlinestatic

Converts the DAS direction from dip/azimuth to a 3D unit vector.

Parameters
[in]dipthe dip of the linear DAS
[in]azimuththe azimuth of the linear DAS
[out]aunit vector pointing in the DAS direction

Definition at line 435 of file WaveSolverUtils.hpp.

◆ computeSeismoTrace()

static void geos::WaveSolverUtils::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 
)
inlinestatic

Compute the seismo traces.

Parameters
[in]time_nCurrent time iteration
[in]dttime-step
[in]timeSeismotime when the seismo is computed
[in]iSeismoi-th seismo trace
[in]receiverNodeIdsindices of the nodes of the element where the receiver is located
[in]receiverConstantsconstant part of the receiver term
[in]receiverIsLocalflag indicating whether the receiver is local or not
[in]var_np1Array containing the variable at time n+1
[in]var_nArray containing the variable at time n
[out]varAtReceiversArray containing the the variable computed at the receivers
[in]coeffsCoefficients array for receivers
[in]addBoolean 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.

◆ dotProduct()

static void geos::WaveSolverUtils::dotProduct ( localIndex const  size,
arrayView1d< real32 > const &  vector1,
arrayView1d< real32 > const &  vector2,
real64 res 
)
inlinestatic

Compute dotProduct between two vectors.

Parameters
numFacesPerElemnumber of face on an element
elemCenterarray containing the center of the elements
faceNormalarray containing the normal of all faces
faceCenterarray containing the center of all faces
elemsToFacesmap to get the global faces from element index and local face index
coordscoordinate of the point
Returns
true if coords is inside the element

Definition at line 411 of file WaveSolverUtils.hpp.

◆ initTrace()

static void geos::WaveSolverUtils::initTrace ( char const *  prefix,
string const &  name,
bool const  outputSeismoTrace,
localIndex const  nReceivers,
arrayView1d< localIndex const > const  receiverIsLocal 
)
inlinestatic

Initialize (clear) the trace file.

Parameters
[in]prefixPrefix of the output file
[in]nameName of the solver on which you write the seismo trace
[in]outputSeismoTraceBoolean equals to 1 if you want to output the seismotrace on a txt file 0 either
[in]nReceiversNumber of receivers
[in]receiverIsLocalArray to check if the receiver is local to the MPI partition

Definition at line 94 of file WaveSolverUtils.hpp.

◆ locateSourceElement()

static GEOS_HOST_DEVICE bool geos::WaveSolverUtils::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] 
)
inlinestatic

Check if the source point is inside an element or not.

Parameters
numFacesPerElemnumber of face on an element
elemCenterarray containing the center of the elements
faceNormalarray containing the normal of all faces
faceCenterarray containing the center of all faces
elemsToFacesmap to get the global faces from element index and local face index
coordscoordinate of the point
Returns
true if coords is inside the element

Definition at line 324 of file WaveSolverUtils.hpp.

◆ writeSeismoTrace()

static void geos::WaveSolverUtils::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 
)
inlinestatic

Write the seismo traces to a file.

Parameters
[in]prefixPrefix of the output file
[in]nameName of the solver on which you write the seismo trace
[in]outputSeismoTraceBoolean equals to 1 if you want to output the seismotrace on a txt file 0 either
[in]nReceiversNumber of receivers
[in]receiverIsLocalArray to check if the receiver is local to the MPI partition
[in]nsamplesSeismoTraceNumber of samples per seismo trace
[in]varAtReceiversArray containing the the variable computed at the receivers

Definition at line 156 of file WaveSolverUtils.hpp.

◆ writeSeismoTraceVector()

static void geos::WaveSolverUtils::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 
)
inlinestatic

Convenient helper for 3D vectors calling 3 times the scalar version with only the sampled variable argument changed.

Parameters
[in]prefixPrefix of the output file
[in]nameName of the solver on which you write the seismo trace
[in]outputSeismoTraceBoolean equals to 1 if you want to output the seismotrace on a txt file 0 either
[in]nReceiversNumber of receivers
[in]receiverIsLocalArray to check if the receiver is local to the MPI partition
[in]nsamplesSeismoTraceNumber of samples per seismo trace
[out]varAtReceiversxArray containing the variable (x-direction) computed at the receivers
[out]varAtReceiversyArray containing the variable (y-direction) computed at the receivers
[out]varAtReceiverszArray containing the variable (z-direction) computed at the receivers

Definition at line 131 of file WaveSolverUtils.hpp.


The documentation for this struct was generated from the following file: