GEOSX
Classes | Public Types | Public Member Functions | Protected Attributes | List of all members
geos::LaplaceFEMKernel< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE > Class Template Reference

Implements kernels for solving Laplace's equation. More...

#include <LaplaceFEMKernels.hpp>

Inheritance diagram for geos::LaplaceFEMKernel< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >:
Inheritance graph
[legend]

Classes

class  StackVariables
 Kernel variables allocated on the stack. More...
 

Public Types

using Base = finiteElement::ImplicitKernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, 1, 1 >
 An alias for the base class.
 
- Public Types inherited from geos::finiteElement::ImplicitKernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, 1, 1 >
using Base = KernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, NUM_DOF_PER_TEST_SP, NUM_DOF_PER_TRIAL_SP >
 Alias for the base class. (i.e. geos::finiteElement::KernelBase)
 

Public Member Functions

 LaplaceFEMKernel (NodeManager const &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, localIndex const targetRegionIndex, SUBREGION_TYPE const &elementSubRegion, FE_TYPE const &finiteElementSpace, CONSTITUTIVE_TYPE &inputConstitutiveType, arrayView1d< globalIndex const > const inputDofNumber, globalIndex const rankOffset, CRSMatrixView< real64, globalIndex const > const inputMatrix, arrayView1d< real64 > const inputRhs, real64 const inputDt, string const fieldName)
 Constructor. More...
 
GEOS_HOST_DEVICE void setup (localIndex const k, StackVariables &stack) const
 Copy global values from primary field to a local stack array. More...
 
GEOS_HOST_DEVICE void quadraturePointKernel (localIndex const k, localIndex const q, StackVariables &stack) const
 Performs a state update at a quadrature point. More...
 
GEOS_HOST_DEVICE real64 complete (localIndex const k, StackVariables &stack) const
 Performs the complete phase for the kernel. More...
 
- Public Member Functions inherited from geos::finiteElement::ImplicitKernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, 1, 1 >
 ImplicitKernelBase (NodeManager const &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, localIndex const targetRegionIndex, SUBREGION_TYPE const &elementSubRegion, FE_TYPE const &finiteElementSpace, CONSTITUTIVE_TYPE &inputConstitutiveType, arrayView1d< globalIndex const > const &inputDofNumber, globalIndex const rankOffset, CRSMatrixView< real64, globalIndex const > const &inputMatrix, arrayView1d< real64 > const &inputRhs, real64 const inputDt)
 Constructor. More...
 
GEOS_HOST_DEVICE void setup (localIndex const k, StackVariables &stack) const
 Performs the setup phase for the kernel. More...
 
- Public Member Functions inherited from geos::finiteElement::KernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, NUM_DOF_PER_TEST_SP, NUM_DOF_PER_TRIAL_SP >
 KernelBase (SUBREGION_TYPE const &elementSubRegion, FE_TYPE const &finiteElementSpace, CONSTITUTIVE_TYPE &inputConstitutiveType)
 Constructor. More...
 
GEOS_HOST_DEVICE void setup (localIndex const k, StackVariables &stack) const
 Performs the setup phase for the kernel. More...
 
GEOS_HOST_DEVICE GEOS_FORCE_INLINE void quadraturePointKernel (localIndex const k, localIndex const q, StackVariables &stack) const
 Performs a state update at a quadrature point. More...
 
GEOS_HOST_DEVICE real64 complete (localIndex const k, StackVariables &stack) const
 Performs the complete phase for the kernel. More...
 

Protected Attributes

const arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > m_X
 The array containing the nodal position array.
 
const arrayView1d< real64 const > m_primaryField
 The global primary field array.
 
- Protected Attributes inherited from geos::finiteElement::ImplicitKernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, 1, 1 >
const arrayView1d< globalIndex const > m_dofNumber
 The global degree of freedom number.
 
const globalIndex m_dofRankOffset
 The global rank offset.
 
const CRSMatrixView< real64, globalIndex const > m_matrix
 The global Jacobian matrix.
 
const arrayView1d< real64m_rhs
 The global residaul vector.
 
FE_TYPE::template MeshData< SUBREGION_TYPE > m_meshData
 Data structure containing mesh data used to setup the finite element.
 
const real64 m_dt
 time increment
 
- Protected Attributes inherited from geos::finiteElement::KernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, NUM_DOF_PER_TEST_SP, NUM_DOF_PER_TRIAL_SP >
const traits::ViewTypeConst< typename SUBREGION_TYPE::NodeMapType::base_type > m_elemsToNodes
 The element to nodes map.
 
const arrayView1d< integer const > m_elemGhostRank
 The element ghost rank array.
 
const CONSTITUTIVE_TYPE::KernelWrapper m_constitutiveUpdate
 
const FE_TYPE & m_finiteElementSpace
 

Additional Inherited Members

- Static Public Member Functions inherited from geos::finiteElement::KernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, NUM_DOF_PER_TEST_SP, NUM_DOF_PER_TRIAL_SP >
template<typename POLICY , typename KERNEL_TYPE >
static real64 kernelLaunch (localIndex const numElems, KERNEL_TYPE const &kernelComponent)
 Kernel Launcher. More...
 
- Public Attributes inherited from geos::finiteElement::ImplicitKernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, 1, 1 >
const traits::ViewTypeConst< typename SUBREGION_TYPE::NodeMapType::base_type > m_elemsToNodes
 The element to nodes map.
 
const FE_TYPE & m_finiteElementSpace
 
- Static Public Attributes inherited from geos::finiteElement::ImplicitKernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, 1, 1 >
static constexpr int maxNumTestSupportPointsPerElem
 
static constexpr int maxNumTrialSupportPointsPerElem
 
static constexpr int numDofPerTestSupportPoint
 
static constexpr int numDofPerTrialSupportPoint
 
- Static Public Attributes inherited from geos::finiteElement::KernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, NUM_DOF_PER_TEST_SP, NUM_DOF_PER_TRIAL_SP >
static constexpr int maxNumTestSupportPointsPerElem = FE_TYPE::maxSupportPoints
 
static constexpr int maxNumTrialSupportPointsPerElem = FE_TYPE::maxSupportPoints
 
static constexpr int numDofPerTestSupportPoint = NUM_DOF_PER_TEST_SP
 
static constexpr int numDofPerTrialSupportPoint = NUM_DOF_PER_TRIAL_SP
 
static constexpr int numQuadraturePointsPerElem = FE_TYPE::numQuadraturePoints
 Compile time value for the number of quadrature points per element.
 

Detailed Description

template<typename SUBREGION_TYPE, typename CONSTITUTIVE_TYPE, typename FE_TYPE>
class geos::LaplaceFEMKernel< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >

Implements kernels for solving Laplace's equation.

Define the base interface for finite element kernels.

Template Parameters
SUBREGION_TYPEThe type of subregion that the kernel will act on.
CONSTITUTIVE_TYPEThe type of constitutive model present in the subregion.
NUM_TEST_SUPPORT_POINTS_PER_ELEMThe number of test space support points per element.
NUM_TRIAL_SUPPORT_POINTS_PER_ELEMThe number of trial space support points per element.
NUM_DOF_PER_TEST_SPThe number of DOF per test support point.
NUM_DOF_PER_TRIAL_SPThe number of DOF per trial support point.

General KernelBase Description

KernelBase defines an interface for implementing finite element kernels that will be callable by the family of kernel launch functions. Specific physics kernels may or may not derive from KernelBase, but must follow the same interface in order to be callable from the generic launching functions.

The template parameters of KernelBase should be duplicated as part of the interface, EXCEPT for NUM_DOF_PER_TEST_SP and NUM_DOF_PER_TRIAL_SP. These values should be set internally by the physics solver since each physics discretization will have a constant intrinsic value for these quantities. For example, when solving or the heat equation with scalar temperature as the primary variable at the support point, these will have a value of 1. In contrast, when solving a solid mechanics problem, with vector displacement as the primary variable at the support point, these will have a value of 3. Note that the interface provided by geos::finiteElement::RegionBasedKernelApplication will construct a kernel assuming only the first 4 template arguments.

Template Parameters
NUM_NODES_PER_ELEMThe number of nodes per element for the SUBREGION_TYPE.
UNUSEDAn unused parameter since we are assuming that the test and trial space have the same number of support points.

LaplaceFEMKernel Description

Implements the KernelBase interface functions required for solving Laplace's equation using on of the finite element kernel application functions such as geos::finiteElement::RegionBasedKernelApplication.

In this implementation, the template parameter NUM_NODES_PER_ELEM is used in place of both NUM_TEST_SUPPORT_POINTS_PER_ELEM and NUM_TRIAL_SUPPORT_POINTS_PER_ELEM, which are assumed to be equal. This results in the UNUSED template parameter as only the NUM_NODES_PER_ELEM is passed to the ImplicitKernelBase template to form the base class.

Additionally, the number of degrees of freedom per support point for both the test and trial spaces are specified as 1 when specifying the base class.

Definition at line 57 of file LaplaceFEMKernels.hpp.

Constructor & Destructor Documentation

◆ LaplaceFEMKernel()

template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
geos::LaplaceFEMKernel< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >::LaplaceFEMKernel ( NodeManager const &  nodeManager,
EdgeManager const &  edgeManager,
FaceManager const &  faceManager,
localIndex const  targetRegionIndex,
SUBREGION_TYPE const &  elementSubRegion,
FE_TYPE const &  finiteElementSpace,
CONSTITUTIVE_TYPE &  inputConstitutiveType,
arrayView1d< globalIndex const > const  inputDofNumber,
globalIndex const  rankOffset,
CRSMatrixView< real64, globalIndex const > const  inputMatrix,
arrayView1d< real64 > const  inputRhs,
real64 const  inputDt,
string const  fieldName 
)
inline

Constructor.

Constructor.

Parameters
nodeManagerReference to the NodeManager object.
edgeManagerReference to the EdgeManager object.
faceManagerReference to the FaceManager object.
targetRegionIndexIndex of the region the subregion belongs to.
inputDofNumberThe dof number for the primary field.
rankOffsetdof index offset of current rank
inputMatrixReference to the Jacobian matrix.
inputRhsReference to the RHS vector.
inputDtThe timestep for the physics update.

Constructor.

Parameters
elementSubRegionReference to the SUBREGION_TYPE(class template parameter) object.
inputConstitutiveTypeThe constitutive object.
finiteElementSpacePlaceholder for the finite element space object, which currently doesn't do much.
fieldNameThe name of the primary field (i.e. Temperature, Pressure, etc.)

Definition at line 90 of file LaplaceFEMKernels.hpp.

Member Function Documentation

◆ complete()

template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
GEOS_HOST_DEVICE real64 geos::LaplaceFEMKernel< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >::complete ( localIndex const  k,
StackVariables stack 
) const
inline

Performs the complete phase for the kernel.

Template Parameters
STACK_VARIABLE_TYPEThe type of StackVariable that holds the stack variables. This is most likely a defined in a type that derives from KernelBase.
Parameters
kThe element index.
stackThe StackVariable object that hold the stack variables.
Returns
The maximum contribution to the residual.

KernelBase::complete() Description

The operations typically found in complete are the mapping of the local Jacobian and Residual into the global Jacobian and Residual.

Form element residual from the fully formed element Jacobian dotted with the primary field and map the element local Jacobian/Residual to the global matrix/vector.

Definition at line 219 of file LaplaceFEMKernels.hpp.

◆ quadraturePointKernel()

template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
GEOS_HOST_DEVICE void geos::LaplaceFEMKernel< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >::quadraturePointKernel ( localIndex const  k,
localIndex const  q,
StackVariables stack 
) const
inline

Performs a state update at a quadrature point.

Template Parameters
STACK_VARIABLE_TYPEThe type of StackVariable that holds the stack variables. This is most likely a defined in a type that derives from KernelBase.
Parameters
kThe element index.
qThe quadrature point index.
stackThe StackVariable object that hold the stack variables.

KernelBase::quadraturePointKernel() Description

The operations found here are the mapping from the support points to the quadrature point, calculation of gradients, etc. From this data the state of the constitutive model is updated if required by the physics package.

Definition at line 194 of file LaplaceFEMKernels.hpp.

◆ setup()

template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
GEOS_HOST_DEVICE void geos::LaplaceFEMKernel< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >::setup ( localIndex const  k,
StackVariables stack 
) const
inline

Copy global values from primary field to a local stack array.

Performs the setup phase for the kernel.

Template Parameters
STACK_VARIABLE_TYPEThe type of StackVariable that holds the stack variables. This is most likely a defined in a type that derives from KernelBase.
Parameters
kThe element index.
stackThe StackVariable object that hold the stack variables.

KernelBase::setup() Description

The operations typically found in setup are thing such as the collection of global data into local stack storage.

ImplicitKernelBase::setup() Description

In this implementation, the element local Row and Column DOF stack arrays are filled for when we fill the global matrix and rhs.

Note
This seems like a waste of register space. We should do this in complete() unless we actually need these dof somewhere else in the kernel.

For the LaplaceFEMKernel implementation, global values from the primaryField, and degree of freedom numbers are placed into element local stack storage.

Definition at line 163 of file LaplaceFEMKernels.hpp.


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