GEOSX
ImplicitKernelBase.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2018-2020 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2020 The Board of Trustees of the Leland Stanford Junior University
7  * Copyright (c) 2018-2020 Total, S.A
8  * Copyright (c) 2019- GEOSX Contributors
9  * All rights reserved
10  *
11  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
12  * ------------------------------------------------------------------------------------------------------------
13  */
14 
15 #include "KernelBase.hpp"
16 
21 #ifndef GEOSX_FINITEELEMENT_IMPLICITKERNELBASE_HPP_
22 #define GEOSX_FINITEELEMENT_IMPLICITKERNELBASE_HPP_
23 
24 
25 
26 namespace geosx
27 {
28 
29 namespace finiteElement
30 {
31 
32 //*****************************************************************************
33 //*****************************************************************************
34 //*****************************************************************************
45 template< typename SUBREGION_TYPE,
46  typename CONSTITUTIVE_TYPE,
47  typename FE_TYPE,
48  int NUM_DOF_PER_TEST_SP,
49  int NUM_DOF_PER_TRIAL_SP >
50 class ImplicitKernelBase : public KernelBase< SUBREGION_TYPE,
51  CONSTITUTIVE_TYPE,
52  FE_TYPE,
53  NUM_DOF_PER_TEST_SP,
54  NUM_DOF_PER_TRIAL_SP >
55 {
56 public:
58  using Base = KernelBase< SUBREGION_TYPE,
59  CONSTITUTIVE_TYPE,
60  FE_TYPE,
61  NUM_DOF_PER_TEST_SP,
62  NUM_DOF_PER_TRIAL_SP >;
63 
69 
70 
82  ImplicitKernelBase( NodeManager const & nodeManager,
83  EdgeManager const & edgeManager,
84  FaceManager const & faceManager,
85  SUBREGION_TYPE const & elementSubRegion,
86  FE_TYPE const & finiteElementSpace,
87  CONSTITUTIVE_TYPE * const inputConstitutiveType,
88  arrayView1d< globalIndex const > const & inputDofNumber,
89  globalIndex const rankOffset,
91  arrayView1d< real64 > const & inputRhs ):
92  Base( elementSubRegion,
93  finiteElementSpace,
94  inputConstitutiveType ),
95  m_dofNumber( inputDofNumber ),
96  m_dofRankOffset( rankOffset ),
97  m_matrix( inputMatrix ),
98  m_rhs( inputRhs )
99  {
100  GEOSX_UNUSED_VAR( nodeManager );
101  GEOSX_UNUSED_VAR( edgeManager );
102  GEOSX_UNUSED_VAR( faceManager );
103  }
104 
105 
106  //***************************************************************************
111  {
114 
117 
123  localRowDofIndex{ 0 },
124  localColDofIndex{ 0 },
125  localResidual{ 0.0 },
126  localJacobian{ {0.0} }
127  {}
128 
131 
134 
137 
140  };
141  //***************************************************************************
142 
156  void setup( localIndex const k,
157  StackVariables & stack ) const
158  {
159  for( localIndex a=0; a<numTestSupportPointsPerElem; ++a )
160  {
161  localIndex const localNodeIndex = m_elemsToNodes[k][a];
162  for( int i=0; i<numDofPerTestSupportPoint; ++i )
163  {
164  stack.localRowDofIndex[a*numDofPerTestSupportPoint+i] = m_dofNumber[localNodeIndex]+i;
165  }
166  }
167 
168  for( localIndex a=0; a<numTrialSupportPointsPerElem; ++a )
169  {
170  localIndex const localNodeIndex = m_elemsToNodes[k][a];
171  for( int i=0; i<numDofPerTrialSupportPoint; ++i )
172  {
173  stack.localColDofIndex[a*numDofPerTrialSupportPoint+i] = m_dofNumber[localNodeIndex]+i;
174  }
175  }
176  }
177 
178 
179 protected:
182 
185 
188 
191 
192 };
193 
194 }
195 }
196 
197 
198 
199 #endif /* GEOSX_FINITEELEMENT_IMPLICITKERNELBASE_HPP_ */
#define GEOSX_HOST_DEVICE
Marks a host-device function.
Definition: GeosxMacros.hpp:48
Define the base interface for implicit finite element kernels.
real64 localResidual[numRows]
C-array storage for the element local residual vector.
Kernel variables allocated on the stack.
Definition: KernelBase.hpp:182
long long int globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:128
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data...
Definition: NodeManager.hpp:47
ImplicitKernelBase(NodeManager const &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, SUBREGION_TYPE const &elementSubRegion, FE_TYPE const &finiteElementSpace, CONSTITUTIVE_TYPE *const inputConstitutiveType, arrayView1d< globalIndex const > const &inputDofNumber, globalIndex const rankOffset, CRSMatrixView< real64, globalIndex const > const &inputMatrix, arrayView1d< real64 > const &inputRhs)
Constructor.
GEOSX_HOST_DEVICE GEOSX_FORCE_INLINE void setup(localIndex const k, StackVariables &stack) const
Performs the setup phase for the kernel.
traits::ViewTypeConst< typename SUBREGION_TYPE::NodeMapType::base_type > const m_elemsToNodes
The element to nodes map.
Definition: KernelBase.hpp:302
This class serves to provide a "view" of a multidimensional array.
Definition: ArrayView.hpp:67
static constexpr int numRows
The number of rows in the element local jacobian matrix.
arrayView1d< globalIndex const > const m_dofNumber
The global degree of freedom number.
globalIndex const m_dofRankOffset
The global rank offset.
double real64
64-bit floating point type.
Definition: DataTypes.hpp:136
Define the base interface for finite element kernels.
Definition: KernelBase.hpp:132
static constexpr int numTrialSupportPointsPerElem
Definition: KernelBase.hpp:141
globalIndex localColDofIndex[numCols]
C-array storage for the element local column degrees of freedom.
static constexpr int numTestSupportPointsPerElem
Definition: KernelBase.hpp:137
static constexpr int numDofPerTestSupportPoint
Definition: KernelBase.hpp:145
This class provides an interface to ObjectManagerBase in order to manage edge data.
Definition: EdgeManager.hpp:42
#define GEOSX_FORCE_INLINE
Marks a function or lambda for inlining.
Definition: GeosxMacros.hpp:50
static constexpr int numCols
The number of columns in the element local jacobian matrix.
#define GEOSX_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:78
static constexpr int numDofPerTrialSupportPoint
Definition: KernelBase.hpp:149
std::ptrdiff_t localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:125
CRSMatrixView< real64, globalIndex const > const m_matrix
The global Jacobian matrix.
The FaceManager class provides an interface to ObjectManagerBase in order to manage face data...
Definition: FaceManager.hpp:40
globalIndex localRowDofIndex[numRows]
C-array storage for the element local row degrees of freedom.
arrayView1d< real64 > const m_rhs
The global residaul vector.
real64 localJacobian[numRows][numCols]
C-array storage for the element local Jacobian matrix.