GEOS
|
Public Types | |
using | EXEC_POLICY = parallelDevicePolicy< > |
Static Public Member Functions | |
static void | LeapFrogWithoutPML (real64 const dt, arrayView1d< real32 > const p_np1, arrayView1d< real32 > const p_n, arrayView1d< real32 > const p_nm1, arrayView1d< real32 const > const mass, arrayView1d< real32 > const stiffnessVector, arrayView1d< real32 const > const damping, arrayView1d< real32 > const rhs, arrayView1d< localIndex const > const freeSurfaceNodeIndicator, SortedArrayView< localIndex const > const solverTargetNodesSet) |
Apply second order Leap-Frog time scheme for isotropic case without PML. More... | |
static void | AttenuationLeapFrogWithoutPML (real64 const dt, arrayView1d< real32 > const p_np1, arrayView1d< real32 > const p_n, arrayView1d< real32 > const p_nm1, arrayView2d< real32 > const divpsi, arrayView1d< real32 const > const mass, arrayView1d< real32 > const stiffnessVector, arrayView1d< real32 > const stiffnessVectorA, arrayView1d< real32 const > const damping, arrayView1d< real32 > const rhs, arrayView1d< localIndex const > const freeSurfaceNodeIndicator, SortedArrayView< localIndex const > const solverTargetNodesSet, arrayView1d< real32 > referenceFrequencies, arrayView1d< real32 > anelasticityCoefficients) |
Apply second order Leap-Frog time scheme for isotropic case without PML, but with attenuation. More... | |
static void | LeapFrogforVTI (localIndex const size, real64 const dt, arrayView1d< real32 > const p_np1, arrayView1d< real32 > const p_n, arrayView1d< real32 > const p_nm1, arrayView1d< real32 > const q_np1, arrayView1d< real32 > const q_n, arrayView1d< real32 > const q_nm1, arrayView1d< real32 const > const mass, arrayView1d< real32 > const stiffnessVector_p, arrayView1d< real32 > const stiffnessVector_q, arrayView1d< real32 const > const damping_p, arrayView1d< real32 const > const damping_pq, arrayView1d< real32 const > const damping_q, arrayView1d< real32 const > const damping_qp, arrayView1d< real32 > const rhs, arrayView1d< localIndex const > const freeSurfaceNodeIndicator, arrayView1d< localIndex const > const lateralSurfaceNodeIndicator, arrayView1d< localIndex const > const bottomSurfaceNodeIndicator) |
Apply second order Leap-Frog time scheme for VTI case without PML. More... | |
static void | AttenuationLeapFrogforVTI (localIndex const size, real64 const dt, arrayView1d< real32 > const p_np1, arrayView1d< real32 > const p_n, arrayView1d< real32 > const p_nm1, arrayView1d< real32 > const q_np1, arrayView1d< real32 > const q_n, arrayView1d< real32 > const q_nm1, arrayView2d< real32 > const divpsi_p, arrayView2d< real32 > const divpsi_q, arrayView1d< real32 const > const mass, arrayView1d< real32 > const stiffnessVector_p, arrayView1d< real32 > const stiffnessVector_q, arrayView1d< real32 > const stiffnessVectorA_p, arrayView1d< real32 > const stiffnessVectorA_q, arrayView1d< real32 const > const damping_p, arrayView1d< real32 const > const damping_pq, arrayView1d< real32 const > const damping_q, arrayView1d< real32 const > const damping_qp, arrayView1d< real32 > const rhs, arrayView1d< localIndex const > const freeSurfaceNodeIndicator, arrayView1d< localIndex const > const lateralSurfaceNodeIndicator, arrayView1d< localIndex const > const bottomSurfaceNodeIndicator, arrayView1d< real32 > referenceFrequencies, arrayView1d< real32 > anelasticityCoefficients) |
Apply second order Leap-Frog time scheme for VTI case without PML, with attenuation. More... | |
Definition at line 26 of file AcousticTimeSchemeSEMKernel.hpp.
|
inlinestatic |
Apply second order Leap-Frog time scheme for VTI case without PML, with attenuation.
[in] | size | The number of nodes in the nodeManager |
[in] | dt | time-step |
[out] | p_np1 | pressure array at time n+1 (updated here) |
[in] | p_n | pressure array at time n |
[in] | p_nm1 | pressure array at time n-1 |
[out] | q_np1 | auxiliary pressure array at time n+1 (updated here) |
[in] | q_n | auxiliary pressure array at time n |
[in] | q_nm1 | auxiliary pressure array at time n-1 |
[in] | divpsi_p | divergence of the p-type acoustc memory variable |
[in] | divpsi_q | divergence of the q-type acoustc memory variable |
[in] | mass | the mass matrix |
[in] | stiffnessVector_p | array containing the product of the stiffness matrix R and the pressure at time n |
[in] | stiffnessVector_q | array containing the product of the stiffness matrix R and the auxiliary pressure at time n |
[in] | stiffnessVectorA_p | array containing the product of the attenuated stiffness matrix R and the pressure at time n |
[in] | stiffnessVectorA_q | array containing the product of the attenuated stiffness matrix R and the auxiliary pressure at time n |
[in] | damping_p | the damping matrix |
[in] | damping_pq | the damping matrix |
[in] | damping_q | the damping matrix |
[in] | damping_qp | the damping matrix |
[in] | rhs | the right-hand-side |
[in] | freeSurfaceNodeIndicator | array which contains indicators to tell if we are on a free-surface boundary or not |
[in] | lateralSurfaceNodeIndicator | array which contains indicators to tell if we are on a lateral boundary or not |
[in] | bottomSurfaceNodeIndicator | array which contains indicators to telle if we are on the bottom boundary or not |
[in] | referenceFrequencies | the reference frequencies for each SLS |
[in] | anelasticityCoefficients | the coefficients of anelasticity for each SLS |
Definition at line 249 of file AcousticTimeSchemeSEMKernel.hpp.
|
inlinestatic |
Apply second order Leap-Frog time scheme for isotropic case without PML, but with attenuation.
[in] | dt | time-step |
[out] | p_np1 | pressure array at time n+1 (updated here) |
[in] | p_n | pressure array at time n |
[in] | p_nm1 | pressure array at time n-1 |
[in] | divpsi | divergence of the acoustic memory variable |
[in] | mass | the mass matrix |
[in] | stiffnessVector | array containing the product of the stiffness matrix R and the pressure at time n |
[in] | stiffnessVectorA | array containing the product of the attenuation stiffness matrix R and the pressure at time n |
[in] | damping | the damping matrix |
[in] | rhs | the right-hand-side |
[in] | freeSurfaceNodeIndicator | array which contains indicators to tell if we are on a free-surface boundary or not |
[in] | solverTargetNodesSet | the targetted nodeset (useful in particular when we do elasto-acoustic simulation ) |
[in] | referenceFrequencies | the reference frequencies for each SLS |
[in] | anelasticityCoefficients | the coefficients of anelasticity for each SLS |
Definition at line 89 of file AcousticTimeSchemeSEMKernel.hpp.
|
inlinestatic |
Apply second order Leap-Frog time scheme for VTI case without PML.
[in] | size | The number of nodes in the nodeManager |
[in] | dt | time-step |
[out] | p_np1 | pressure array at time n+1 (updated here) |
[in] | p_n | pressure array at time n |
[in] | p_nm1 | pressure array at time n-1 |
[out] | q_np1 | auxiliary pressure array at time n+1 (updated here) |
[in] | q_n | auxiliary pressure array at time n |
[in] | q_nm1 | auxiliary pressure array at time n-1 |
[in] | mass | the mass matrix |
[in] | stiffnessVector_p | array containing the product of the stiffness matrix R and the pressure at time n |
[in] | stiffnessVector_q | array containing the product of the stiffness matrix R and the auxiliary pressure at time n |
[in] | damping_p | the damping matrix |
[in] | damping_pq | the damping matrix |
[in] | damping_q | the damping matrix |
[in] | damping_qp | the damping matrix |
[in] | rhs | the right-hand-side |
[in] | freeSurfaceNodeIndicator | array which contains indicators to tell if we are on a free-surface boundary or not |
[in] | lateralSurfaceNodeIndicator | array which contains indicators to tell if we are on a lateral boundary or not |
[in] | bottomSurfaceNodeIndicator | array which contains indicators to telle if we are on the bottom boundary or not |
Definition at line 152 of file AcousticTimeSchemeSEMKernel.hpp.
|
inlinestatic |
Apply second order Leap-Frog time scheme for isotropic case without PML.
[in] | dt | time-step |
[out] | p_np1 | pressure array at time n+1 (updated here) |
[in] | p_n | pressure array at time n |
[in] | p_nm1 | pressure array at time n-1 |
[in] | mass | the mass matrix |
[in] | stiffnessVector | array containing the product of the stiffness matrix R and the pressure at time n |
[in] | damping | the damping matrix |
[in] | rhs | the right-hand-side |
[in] | freeSurfaceNodeIndicator | array which contains indicators to tell if we are on a free-surface boundary or not |
[in] | solverTargetNodesSet | the targetted nodeset (useful in particular when we do elasto-acoustic simulation ) |
Definition at line 45 of file AcousticTimeSchemeSEMKernel.hpp.