GEOS
ImplicitKernelBase.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 #include "KernelBase.hpp"
17 
22 #ifndef GEOS_FINITEELEMENT_IMPLICITKERNELBASE_HPP_
23 #define GEOS_FINITEELEMENT_IMPLICITKERNELBASE_HPP_
24 
25 
26 
27 namespace geos
28 {
29 
30 namespace finiteElement
31 {
32 
33 //*****************************************************************************
34 //*****************************************************************************
35 //*****************************************************************************
46 template< typename SUBREGION_TYPE,
47  typename CONSTITUTIVE_TYPE,
48  typename FE_TYPE,
49  int NUM_DOF_PER_TEST_SP,
50  int NUM_DOF_PER_TRIAL_SP >
51 class ImplicitKernelBase : public KernelBase< SUBREGION_TYPE,
52  CONSTITUTIVE_TYPE,
53  FE_TYPE,
54  NUM_DOF_PER_TEST_SP,
55  NUM_DOF_PER_TRIAL_SP >
56 {
57 public:
59  using Base = KernelBase< SUBREGION_TYPE,
60  CONSTITUTIVE_TYPE,
61  FE_TYPE,
62  NUM_DOF_PER_TEST_SP,
63  NUM_DOF_PER_TRIAL_SP >;
64 
71 
72 
86  ImplicitKernelBase( NodeManager const & nodeManager,
87  EdgeManager const & edgeManager,
88  FaceManager const & faceManager,
89  localIndex const targetRegionIndex,
90  SUBREGION_TYPE const & elementSubRegion,
91  FE_TYPE const & finiteElementSpace,
92  CONSTITUTIVE_TYPE & inputConstitutiveType,
93  arrayView1d< globalIndex const > const & inputDofNumber,
94  globalIndex const rankOffset,
96  arrayView1d< real64 > const & inputRhs,
97  real64 const inputDt ):
98  Base( elementSubRegion,
99  finiteElementSpace,
100  inputConstitutiveType ),
101  m_dofNumber( inputDofNumber ),
102  m_dofRankOffset( rankOffset ),
103  m_matrix( inputMatrix ),
104  m_rhs( inputRhs ),
105  m_dt( inputDt )
106  {
107  FiniteElementBase::initialize< FE_TYPE >( nodeManager,
108  edgeManager,
109  faceManager,
110  elementSubRegion,
111  m_meshData );
112  GEOS_UNUSED_VAR( targetRegionIndex );
113  }
114 
115 
116  //***************************************************************************
121  {
124 
127 
133  localRowDofIndex{ 0 },
134  localColDofIndex{ 0 },
135  localResidual{ 0.0 },
136  localJacobian{ {0.0} }
137  {
138  for( int ii = 0; ii < maxNumRows; ++ii )
139  {
140  for( int jj = 0; jj < maxNumCols; ++jj )
141  {
142  localJacobian[ii][jj] = 0.0;
143  }
144  }
145  }
146 
149 
152 
154  globalIndex localRowDofIndex[maxNumRows];
155 
157  globalIndex localColDofIndex[maxNumCols];
158 
160  real64 localResidual[maxNumRows];
161 
163  real64 localJacobian[maxNumRows][maxNumCols];
164 
166  typename FE_TYPE::StackVariables feStack;
167  };
168  //***************************************************************************
169 
182  inline
183  void setup( localIndex const k,
184  StackVariables & stack ) const
185  {
186  m_finiteElementSpace.template setup< FE_TYPE >( k, m_meshData, stack.feStack );
187  localIndex numTestSupportPoints = m_finiteElementSpace.template numSupportPoints< FE_TYPE >( stack.feStack );
188  localIndex numTrialSupportPoints = numTestSupportPoints;
189  stack.numRows = numTestSupportPoints * numDofPerTestSupportPoint;
190  stack.numCols = numTrialSupportPoints * numDofPerTrialSupportPoint;
191  for( localIndex a=0; a<numTestSupportPoints; ++a )
192  {
193  localIndex const localNodeIndex = m_elemsToNodes[k][a];
194  for( int i=0; i<numDofPerTestSupportPoint; ++i )
195  {
196  stack.localRowDofIndex[a*numDofPerTestSupportPoint+i] = m_dofNumber[localNodeIndex]+i;
197  }
198  }
199 
200  for( localIndex a=0; a<numTrialSupportPoints; ++a )
201  {
202  localIndex const localNodeIndex = m_elemsToNodes[k][a];
203  for( int i=0; i<numDofPerTrialSupportPoint; ++i )
204  {
205  stack.localColDofIndex[a*numDofPerTrialSupportPoint+i] = m_dofNumber[localNodeIndex]+i;
206  }
207  }
208  }
209 
210 
211 protected:
214 
217 
220 
223 
225  typename FE_TYPE::template MeshData< SUBREGION_TYPE > m_meshData;
226 
228  real64 const m_dt;
229 
230 };
231 
232 }
233 }
234 
235 
236 
237 #endif /* GEOS_FINITEELEMENT_IMPLICITKERNELBASE_HPP_ */
#define GEOS_HOST_DEVICE
Marks a host-device function.
Definition: GeosxMacros.hpp:49
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:84
This class provides an interface to ObjectManagerBase in order to manage edge data.
Definition: EdgeManager.hpp:43
The FaceManager class provides an interface to ObjectManagerBase in order to manage face data.
Definition: FaceManager.hpp:44
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
Definition: NodeManager.hpp:46
Define the base interface for implicit finite element kernels.
globalIndex const m_dofRankOffset
The global rank offset.
CRSMatrixView< real64, globalIndex const > const m_matrix
The global Jacobian matrix.
FE_TYPE::template MeshData< SUBREGION_TYPE > m_meshData
Data structure containing mesh data used to setup the finite element.
static constexpr int numDofPerTestSupportPoint
Definition: KernelBase.hpp:99
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.
static constexpr int numDofPerTrialSupportPoint
Definition: KernelBase.hpp:103
GEOS_HOST_DEVICE void setup(localIndex const k, StackVariables &stack) const
Performs the setup phase for the kernel.
arrayView1d< globalIndex const > const m_dofNumber
The global degree of freedom number.
static constexpr int maxNumTrialSupportPointsPerElem
Definition: KernelBase.hpp:95
static constexpr int maxNumTestSupportPointsPerElem
Definition: KernelBase.hpp:91
arrayView1d< real64 > const m_rhs
The global residaul vector.
Define the base interface for finite element kernels.
Definition: KernelBase.hpp:87
traits::ViewTypeConst< typename SUBREGION_TYPE::NodeMapType::base_type > const m_elemsToNodes
The element to nodes map.
Definition: KernelBase.hpp:257
static constexpr int numDofPerTestSupportPoint
Definition: KernelBase.hpp:99
static constexpr int numDofPerTrialSupportPoint
Definition: KernelBase.hpp:103
static constexpr int maxNumTrialSupportPointsPerElem
Definition: KernelBase.hpp:95
static constexpr int maxNumTestSupportPointsPerElem
Definition: KernelBase.hpp:91
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:180
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:310
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:88
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
localIndex numCols
The actual number of columns in the element local jacobian matrix (<= maxNumCols).
real64 localJacobian[maxNumRows][maxNumCols]
C-array storage for the element local Jacobian matrix.
globalIndex localRowDofIndex[maxNumRows]
C-array storage for the element local row degrees of freedom.
globalIndex localColDofIndex[maxNumCols]
C-array storage for the element local column degrees of freedom.
FE_TYPE::StackVariables feStack
Stack variables needed for the underlying FEM type.
static constexpr int maxNumRows
The number of rows in the pre-allocated element local jacobian matrix (upper bound for numRows).
real64 localResidual[maxNumRows]
C-array storage for the element local residual vector.
static constexpr int maxNumCols
The number of columns in the pre-allocated element local jacobian matrix (upper bound for numCols).
localIndex numRows
The actual number of rows in the element local jacobian matrix (<= maxNumRows).
Kernel variables allocated on the stack.
Definition: KernelBase.hpp:137