GEOSX
LaplaceFEMKernels.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 TotalEnergies
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 
20 #ifndef GEOS_PHYSICSSOLVERS_SIMPLEPDE_LAPLACEFEMKERNELS_HPP_
21 #define GEOS_PHYSICSSOLVERS_SIMPLEPDE_LAPLACEFEMKERNELS_HPP_
22 
23 #define GEOSX_DISPATCH_VEM
24 
26 
27 namespace geos
28 {
29 
30 //*****************************************************************************
54 template< typename SUBREGION_TYPE,
55  typename CONSTITUTIVE_TYPE,
56  typename FE_TYPE >
58  public finiteElement::ImplicitKernelBase< SUBREGION_TYPE,
59  CONSTITUTIVE_TYPE,
60  FE_TYPE,
61  1,
62  1 >
63 {
64 public:
66  using Base = finiteElement::ImplicitKernelBase< SUBREGION_TYPE,
67  CONSTITUTIVE_TYPE,
68  FE_TYPE,
69  1,
70  1 >;
71 
75  using Base::m_dofNumber;
77  using Base::m_matrix;
78  using Base::m_rhs;
81  using Base::m_meshData;
82  using Base::m_dt;
83 
90  LaplaceFEMKernel( NodeManager const & nodeManager,
91  EdgeManager const & edgeManager,
92  FaceManager const & faceManager,
93  localIndex const targetRegionIndex,
94  SUBREGION_TYPE const & elementSubRegion,
95  FE_TYPE const & finiteElementSpace,
96  CONSTITUTIVE_TYPE & inputConstitutiveType,
97  arrayView1d< globalIndex const > const inputDofNumber,
98  globalIndex const rankOffset,
100  arrayView1d< real64 > const inputRhs,
101  real64 const inputDt,
102  string const fieldName ):
103  Base( nodeManager,
104  edgeManager,
105  faceManager,
106  targetRegionIndex,
107  elementSubRegion,
108  finiteElementSpace,
109  inputConstitutiveType,
110  inputDofNumber,
111  rankOffset,
112  inputMatrix,
113  inputRhs,
114  inputDt ),
115  m_X( nodeManager.referencePosition() ),
116  m_primaryField( nodeManager.template getReference< array1d< real64 > >( fieldName ))
117  {}
118 
119  //***************************************************************************
127  {
128 public:
129 
135  Base::StackVariables(),
136  xLocal(),
137  primaryField_local{ 0.0 }
138  {}
139 
140 #if !defined(CALC_FEM_SHAPE_IN_KERNEL)
142  int xLocal;
143 #else
146 #endif
147 
150  };
151 
152 
162  inline
163  void setup( localIndex const k,
164  StackVariables & stack ) const
165  {
166  m_finiteElementSpace.template setup< FE_TYPE >( k, m_meshData, stack.feStack );
167  stack.numRows = m_finiteElementSpace.template numSupportPoints< FE_TYPE >( stack.feStack );
168  stack.numCols = stack.numRows;
169  for( localIndex a = 0; a < stack.numRows; ++a )
170  {
171  localIndex const localNodeIndex = m_elemsToNodes( k, a );
172 
173 #if defined(CALC_FEM_SHAPE_IN_KERNEL)
174  for( int i=0; i<3; ++i )
175  {
176  stack.xLocal[ a ][ i ] = m_X[ localNodeIndex ][ i ];
177  }
178 #endif
179 
180  stack.primaryField_local[ a ] = m_primaryField[ localNodeIndex ];
181  stack.localRowDofIndex[a] = m_dofNumber[localNodeIndex];
182  stack.localColDofIndex[a] = m_dofNumber[localNodeIndex];
183  }
184  m_finiteElementSpace.template
185  addGradGradStabilizationMatrix< FE_TYPE, numDofPerTrialSupportPoint >( stack.feStack,
186  stack.localJacobian );
187  }
188 
193  inline
195  localIndex const q,
196  StackVariables & stack ) const
197  {
199  real64 const detJ = m_finiteElementSpace.template getGradN< FE_TYPE >( k, q, stack.xLocal,
200  stack.feStack, dNdX );
201  for( localIndex a = 0; a < stack.numRows; ++a )
202  {
203  for( localIndex b = 0; b < stack.numCols; ++b )
204  {
205  stack.localJacobian[ a ][ b ] += LvArray::tensorOps::AiBi< 3 >( dNdX[a], dNdX[b] ) * detJ;
206  }
207  }
208  }
209 
218  inline
220  StackVariables & stack ) const
221  {
222  GEOS_UNUSED_VAR( k );
223  real64 maxForce = 0;
224 
225  for( localIndex a = 0; a < stack.numRows; ++a )
226  {
227  for( localIndex b = 0; b < stack.numCols; ++b )
228  {
229  stack.localResidual[ a ] += stack.localJacobian[ a ][ b ] * stack.primaryField_local[ b ];
230  }
231  }
232 
233  for( int a = 0; a < stack.numRows; ++a )
234  {
235  localIndex const dof = LvArray::integerConversion< localIndex >( stack.localRowDofIndex[ a ] - m_dofRankOffset );
236  if( dof < 0 || dof >= m_matrix.numRows() ) continue;
237  m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
238  stack.localColDofIndex,
239  stack.localJacobian[ a ],
240  stack.numCols );
241 
242  RAJA::atomicAdd< parallelDeviceAtomic >( &m_rhs[ dof ], stack.localResidual[ a ] );
243  maxForce = fmax( maxForce, fabs( stack.localResidual[ a ] ) );
244  }
245 
246  return maxForce;
247  }
248 
249 protected:
252 
255 
256 };
257 
261  globalIndex const,
263  arrayView1d< real64 > const,
264  real64 const,
265  string const >;
266 
267 } // namespace geos
268 
270 
271 #endif // GEOS_PHYSICSSOLVERS_SIMPLEPDE_LAPLACEFEMKERNELS_HPP_
#define GEOS_HOST_DEVICE
Marks a host-device function.
Definition: GeosxMacros.hpp:48
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:83
This class provides an interface to ObjectManagerBase in order to manage edge data.
Definition: EdgeManager.hpp:42
The FaceManager class provides an interface to ObjectManagerBase in order to manage face data.
Definition: FaceManager.hpp:43
Implements kernels for solving Laplace's equation.
GEOS_HOST_DEVICE real64 complete(localIndex const k, StackVariables &stack) const
Performs the complete phase for the kernel.
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.
GEOS_HOST_DEVICE void quadraturePointKernel(localIndex const k, localIndex const q, StackVariables &stack) const
Performs a state update at a quadrature point.
arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const m_X
The array containing the nodal position array.
GEOS_HOST_DEVICE void setup(localIndex const k, StackVariables &stack) const
Copy global values from primary field to a local stack array.
arrayView1d< real64 const > const m_primaryField
The global primary field array.
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
Definition: NodeManager.hpp:45
Define the base interface for implicit finite element kernels.
traits::ViewTypeConst< typename SUBREGION_TYPE::NodeMapType::base_type > const m_elemsToNodes
The element to nodes map.
Definition: KernelBase.hpp:256
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.
arrayView1d< globalIndex const > const m_dofNumber
The global degree of freedom number.
Used to forward arguments to a class that implements the KernelBase interface.
Definition: KernelBase.hpp:281
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:220
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:350
double real64
64-bit floating point type.
Definition: DataTypes.hpp:139
GEOSX_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:128
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Definition: DataTypes.hpp:236
GEOSX_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:125
Array< T, 1 > array1d
Alias for 1D array.
Definition: DataTypes.hpp:216
Kernel variables allocated on the stack.
GEOS_HOST_DEVICE StackVariables()
Constructor.
real64 primaryField_local[maxNumTestSupportPointsPerElem]
C-array storage for the element local primary field variable.
Kernel variables allocated on the stack.
Definition: KernelBase.hpp:136