https://mooseframework.inl.gov
Classes | Public Member Functions | Private Member Functions | Private Attributes | List of all members
Moose::Kokkos::System Class Reference

The Kokkos system class. More...

#include <KokkosSystem.h>

Inheritance diagram for Moose::Kokkos::System:
[legend]

Classes

struct  Sparsity
 CSR format sparsity data. More...
 

Public Member Functions

 System (SystemBase &system)
 Constructor. More...
 
void sync (const MemcpyKind dir)
 Synchronize the active tagged vectors and matrices between host and device. More...
 
void sync (const std::set< TagID > &tags, const MemcpyKind dir)
 Synchronize the specified tagged vectors between host and device. More...
 
void reinit ()
 Allocate the quadrature point solution vectors for active variable and tags and cache quadrature point values. More...
 
void setActiveVariables (const std::set< MooseVariableFieldBase *> &vars)
 Set the active variables. More...
 
void setActiveVariableTags (const std::set< TagID > &tags)
 Set the active solution tags. More...
 
void setActiveResidualTags (const std::set< TagID > &tags)
 Set the active residual tags. More...
 
void setActiveMatrixTags (const std::set< TagID > &tags)
 Set the active matrix tags. More...
 
void clearActiveVariables ()
 Clear the cached active variables. More...
 
void clearActiveVariableTags ()
 Clear the cached active solution tags. More...
 
void clearActiveResidualTags ()
 Clear the cached active residual tags. More...
 
void clearActiveMatrixTags ()
 Clear the cached active matrix tags. More...
 
const auto & getSystem () const
 Get the MOOSE system. More...
 
const auto & getDofMap () const
 Get the libMesh DOF map. More...
 
const auto & getComm () const
 Get the libMesh communicator. More...
 
const auto & getLocalCommList () const
 Get the list of local DOF indices to communicate. More...
 
const auto & getGhostCommList () const
 Get the list of ghost DOF indices to communicate. More...
 
const auto & getSparsity () const
 Get the sparisty pattern data. More...
 
KOKKOS_FUNCTION const auto & getCoupling (unsigned int var) const
 Get the list of off-diagonal coupled variable numbers of a variable. More...
 
KOKKOS_FUNCTION bool isVariableActive (unsigned int var, ContiguousSubdomainID subdomain) const
 Check whether a variable is active on a subdomain. More...
 
KOKKOS_FUNCTION bool isResidualTagActive (TagID tag) const
 Check whether a residual tag is active. More...
 
KOKKOS_FUNCTION bool isMatrixTagActive (TagID tag) const
 Check whether a matrix tag is active. More...
 
KOKKOS_FUNCTION bool hasNodalBC (dof_id_type dof) const
 Check whether a local DOF index is covered by a nodal BC. More...
 
KOKKOS_FUNCTION bool hasNodalBCResidualTag (dof_id_type dof, TagID tag) const
 Check whether a local DOF index is covered by a nodal BC for an extra residual tag. More...
 
KOKKOS_FUNCTION bool hasNodalBCMatrixTag (dof_id_type dof, TagID tag) const
 Check whether a local DOF index is associated with a nodal BC for an extra matrix tag. More...
 
KOKKOS_FUNCTION unsigned int getFETypeID (unsigned int var) const
 Get the FE type ID of a variable. More...
 
KOKKOS_FUNCTION dof_id_type getNumLocalDofs () const
 Get the number of local DOFs. More...
 
KOKKOS_FUNCTION dof_id_type getNumGhostDofs () const
 Get the number of ghost DOFs. More...
 
KOKKOS_FUNCTION unsigned int getMaxDofsPerElem (unsigned int var) const
 Get the maximum number of DOFs per element of a variable This number is local to each process. More...
 
KOKKOS_FUNCTION dof_id_type getElemLocalDofIndex (ContiguousElementID elem, unsigned int i, unsigned int var) const
 Get the local DOF index of a variable for an element. More...
 
KOKKOS_FUNCTION dof_id_type getNodeLocalDofIndex (ContiguousNodeID node, unsigned int var) const
 Get the local DOF index of a variable for a node. More...
 
KOKKOS_FUNCTION dof_id_type getElemGlobalDofIndex (ContiguousElementID elem, unsigned int i, unsigned int var) const
 Get the global DOF index of a variable for an element. More...
 
KOKKOS_FUNCTION dof_id_type getNodeGlobalDofIndex (ContiguousNodeID node, unsigned int var) const
 Get the global DOF index of a variable for a node. More...
 
KOKKOS_FUNCTION dof_id_type localToGlobalDofIndex (dof_id_type dof) const
 Get the global DOF index of a local DOF index. More...
 
KOKKOS_FUNCTION bool isNodalDefined (ContiguousNodeID node, unsigned int var) const
 Get whether a variable is defined on a node. More...
 
KOKKOS_FUNCTION auto & getVector (TagID tag) const
 Get a tagged Kokkos vector. More...
 
KOKKOS_FUNCTION auto & getMatrix (TagID tag) const
 Get a tagged Kokkos matrix. More...
 
KOKKOS_FUNCTION RealgetVectorDofValue (dof_id_type dof, TagID tag) const
 Get the DOF value of a tagged vector. More...
 
KOKKOS_FUNCTION RealgetVectorQpValue (ElementInfo info, dof_id_type qp, unsigned int var, TagID tag) const
 Get the quadrature value of a variable from a tagged vector. More...
 
KOKKOS_FUNCTION Real3getVectorQpGrad (ElementInfo info, dof_id_type qp, unsigned int var, TagID tag) const
 Get the quadrature gradient of a variable from a tagged vector. More...
 
KOKKOS_FUNCTION Real getVectorQpValueFace (const ElementInfo info, const unsigned int side, const unsigned int qp, const unsigned int var, const TagID tag) const
 Get the face quadrature value of a variable from a tagged vector. More...
 
KOKKOS_FUNCTION Real3 getVectorQpGradFace (const ElementInfo info, const unsigned int side, const Real33 jacobian, const unsigned int qp, const unsigned int var, const TagID tag) const
 Get the face quadrature gradient of a variable from a tagged vector. More...
 
KOKKOS_FUNCTION RealgetMatrixValue (dof_id_type row, dof_id_type col, TagID tag) const
 Get an entry from a tagged matrix. More...
 
KOKKOS_FUNCTION void operator() (const ThreadID tid) const
 Kokkos function for caching variable values on element quadrature points. More...
 
KOKKOS_FUNCTION const MeshkokkosMesh () const
 Get the const reference of the Kokkos mesh. More...
 
KOKKOS_FUNCTION const AssemblykokkosAssembly () const
 Get the const reference of the Kokkos assembly. More...
 

Private Member Functions

void setupVariables ()
 Setup variable data. More...
 
void setupDofs ()
 Setup DOF data. More...
 
void setupSparsity ()
 Setup sparsity data. More...
 
void checkNodalBCs ()
 Check if the DOFs are covered by nodal BCs. More...
 
void getNodalBCDofs (const NodalBCBase *nbc, Array< bool > &dofs)
 Get the list of DOFs covered by a nodal BC. More...
 

Private Attributes

Thread _thread
 Kokkos thread object. More...
 
SystemBase_system
 Reference of the MOOSE system. More...
 
const MooseMesh_mesh
 Reference of the MOOSE mesh. More...
 
const libMesh::DofMap_dof_map
 Reference of the libMesh DOF map. More...
 
const Parallel::Communicator_comm
 Reference of the libMesh communicator. More...
 
const unsigned int _num_vars
 Number of variables. More...
 
const dof_id_type _num_local_dofs
 Number of local DOFs. More...
 
const dof_id_type _num_ghost_dofs
 Number of ghost DOFs. More...
 
Array< dof_id_type_local_to_global_dof_index
 Map from local DOF index to global DOF index. More...
 
Array< unsigned int_max_dofs_per_elem
 Maximum number of DOFs per element for each variable. More...
 
Array< unsigned int_var_fe_types
 FE type ID of each variable. More...
 
Array2D< bool > _var_subdomain_active
 Whether each variable is active on subdomains. More...
 
Array< Array< unsigned int > > _coupling
 Off-diagonal coupled variable numbers of each variable. More...
 
Array< unsigned int_active_variables
 List of active variable numbers. More...
 
Array< bool > _nbc_dof
 Flag whether each DOF is covered by a nodal BC. More...
 
Sparsity _sparsity
 Matrix sparsity pattern data. More...
 
Array< Vector_vectors
 Kokkos vectors and matrices on device. More...
 
Array< Matrix_matrices
 
Array< Array2D< Array< Real > > > _qp_solutions
 Cached elemental quadrature values and gradients. More...
 
Array< Array2D< Array< Real3 > > > _qp_solutions_grad
 
Array< Array2D< dof_id_type > > _local_elem_dof_index
 Local DOF index of each variable. More...
 
Array< Array< dof_id_type > > _local_node_dof_index
 
Array< TagID_active_variable_tags
 List of active tags. More...
 
Array< TagID_active_residual_tags
 
Array< TagID_active_matrix_tags
 
Array< bool > _residual_tag_active
 Flag whether each tag is active. More...
 
Array< bool > _matrix_tag_active
 
Array< Array< bool > > _nbc_residual_tag_dof
 Flag whether each DOF is covered by a nodal BC for extra residual tags. More...
 
Array< Array< bool > > _nbc_matrix_tag_dof
 
Array< Array< dof_id_type > > _local_comm_list
 List of DOFs to send and receive. More...
 
Array< Array< dof_id_type > > _ghost_comm_list
 

Detailed Description

The Kokkos system class.

Each nonlinear and auxiliary system in MOOSE has a corresponding Kokkos system.

Definition at line 35 of file KokkosSystem.h.

Constructor & Destructor Documentation

◆ System()

Moose::Kokkos::System::System ( SystemBase system)

Constructor.

Parameters
systemThe associated MOOSE system

Member Function Documentation

◆ checkNodalBCs()

void Moose::Kokkos::System::checkNodalBCs ( )
private

Check if the DOFs are covered by nodal BCs.

◆ clearActiveMatrixTags()

void Moose::Kokkos::System::clearActiveMatrixTags ( )
inline

Clear the cached active matrix tags.

Definition at line 100 of file KokkosSystem.h.

101  {
102  _active_matrix_tags.destroy();
103  _matrix_tag_active = false;
104  }
Array< TagID > _active_matrix_tags
Definition: KokkosSystem.h:498
Array< bool > _matrix_tag_active
Definition: KokkosSystem.h:505

◆ clearActiveResidualTags()

void Moose::Kokkos::System::clearActiveResidualTags ( )
inline

Clear the cached active residual tags.

Definition at line 92 of file KokkosSystem.h.

93  {
94  _active_residual_tags.destroy();
95  _residual_tag_active = false;
96  }
Array< bool > _residual_tag_active
Flag whether each tag is active.
Definition: KokkosSystem.h:504
Array< TagID > _active_residual_tags
Definition: KokkosSystem.h:497

◆ clearActiveVariables()

void Moose::Kokkos::System::clearActiveVariables ( )
inline

Clear the cached active variables.

Definition at line 84 of file KokkosSystem.h.

84 { _active_variables.destroy(); }
Array< unsigned int > _active_variables
List of active variable numbers.
Definition: KokkosSystem.h:491

◆ clearActiveVariableTags()

void Moose::Kokkos::System::clearActiveVariableTags ( )
inline

Clear the cached active solution tags.

Definition at line 88 of file KokkosSystem.h.

88 { _active_variable_tags.destroy(); }
Array< TagID > _active_variable_tags
List of active tags.
Definition: KokkosSystem.h:496

◆ getComm()

const auto& Moose::Kokkos::System::getComm ( ) const
inline

Get the libMesh communicator.

Returns
The libMesh communicator

Definition at line 119 of file KokkosSystem.h.

119 { return _comm; }
const Parallel::Communicator & _comm
Reference of the libMesh communicator.
Definition: KokkosSystem.h:432

◆ getCoupling()

KOKKOS_FUNCTION const auto& Moose::Kokkos::System::getCoupling ( unsigned int  var) const
inline

Get the list of off-diagonal coupled variable numbers of a variable.

Parameters
varThe variable number
Returns
The list of off-diagonal coupled variable numbers

Definition at line 140 of file KokkosSystem.h.

140 { return _coupling[var]; }
Array< Array< unsigned int > > _coupling
Off-diagonal coupled variable numbers of each variable.
Definition: KokkosSystem.h:486

◆ getDofMap()

const auto& Moose::Kokkos::System::getDofMap ( ) const
inline

Get the libMesh DOF map.

Returns
The libMesh DOF map

Definition at line 114 of file KokkosSystem.h.

114 { return _system.dofMap(); }
SystemBase & _system
Reference of the MOOSE system.
Definition: KokkosSystem.h:420
virtual libMesh::DofMap & dofMap()
Gets writeable reference to the dof map.
Definition: SystemBase.C:1163

◆ getElemGlobalDofIndex()

KOKKOS_FUNCTION dof_id_type Moose::Kokkos::System::getElemGlobalDofIndex ( ContiguousElementID  elem,
unsigned int  i,
unsigned int  var 
) const
inline

Get the global DOF index of a variable for an element.

Parameters
elemThe contiguous element ID
iThe element-local DOF index
varThe variable number
Returns
The global DOF index

Definition at line 248 of file KokkosSystem.h.

251  {
253  }
Array< dof_id_type > _local_to_global_dof_index
Map from local DOF index to global DOF index.
Definition: KokkosSystem.h:470
Array< Array2D< dof_id_type > > _local_elem_dof_index
Local DOF index of each variable.
Definition: KokkosSystem.h:464

◆ getElemLocalDofIndex()

KOKKOS_FUNCTION dof_id_type Moose::Kokkos::System::getElemLocalDofIndex ( ContiguousElementID  elem,
unsigned int  i,
unsigned int  var 
) const
inline

Get the local DOF index of a variable for an element.

Parameters
elemThe contiguous element ID
iThe element-local DOF index
varThe variable number
Returns
The local DOF index

Definition at line 225 of file KokkosSystem.h.

Referenced by Moose::Kokkos::ResidualObject::accumulateTaggedElementalMatrix(), Moose::Kokkos::ResidualObject::accumulateTaggedElementalResidual(), getVectorQpGradFace(), and getVectorQpValueFace().

228  {
229  return _local_elem_dof_index[var](elem, i);
230  }
Array< Array2D< dof_id_type > > _local_elem_dof_index
Local DOF index of each variable.
Definition: KokkosSystem.h:464

◆ getFETypeID()

KOKKOS_FUNCTION unsigned int Moose::Kokkos::System::getFETypeID ( unsigned int  var) const
inline

Get the FE type ID of a variable.

Parameters
varThe variable number
Returns
The FE type ID

Definition at line 197 of file KokkosSystem.h.

197 { return _var_fe_types[var]; }
Array< unsigned int > _var_fe_types
FE type ID of each variable.
Definition: KokkosSystem.h:478

◆ getGhostCommList()

const auto& Moose::Kokkos::System::getGhostCommList ( ) const
inline

Get the list of ghost DOF indices to communicate.

Returns
The list of ghost DOF indices to communicate

Definition at line 129 of file KokkosSystem.h.

129 { return _ghost_comm_list; }
Array< Array< dof_id_type > > _ghost_comm_list
Definition: KokkosSystem.h:523

◆ getLocalCommList()

const auto& Moose::Kokkos::System::getLocalCommList ( ) const
inline

Get the list of local DOF indices to communicate.

Returns
The list of local DOF indices to communicate

Definition at line 124 of file KokkosSystem.h.

124 { return _local_comm_list; }
Array< Array< dof_id_type > > _local_comm_list
List of DOFs to send and receive.
Definition: KokkosSystem.h:522

◆ getMatrix()

KOKKOS_FUNCTION auto& Moose::Kokkos::System::getMatrix ( TagID  tag) const
inline

Get a tagged Kokkos matrix.

Parameters
tagThe matrix tag
Returns
The Kokkos matrix

Definition at line 294 of file KokkosSystem.h.

294 { return _matrices[tag]; }
Array< Matrix > _matrices
Definition: KokkosSystem.h:451

◆ getMatrixValue()

KOKKOS_FUNCTION Real& Moose::Kokkos::System::getMatrixValue ( dof_id_type  row,
dof_id_type  col,
TagID  tag 
) const
inline

Get an entry from a tagged matrix.

Parameters
rowThe local row index
colThe global column index
tagThe matrix tag
Returns
The entry from the tagged matrix

Definition at line 368 of file KokkosSystem.h.

369  {
370  return _matrices[tag](row, col);
371  }
Array< Matrix > _matrices
Definition: KokkosSystem.h:451

◆ getMaxDofsPerElem()

KOKKOS_FUNCTION unsigned int Moose::Kokkos::System::getMaxDofsPerElem ( unsigned int  var) const
inline

Get the maximum number of DOFs per element of a variable This number is local to each process.

Parameters
varThe variable number
Returns
The maximum number of DOFs per element

Definition at line 214 of file KokkosSystem.h.

215  {
216  return _max_dofs_per_elem[var];
217  }
Array< unsigned int > _max_dofs_per_elem
Maximum number of DOFs per element for each variable.
Definition: KokkosSystem.h:474

◆ getNodalBCDofs()

void Moose::Kokkos::System::getNodalBCDofs ( const NodalBCBase nbc,
Array< bool > &  dofs 
)
private

Get the list of DOFs covered by a nodal BC.

Parameters
nbcThe Kokkos nodal BC object
dofsThe flag whether each DOF is covered by the nodal BC

◆ getNodeGlobalDofIndex()

KOKKOS_FUNCTION dof_id_type Moose::Kokkos::System::getNodeGlobalDofIndex ( ContiguousNodeID  node,
unsigned int  var 
) const
inline

Get the global DOF index of a variable for a node.

Parameters
nodeThe contiguous node ID
varThe variable number
Returns
The global DOF index

Definition at line 260 of file KokkosSystem.h.

261  {
263  }
Array< dof_id_type > _local_to_global_dof_index
Map from local DOF index to global DOF index.
Definition: KokkosSystem.h:470
Array< Array< dof_id_type > > _local_node_dof_index
Definition: KokkosSystem.h:465

◆ getNodeLocalDofIndex()

KOKKOS_FUNCTION dof_id_type Moose::Kokkos::System::getNodeLocalDofIndex ( ContiguousNodeID  node,
unsigned int  var 
) const
inline

Get the local DOF index of a variable for a node.

Parameters
nodeThe contiguous node ID
varThe variable number
Returns
The local DOF index

Definition at line 237 of file KokkosSystem.h.

Referenced by Moose::Kokkos::ResidualObject::accumulateTaggedNodalMatrix(), and Moose::Kokkos::ResidualObject::accumulateTaggedNodalResidual().

238  {
239  return _local_node_dof_index[var][node];
240  }
Array< Array< dof_id_type > > _local_node_dof_index
Definition: KokkosSystem.h:465

◆ getNumGhostDofs()

KOKKOS_FUNCTION dof_id_type Moose::Kokkos::System::getNumGhostDofs ( ) const
inline

Get the number of ghost DOFs.

Returns
The number of ghost DOFs

Definition at line 207 of file KokkosSystem.h.

207 { return _num_ghost_dofs; }
const dof_id_type _num_ghost_dofs
Number of ghost DOFs.
Definition: KokkosSystem.h:444

◆ getNumLocalDofs()

KOKKOS_FUNCTION dof_id_type Moose::Kokkos::System::getNumLocalDofs ( ) const
inline

Get the number of local DOFs.

Returns
The number of local DOFs

Definition at line 202 of file KokkosSystem.h.

202 { return _num_local_dofs; }
const dof_id_type _num_local_dofs
Number of local DOFs.
Definition: KokkosSystem.h:440

◆ getSparsity()

const auto& Moose::Kokkos::System::getSparsity ( ) const
inline

Get the sparisty pattern data.

Returns
The sparisty pattern data

Definition at line 134 of file KokkosSystem.h.

134 { return _sparsity; }
Sparsity _sparsity
Matrix sparsity pattern data.
Definition: KokkosSystem.h:529

◆ getSystem()

const auto& Moose::Kokkos::System::getSystem ( ) const
inline

Get the MOOSE system.

Returns
The MOOSE system

Definition at line 109 of file KokkosSystem.h.

109 { return _system; }
SystemBase & _system
Reference of the MOOSE system.
Definition: KokkosSystem.h:420

◆ getVector()

KOKKOS_FUNCTION auto& Moose::Kokkos::System::getVector ( TagID  tag) const
inline

Get a tagged Kokkos vector.

Parameters
tagThe vector tag
Returns
The Kokkos vector

Definition at line 288 of file KokkosSystem.h.

288 { return _vectors[tag]; }
Array< Vector > _vectors
Kokkos vectors and matrices on device.
Definition: KokkosSystem.h:450

◆ getVectorDofValue()

KOKKOS_FUNCTION Real& Moose::Kokkos::System::getVectorDofValue ( dof_id_type  dof,
TagID  tag 
) const
inline

Get the DOF value of a tagged vector.

Parameters
dofThe local DOF index
tagThe vector tag
Returns
The DOF value

Definition at line 301 of file KokkosSystem.h.

Referenced by getVectorQpGradFace(), and getVectorQpValueFace().

302  {
303  return _vectors[tag][dof];
304  }
Array< Vector > _vectors
Kokkos vectors and matrices on device.
Definition: KokkosSystem.h:450

◆ getVectorQpGrad()

KOKKOS_FUNCTION Real3& Moose::Kokkos::System::getVectorQpGrad ( ElementInfo  info,
dof_id_type  qp,
unsigned int  var,
TagID  tag 
) const
inline

Get the quadrature gradient of a variable from a tagged vector.

Parameters
infoThe element information object
qpThe global quadrature point index
varThe variable number
tagThe vector tag
Returns
The quadrature gradient

Definition at line 327 of file KokkosSystem.h.

328  {
329  return _qp_solutions_grad[tag](info.subdomain, var)[qp];
330  }
MPI_Info info
Array< Array2D< Array< Real3 > > > _qp_solutions_grad
Definition: KokkosSystem.h:458

◆ getVectorQpGradFace()

KOKKOS_FUNCTION Real3 Moose::Kokkos::System::getVectorQpGradFace ( const ElementInfo  info,
const unsigned int  side,
const Real33  jacobian,
const unsigned int  qp,
const unsigned int  var,
const TagID  tag 
) const
inline

Get the face quadrature gradient of a variable from a tagged vector.

Parameters
infoThe element information object
sideThe side index
jacobianThe inverse Jacobian matrix
qpThe local quadrature point index
varThe variable number
tagThe vector tag
Returns
The face quadrature gradient

Definition at line 549 of file KokkosSystem.h.

555 {
556  auto fe = _var_fe_types[var];
557  auto n_dofs = kokkosAssembly().getNumDofs(info.type, fe);
558  auto & grad_phi = kokkosAssembly().getGradPhiFace(info.subdomain, info.type, fe)(side);
559 
560  Real3 grad = 0;
561 
562  for (unsigned int i = 0; i < n_dofs; ++i)
563  grad += getVectorDofValue(getElemLocalDofIndex(info.id, i, var), tag) *
564  (jacobian * grad_phi(i, qp));
565 
566  return grad;
567 }
KOKKOS_FUNCTION dof_id_type getElemLocalDofIndex(ContiguousElementID elem, unsigned int i, unsigned int var) const
Get the local DOF index of a variable for an element.
Definition: KokkosSystem.h:225
KOKKOS_FUNCTION const Assembly & kokkosAssembly() const
Get the const reference of the Kokkos assembly.
MPI_Info info
KOKKOS_FUNCTION const auto & getGradPhiFace(ContiguousSubdomainID subdomain, unsigned int elem_type, unsigned int fe_type) const
Get the gradient of face shape functions of a FE type for an element type and subdomain.
KOKKOS_FUNCTION unsigned int getNumDofs(unsigned int elem_type, unsigned int fe_type) const
Get the number of DOFs of a FE type for an element type.
Array< unsigned int > _var_fe_types
FE type ID of each variable.
Definition: KokkosSystem.h:478
KOKKOS_FUNCTION Real & getVectorDofValue(dof_id_type dof, TagID tag) const
Get the DOF value of a tagged vector.
Definition: KokkosSystem.h:301

◆ getVectorQpValue()

KOKKOS_FUNCTION Real& Moose::Kokkos::System::getVectorQpValue ( ElementInfo  info,
dof_id_type  qp,
unsigned int  var,
TagID  tag 
) const
inline

Get the quadrature value of a variable from a tagged vector.

Parameters
infoThe element information object
qpThe global quadrature point index
varThe variable number
tagThe vector tag
Returns
The quadrature value

Definition at line 314 of file KokkosSystem.h.

315  {
316  return _qp_solutions[tag](info.subdomain, var)[qp];
317  }
Array< Array2D< Array< Real > > > _qp_solutions
Cached elemental quadrature values and gradients.
Definition: KokkosSystem.h:457
MPI_Info info

◆ getVectorQpValueFace()

KOKKOS_FUNCTION Real Moose::Kokkos::System::getVectorQpValueFace ( const ElementInfo  info,
const unsigned int  side,
const unsigned int  qp,
const unsigned int  var,
const TagID  tag 
) const
inline

Get the face quadrature value of a variable from a tagged vector.

Parameters
infoThe element information object
sideThe side index
qpThe local quadrature point index
varThe vriable number
tagThe vector tag
Returns
The face quadrature value

Definition at line 534 of file KokkosSystem.h.

536 {
537  auto fe = _var_fe_types[var];
538  auto n_dofs = kokkosAssembly().getNumDofs(info.type, fe);
539  auto & phi = kokkosAssembly().getPhiFace(info.subdomain, info.type, fe)(side);
540 
541  Real value = 0;
542 
543  for (unsigned int i = 0; i < n_dofs; ++i)
544  value += getVectorDofValue(getElemLocalDofIndex(info.id, i, var), tag) * phi(i, qp);
545 
546  return value;
547 }
KOKKOS_FUNCTION dof_id_type getElemLocalDofIndex(ContiguousElementID elem, unsigned int i, unsigned int var) const
Get the local DOF index of a variable for an element.
Definition: KokkosSystem.h:225
KOKKOS_FUNCTION const Assembly & kokkosAssembly() const
Get the const reference of the Kokkos assembly.
MPI_Info info
KOKKOS_FUNCTION unsigned int getNumDofs(unsigned int elem_type, unsigned int fe_type) const
Get the number of DOFs of a FE type for an element type.
KOKKOS_FUNCTION const auto & getPhiFace(ContiguousSubdomainID subdomain, unsigned int elem_type, unsigned int fe_type) const
Get the face shape functions of a FE type for an element type and subdomain.
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
Array< unsigned int > _var_fe_types
FE type ID of each variable.
Definition: KokkosSystem.h:478
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
KOKKOS_FUNCTION Real & getVectorDofValue(dof_id_type dof, TagID tag) const
Get the DOF value of a tagged vector.
Definition: KokkosSystem.h:301

◆ hasNodalBC()

KOKKOS_FUNCTION bool Moose::Kokkos::System::hasNodalBC ( dof_id_type  dof) const
inline

Check whether a local DOF index is covered by a nodal BC.

Parameters
dofThe local DOF index
Returns
Whether the local DOF index is covered by a nodal BC

Definition at line 168 of file KokkosSystem.h.

169  {
170  return _nbc_dof.isAlloc() && _nbc_dof[dof];
171  }
Array< bool > _nbc_dof
Flag whether each DOF is covered by a nodal BC.
Definition: KokkosSystem.h:510

◆ hasNodalBCMatrixTag()

KOKKOS_FUNCTION bool Moose::Kokkos::System::hasNodalBCMatrixTag ( dof_id_type  dof,
TagID  tag 
) const
inline

Check whether a local DOF index is associated with a nodal BC for an extra matrix tag.

Parameters
dofThe local DOF index
tagThe extra matrix tag
Returns
Whether the local DOF index is covered by a nodal BC

Definition at line 188 of file KokkosSystem.h.

189  {
190  return _nbc_matrix_tag_dof[tag].isAlloc() && _nbc_matrix_tag_dof[tag][dof];
191  }
Array< Array< bool > > _nbc_matrix_tag_dof
Definition: KokkosSystem.h:516

◆ hasNodalBCResidualTag()

KOKKOS_FUNCTION bool Moose::Kokkos::System::hasNodalBCResidualTag ( dof_id_type  dof,
TagID  tag 
) const
inline

Check whether a local DOF index is covered by a nodal BC for an extra residual tag.

Parameters
dofThe local DOF index
tagThe extra residual tag
Returns
Whether the local DOF index is covered by a nodal BC

Definition at line 178 of file KokkosSystem.h.

179  {
180  return _nbc_residual_tag_dof[tag].isAlloc() && _nbc_residual_tag_dof[tag][dof];
181  }
Array< Array< bool > > _nbc_residual_tag_dof
Flag whether each DOF is covered by a nodal BC for extra residual tags.
Definition: KokkosSystem.h:515

◆ isMatrixTagActive()

KOKKOS_FUNCTION bool Moose::Kokkos::System::isMatrixTagActive ( TagID  tag) const
inline

Check whether a matrix tag is active.

Parameters
tagThe matrix tag
Returns
Whether the matrix tag is active

Definition at line 162 of file KokkosSystem.h.

162 { return _matrix_tag_active[tag]; }
Array< bool > _matrix_tag_active
Definition: KokkosSystem.h:505

◆ isNodalDefined()

KOKKOS_FUNCTION bool Moose::Kokkos::System::isNodalDefined ( ContiguousNodeID  node,
unsigned int  var 
) const
inline

Get whether a variable is defined on a node.

Parameters
nodeThe contiguous node ID
varThe variable number
Returns
Whether the variable is defined on the node

Definition at line 279 of file KokkosSystem.h.

280  {
282  }
Array< Array< dof_id_type > > _local_node_dof_index
Definition: KokkosSystem.h:465
static const dof_id_type invalid_id

◆ isResidualTagActive()

KOKKOS_FUNCTION bool Moose::Kokkos::System::isResidualTagActive ( TagID  tag) const
inline

Check whether a residual tag is active.

Parameters
tagThe residual tag
Returns
Whether the residual tag is active

Definition at line 156 of file KokkosSystem.h.

156 { return _residual_tag_active[tag]; }
Array< bool > _residual_tag_active
Flag whether each tag is active.
Definition: KokkosSystem.h:504

◆ isVariableActive()

KOKKOS_FUNCTION bool Moose::Kokkos::System::isVariableActive ( unsigned int  var,
ContiguousSubdomainID  subdomain 
) const
inline

Check whether a variable is active on a subdomain.

Parameters
varThe variable number
subdomainThe contiguous subdomain ID
Returns
Whether the variable is active

Definition at line 147 of file KokkosSystem.h.

148  {
149  return _var_subdomain_active(var, subdomain);
150  }
Array2D< bool > _var_subdomain_active
Whether each variable is active on subdomains.
Definition: KokkosSystem.h:482

◆ kokkosAssembly()

KOKKOS_FUNCTION const Assembly& Moose::Kokkos::AssemblyHolder::kokkosAssembly ( ) const
inlineinherited

Get the const reference of the Kokkos assembly.

Returns
The const reference of the Kokkos assembly depending on the architecture this function is being called on

Definition at line 499 of file KokkosAssembly.h.

Referenced by getVectorQpGradFace(), and getVectorQpValueFace().

500  {
501  KOKKOS_IF_ON_HOST(return _assembly_host;)
502 
503  return _assembly_device;
504  }
const Assembly _assembly_device
Device copy of the Kokkos assembly.
const Assembly & _assembly_host
Host reference of the Kokkos assembly.

◆ kokkosMesh()

KOKKOS_FUNCTION const Mesh& Moose::Kokkos::MeshHolder::kokkosMesh ( ) const
inlineinherited

Get the const reference of the Kokkos mesh.

Returns
The const reference of the Kokkos mesh depending on the architecture this function is being called on

Definition at line 360 of file KokkosMesh.h.

Referenced by Moose::Kokkos::Assembly::computePhysicalMap().

361  {
362  KOKKOS_IF_ON_HOST(return _mesh_host;)
363 
364  return _mesh_device;
365  }
const Mesh _mesh_device
Device copy of the Kokkos mesh.
Definition: KokkosMesh.h:376
const Mesh & _mesh_host
Host reference of the Kokkos mesh.
Definition: KokkosMesh.h:372

◆ localToGlobalDofIndex()

KOKKOS_FUNCTION dof_id_type Moose::Kokkos::System::localToGlobalDofIndex ( dof_id_type  dof) const
inline

Get the global DOF index of a local DOF index.

Parameters
dofThe local DOF index
Returns
The global DOF index

Definition at line 269 of file KokkosSystem.h.

270  {
271  return _local_to_global_dof_index[dof];
272  }
Array< dof_id_type > _local_to_global_dof_index
Map from local DOF index to global DOF index.
Definition: KokkosSystem.h:470

◆ operator()()

KOKKOS_FUNCTION void Moose::Kokkos::System::operator() ( const ThreadID  tid) const

Kokkos function for caching variable values on element quadrature points.

◆ reinit()

void Moose::Kokkos::System::reinit ( )

Allocate the quadrature point solution vectors for active variable and tags and cache quadrature point values.

◆ setActiveMatrixTags()

void Moose::Kokkos::System::setActiveMatrixTags ( const std::set< TagID > &  tags)

Set the active matrix tags.

Parameters
varsThe active matrix tags

◆ setActiveResidualTags()

void Moose::Kokkos::System::setActiveResidualTags ( const std::set< TagID > &  tags)

Set the active residual tags.

Parameters
tagsThe active residual tags

◆ setActiveVariables()

void Moose::Kokkos::System::setActiveVariables ( const std::set< MooseVariableFieldBase *> &  vars)

Set the active variables.

Parameters
varsThe active MOOSE variables

◆ setActiveVariableTags()

void Moose::Kokkos::System::setActiveVariableTags ( const std::set< TagID > &  tags)

Set the active solution tags.

Parameters
tagsThe active solution tags

◆ setupDofs()

void Moose::Kokkos::System::setupDofs ( )
private

Setup DOF data.

◆ setupSparsity()

void Moose::Kokkos::System::setupSparsity ( )
private

Setup sparsity data.

◆ setupVariables()

void Moose::Kokkos::System::setupVariables ( )
private

Setup variable data.

◆ sync() [1/2]

void Moose::Kokkos::System::sync ( const MemcpyKind  dir)

Synchronize the active tagged vectors and matrices between host and device.

Parameters
dirCopy direction

◆ sync() [2/2]

void Moose::Kokkos::System::sync ( const std::set< TagID > &  tags,
const MemcpyKind  dir 
)

Synchronize the specified tagged vectors between host and device.

Parameters
tagsThe vector tags
dirCopy direction

Member Data Documentation

◆ _active_matrix_tags

Array<TagID> Moose::Kokkos::System::_active_matrix_tags
private

Definition at line 498 of file KokkosSystem.h.

Referenced by clearActiveMatrixTags().

◆ _active_residual_tags

Array<TagID> Moose::Kokkos::System::_active_residual_tags
private

Definition at line 497 of file KokkosSystem.h.

Referenced by clearActiveResidualTags().

◆ _active_variable_tags

Array<TagID> Moose::Kokkos::System::_active_variable_tags
private

List of active tags.

Definition at line 496 of file KokkosSystem.h.

Referenced by clearActiveVariableTags().

◆ _active_variables

Array<unsigned int> Moose::Kokkos::System::_active_variables
private

List of active variable numbers.

Definition at line 491 of file KokkosSystem.h.

Referenced by clearActiveVariables().

◆ _comm

const Parallel::Communicator& Moose::Kokkos::System::_comm
private

Reference of the libMesh communicator.

Definition at line 432 of file KokkosSystem.h.

Referenced by getComm().

◆ _coupling

Array<Array<unsigned int> > Moose::Kokkos::System::_coupling
private

Off-diagonal coupled variable numbers of each variable.

Definition at line 486 of file KokkosSystem.h.

Referenced by getCoupling().

◆ _dof_map

const libMesh::DofMap& Moose::Kokkos::System::_dof_map
private

Reference of the libMesh DOF map.

Definition at line 428 of file KokkosSystem.h.

◆ _ghost_comm_list

Array<Array<dof_id_type> > Moose::Kokkos::System::_ghost_comm_list
private

Definition at line 523 of file KokkosSystem.h.

Referenced by getGhostCommList().

◆ _local_comm_list

Array<Array<dof_id_type> > Moose::Kokkos::System::_local_comm_list
private

List of DOFs to send and receive.

Definition at line 522 of file KokkosSystem.h.

Referenced by getLocalCommList().

◆ _local_elem_dof_index

Array<Array2D<dof_id_type> > Moose::Kokkos::System::_local_elem_dof_index
private

Local DOF index of each variable.

Definition at line 464 of file KokkosSystem.h.

Referenced by getElemGlobalDofIndex(), and getElemLocalDofIndex().

◆ _local_node_dof_index

Array<Array<dof_id_type> > Moose::Kokkos::System::_local_node_dof_index
private

Definition at line 465 of file KokkosSystem.h.

Referenced by getNodeGlobalDofIndex(), getNodeLocalDofIndex(), and isNodalDefined().

◆ _local_to_global_dof_index

Array<dof_id_type> Moose::Kokkos::System::_local_to_global_dof_index
private

Map from local DOF index to global DOF index.

Definition at line 470 of file KokkosSystem.h.

Referenced by getElemGlobalDofIndex(), getNodeGlobalDofIndex(), and localToGlobalDofIndex().

◆ _matrices

Array<Matrix> Moose::Kokkos::System::_matrices
private

Definition at line 451 of file KokkosSystem.h.

Referenced by getMatrix(), and getMatrixValue().

◆ _matrix_tag_active

Array<bool> Moose::Kokkos::System::_matrix_tag_active
private

Definition at line 505 of file KokkosSystem.h.

Referenced by clearActiveMatrixTags(), and isMatrixTagActive().

◆ _max_dofs_per_elem

Array<unsigned int> Moose::Kokkos::System::_max_dofs_per_elem
private

Maximum number of DOFs per element for each variable.

Definition at line 474 of file KokkosSystem.h.

Referenced by getMaxDofsPerElem().

◆ _mesh

const MooseMesh& Moose::Kokkos::System::_mesh
private

Reference of the MOOSE mesh.

Definition at line 424 of file KokkosSystem.h.

◆ _nbc_dof

Array<bool> Moose::Kokkos::System::_nbc_dof
private

Flag whether each DOF is covered by a nodal BC.

Definition at line 510 of file KokkosSystem.h.

Referenced by hasNodalBC().

◆ _nbc_matrix_tag_dof

Array<Array<bool> > Moose::Kokkos::System::_nbc_matrix_tag_dof
private

Definition at line 516 of file KokkosSystem.h.

Referenced by hasNodalBCMatrixTag().

◆ _nbc_residual_tag_dof

Array<Array<bool> > Moose::Kokkos::System::_nbc_residual_tag_dof
private

Flag whether each DOF is covered by a nodal BC for extra residual tags.

Definition at line 515 of file KokkosSystem.h.

Referenced by hasNodalBCResidualTag().

◆ _num_ghost_dofs

const dof_id_type Moose::Kokkos::System::_num_ghost_dofs
private

Number of ghost DOFs.

Definition at line 444 of file KokkosSystem.h.

Referenced by getNumGhostDofs().

◆ _num_local_dofs

const dof_id_type Moose::Kokkos::System::_num_local_dofs
private

Number of local DOFs.

Definition at line 440 of file KokkosSystem.h.

Referenced by getNumLocalDofs().

◆ _num_vars

const unsigned int Moose::Kokkos::System::_num_vars
private

Number of variables.

Definition at line 436 of file KokkosSystem.h.

◆ _qp_solutions

Array<Array2D<Array<Real> > > Moose::Kokkos::System::_qp_solutions
private

Cached elemental quadrature values and gradients.

Definition at line 457 of file KokkosSystem.h.

Referenced by getVectorQpValue().

◆ _qp_solutions_grad

Array<Array2D<Array<Real3> > > Moose::Kokkos::System::_qp_solutions_grad
private

Definition at line 458 of file KokkosSystem.h.

Referenced by getVectorQpGrad().

◆ _residual_tag_active

Array<bool> Moose::Kokkos::System::_residual_tag_active
private

Flag whether each tag is active.

Definition at line 504 of file KokkosSystem.h.

Referenced by clearActiveResidualTags(), and isResidualTagActive().

◆ _sparsity

Sparsity Moose::Kokkos::System::_sparsity
private

Matrix sparsity pattern data.

Definition at line 529 of file KokkosSystem.h.

Referenced by getSparsity().

◆ _system

SystemBase& Moose::Kokkos::System::_system
private

Reference of the MOOSE system.

Definition at line 420 of file KokkosSystem.h.

Referenced by getDofMap(), and getSystem().

◆ _thread

Thread Moose::Kokkos::System::_thread
private

Kokkos thread object.

Definition at line 416 of file KokkosSystem.h.

◆ _var_fe_types

Array<unsigned int> Moose::Kokkos::System::_var_fe_types
private

FE type ID of each variable.

Definition at line 478 of file KokkosSystem.h.

Referenced by getFETypeID(), getVectorQpGradFace(), and getVectorQpValueFace().

◆ _var_subdomain_active

Array2D<bool> Moose::Kokkos::System::_var_subdomain_active
private

Whether each variable is active on subdomains.

Definition at line 482 of file KokkosSystem.h.

Referenced by isVariableActive().

◆ _vectors

Array<Vector> Moose::Kokkos::System::_vectors
private

Kokkos vectors and matrices on device.

Definition at line 450 of file KokkosSystem.h.

Referenced by getVector(), and getVectorDofValue().


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