https://mooseframework.inl.gov
KokkosSystem.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 #include "KokkosTypes.h"
13 #include "KokkosVector.h"
14 #include "KokkosMatrix.h"
15 #include "KokkosAssembly.h"
16 
17 #include "SystemBase.h"
18 
19 #include "libmesh/communicator.h"
20 
21 class MooseMesh;
22 class SystemBase;
23 
24 namespace Moose
25 {
26 namespace Kokkos
27 {
28 
29 class NodalBCBase;
30 
35 class System : public MeshHolder, public AssemblyHolder
36 {
37 public:
42  System(SystemBase & system);
43 
44 #ifdef MOOSE_KOKKOS_SCOPE
45 
49  void sync(const MemcpyKind dir);
55  void sync(const std::set<TagID> & tags, const MemcpyKind dir);
60  void reinit();
65  void setActiveVariables(const std::set<MooseVariableFieldBase *> & vars);
70  void setActiveVariableTags(const std::set<TagID> & tags);
75  void setActiveResidualTags(const std::set<TagID> & tags);
80  void setActiveMatrixTags(const std::set<TagID> & tags);
93  {
94  _active_residual_tags.destroy();
95  _residual_tag_active = false;
96  }
101  {
102  _active_matrix_tags.destroy();
103  _matrix_tag_active = false;
104  }
109  const auto & getSystem() const { return _system; }
114  const auto & getDofMap() const { return _system.dofMap(); }
119  const auto & getComm() const { return _comm; }
124  const auto & getLocalCommList() const { return _local_comm_list; }
129  const auto & getGhostCommList() const { return _ghost_comm_list; }
134  const auto & getSparsity() const { return _sparsity; }
140  KOKKOS_FUNCTION const auto & getCoupling(unsigned int var) const { return _coupling[var]; }
147  KOKKOS_FUNCTION bool isVariableActive(unsigned int var, ContiguousSubdomainID subdomain) const
148  {
149  return _var_subdomain_active(var, subdomain);
150  }
156  KOKKOS_FUNCTION bool isResidualTagActive(TagID tag) const { return _residual_tag_active[tag]; }
162  KOKKOS_FUNCTION bool isMatrixTagActive(TagID tag) const { return _matrix_tag_active[tag]; }
168  KOKKOS_FUNCTION bool hasNodalBC(dof_id_type dof) const
169  {
170  return _nbc_dof.isAlloc() && _nbc_dof[dof];
171  }
178  KOKKOS_FUNCTION bool hasNodalBCResidualTag(dof_id_type dof, TagID tag) const
179  {
180  return _nbc_residual_tag_dof[tag].isAlloc() && _nbc_residual_tag_dof[tag][dof];
181  }
188  KOKKOS_FUNCTION bool hasNodalBCMatrixTag(dof_id_type dof, TagID tag) const
189  {
190  return _nbc_matrix_tag_dof[tag].isAlloc() && _nbc_matrix_tag_dof[tag][dof];
191  }
197  KOKKOS_FUNCTION unsigned int getFETypeID(unsigned int var) const { return _var_fe_types[var]; }
202  KOKKOS_FUNCTION dof_id_type getNumLocalDofs() const { return _num_local_dofs; }
207  KOKKOS_FUNCTION dof_id_type getNumGhostDofs() const { return _num_ghost_dofs; }
214  KOKKOS_FUNCTION unsigned int getMaxDofsPerElem(unsigned int var) const
215  {
216  return _max_dofs_per_elem[var];
217  }
226  unsigned int i,
227  unsigned int var) const
228  {
229  return _local_elem_dof_index[var](elem, i);
230  }
237  KOKKOS_FUNCTION dof_id_type getNodeLocalDofIndex(ContiguousNodeID node, unsigned int var) const
238  {
239  return _local_node_dof_index[var][node];
240  }
249  unsigned int i,
250  unsigned int var) const
251  {
253  }
260  KOKKOS_FUNCTION dof_id_type getNodeGlobalDofIndex(ContiguousNodeID node, unsigned int var) const
261  {
263  }
269  KOKKOS_FUNCTION dof_id_type localToGlobalDofIndex(dof_id_type dof) const
270  {
271  return _local_to_global_dof_index[dof];
272  }
279  KOKKOS_FUNCTION bool isNodalDefined(ContiguousNodeID node, unsigned int var) const
280  {
282  }
288  KOKKOS_FUNCTION auto & getVector(TagID tag) const { return _vectors[tag]; }
294  KOKKOS_FUNCTION auto & getMatrix(TagID tag) const { return _matrices[tag]; }
301  KOKKOS_FUNCTION Real & getVectorDofValue(dof_id_type dof, TagID tag) const
302  {
303  return _vectors[tag][dof];
304  }
313  KOKKOS_FUNCTION Real &
314  getVectorQpValue(ElementInfo info, dof_id_type qp, unsigned int var, TagID tag) const
315  {
316  return _qp_solutions[tag](info.subdomain, var)[qp];
317  }
326  KOKKOS_FUNCTION Real3 &
327  getVectorQpGrad(ElementInfo info, dof_id_type qp, unsigned int var, TagID tag) const
328  {
329  return _qp_solutions_grad[tag](info.subdomain, var)[qp];
330  }
340  KOKKOS_FUNCTION Real getVectorQpValueFace(const ElementInfo info,
341  const unsigned int side,
342  const unsigned int qp,
343  const unsigned int var,
344  const TagID tag) const;
355  KOKKOS_FUNCTION Real3 getVectorQpGradFace(const ElementInfo info,
356  const unsigned int side,
357  const Real33 jacobian,
358  const unsigned int qp,
359  const unsigned int var,
360  const TagID tag) const;
368  KOKKOS_FUNCTION Real & getMatrixValue(dof_id_type row, dof_id_type col, TagID tag) const
369  {
370  return _matrices[tag](row, col);
371  }
372 
376  KOKKOS_FUNCTION void operator()(const ThreadID tid) const;
377 #endif
378 
382  struct Sparsity
383  {
387  };
388 
389 private:
393  void setupVariables();
397  void setupDofs();
401  void setupSparsity();
405  void checkNodalBCs();
411  void getNodalBCDofs(const NodalBCBase * nbc, Array<bool> & dofs);
412 
424  const MooseMesh & _mesh;
436  const unsigned int _num_vars;
445 
453 
460 
467 
487 
500 
507 
518 
525 
530 };
531 
532 #ifdef MOOSE_KOKKOS_SCOPE
533 KOKKOS_FUNCTION inline Real
535  ElementInfo info, unsigned int side, unsigned int qp, unsigned int var, TagID tag) const
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 }
548 KOKKOS_FUNCTION inline Real3
550  unsigned int side,
551  Real33 jacobian,
552  unsigned int qp,
553  unsigned int var,
554  TagID tag) const
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 }
568 #endif
569 
577 {
578 public:
583  SystemHolder(Array<System> & systems) : _systems_host(systems), _systems_device(systems) {}
587  SystemHolder(const SystemHolder & holder)
589  {
590  }
591 
592 #ifdef MOOSE_KOKKOS_SCOPE
593 
598  KOKKOS_FUNCTION const Array<System> & kokkosSystems() const
599  {
600  KOKKOS_IF_ON_HOST(return _systems_host;)
601 
602  return _systems_device;
603  }
615  KOKKOS_FUNCTION const System & kokkosSystem(unsigned int sys) const
616  {
617  KOKKOS_IF_ON_HOST(return _systems_host[sys];)
618 
619  return _systems_device[sys];
620  }
626  System & kokkosSystem(unsigned int sys) { return _systems_host[sys]; }
627 #endif
628 
629 private:
638 };
639 
640 } // namespace Kokkos
641 } // namespace Moose
SystemHolder(Array< System > &systems)
Constructor.
Definition: KokkosSystem.h:583
KOKKOS_FUNCTION auto & getVector(TagID tag) const
Get a tagged Kokkos vector.
Definition: KokkosSystem.h:288
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
Array< Array2D< Array< Real > > > _qp_solutions
Cached elemental quadrature values and gradients.
Definition: KokkosSystem.h:457
const auto & getLocalCommList() const
Get the list of local DOF indices to communicate.
Definition: KokkosSystem.h:124
The Kokkos object that contains the information of an element The IDs used in Kokkos are different fr...
Definition: KokkosMesh.h:35
KOKKOS_FUNCTION Real & getMatrixValue(dof_id_type row, dof_id_type col, TagID tag) const
Get an entry from a tagged matrix.
Definition: KokkosSystem.h:368
void clearActiveResidualTags()
Clear the cached active residual tags.
Definition: KokkosSystem.h:92
Array< Array< bool > > _nbc_matrix_tag_dof
Definition: KokkosSystem.h:516
Array< dof_id_type > _local_to_global_dof_index
Map from local DOF index to global DOF index.
Definition: KokkosSystem.h:470
Array< bool > _residual_tag_active
Flag whether each tag is active.
Definition: KokkosSystem.h:504
Array< bool > _nbc_dof
Flag whether each DOF is covered by a nodal BC.
Definition: KokkosSystem.h:510
void clearActiveVariableTags()
Clear the cached active solution tags.
Definition: KokkosSystem.h:88
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
KOKKOS_FUNCTION const Assembly & kokkosAssembly() const
Get the const reference of the Kokkos assembly.
unsigned int TagID
Definition: MooseTypes.h:210
Thread _thread
Kokkos thread object.
Definition: KokkosSystem.h:416
MPI_Info info
dof_id_type ThreadID
Definition: KokkosThread.h:18
dof_id_type ContiguousElementID
Definition: KokkosMesh.h:20
const auto & getGhostCommList() const
Get the list of ghost DOF indices to communicate.
Definition: KokkosSystem.h:129
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...
Definition: KokkosSystem.h:214
KOKKOS_FUNCTION dof_id_type getNumGhostDofs() const
Get the number of ghost DOFs.
Definition: KokkosSystem.h:207
KOKKOS_FUNCTION Real3 & getVectorQpGrad(ElementInfo info, dof_id_type qp, unsigned int var, TagID tag) const
Get the quadrature gradient of a variable from a tagged vector.
Definition: KokkosSystem.h:327
KOKKOS_FUNCTION dof_id_type getNumLocalDofs() const
Get the number of local DOFs.
Definition: KokkosSystem.h:202
char ** vars
KOKKOS_FUNCTION dof_id_type getNodeLocalDofIndex(ContiguousNodeID node, unsigned int var) const
Get the local DOF index of a variable for a node.
Definition: KokkosSystem.h:237
const auto & getComm() const
Get the libMesh communicator.
Definition: KokkosSystem.h:119
const MooseMesh & _mesh
Reference of the MOOSE mesh.
Definition: KokkosSystem.h:424
MemcpyKind
The enumerator that dictates the memory copy direction.
Definition: KokkosArray.h:44
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.
void setActiveVariableTags(const std::set< TagID > &tags)
Set the active solution tags.
Array< Array< dof_id_type > > _local_node_dof_index
Definition: KokkosSystem.h:465
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.
The Kokkos interface that holds the host reference of the Kokkos systems and copies it to device duri...
Definition: KokkosSystem.h:576
const auto & getSparsity() const
Get the sparisty pattern data.
Definition: KokkosSystem.h:134
void setupSparsity()
Setup sparsity data.
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.
Definition: KokkosSystem.h:534
Base class for a system (of equations)
Definition: SystemBase.h:84
The Kokkos system class.
Definition: KokkosSystem.h:35
The Kokkos interface that holds the host reference of the Kokkos mesh and copies it to device during ...
Definition: KokkosMesh.h:338
Array< TagID > _active_matrix_tags
Definition: KokkosSystem.h:498
const auto & getSystem() const
Get the MOOSE system.
Definition: KokkosSystem.h:109
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.
SystemBase & _system
Reference of the MOOSE system.
Definition: KokkosSystem.h:420
void clearActiveMatrixTags()
Clear the cached active matrix tags.
Definition: KokkosSystem.h:100
Array< TagID > _active_variable_tags
List of active tags.
Definition: KokkosSystem.h:496
void setupVariables()
Setup variable data.
Array< Vector > _vectors
Kokkos vectors and matrices on device.
Definition: KokkosSystem.h:450
virtual libMesh::DofMap & dofMap()
Gets writeable reference to the dof map.
Definition: SystemBase.C:1163
KOKKOS_FUNCTION dof_id_type localToGlobalDofIndex(dof_id_type dof) const
Get the global DOF index of a local DOF index.
Definition: KokkosSystem.h:269
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
KOKKOS_FUNCTION bool isVariableActive(unsigned int var, ContiguousSubdomainID subdomain) const
Check whether a variable is active on a subdomain.
Definition: KokkosSystem.h:147
dof_id_type ContiguousNodeID
Definition: KokkosMesh.h:21
KOKKOS_FUNCTION bool isResidualTagActive(TagID tag) const
Check whether a residual tag is active.
Definition: KokkosSystem.h:156
const auto & getDofMap() const
Get the libMesh DOF map.
Definition: KokkosSystem.h:114
void getNodalBCDofs(const NodalBCBase *nbc, Array< bool > &dofs)
Get the list of DOFs covered by a nodal BC.
System(SystemBase &system)
Constructor.
const unsigned int _num_vars
Number of variables.
Definition: KokkosSystem.h:436
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.
Definition: KokkosSystem.h:549
void checkNodalBCs()
Check if the DOFs are covered by nodal BCs.
KOKKOS_FUNCTION auto & getMatrix(TagID tag) const
Get a tagged Kokkos matrix.
Definition: KokkosSystem.h:294
KOKKOS_FUNCTION bool isNodalDefined(ContiguousNodeID node, unsigned int var) const
Get whether a variable is defined on a node.
Definition: KokkosSystem.h:279
void clearActiveVariables()
Clear the cached active variables.
Definition: KokkosSystem.h:84
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
Definition: MooseMesh.h:92
KOKKOS_FUNCTION const System & kokkosSystem(unsigned int sys) const
Get the const reference of a Kokkos system.
Definition: KokkosSystem.h:615
The Kokkos interface that holds the host reference of the Kokkos assembly and copies it to device dur...
Array< System > & kokkosSystems()
Get the writeable host reference of the Kokkos systems.
Definition: KokkosSystem.h:608
void setupDofs()
Setup DOF data.
static const dof_id_type invalid_id
void setActiveVariables(const std::set< MooseVariableFieldBase *> &vars)
Set the active variables.
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.
Definition: KokkosSystem.h:178
KOKKOS_FUNCTION Real & getVectorQpValue(ElementInfo info, dof_id_type qp, unsigned int var, TagID tag) const
Get the quadrature value of a variable from a tagged vector.
Definition: KokkosSystem.h:314
Array< TagID > _active_residual_tags
Definition: KokkosSystem.h:497
void setActiveResidualTags(const std::set< TagID > &tags)
Set the active residual tags.
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...
Definition: KokkosSystem.h:188
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.
Definition: KokkosSystem.h:248
SystemHolder(const SystemHolder &holder)
Copy constructor.
Definition: KokkosSystem.h:587
const dof_id_type _num_local_dofs
Number of local DOFs.
Definition: KokkosSystem.h:440
KOKKOS_FUNCTION bool isMatrixTagActive(TagID tag) const
Check whether a matrix tag is active.
Definition: KokkosSystem.h:162
void setActiveMatrixTags(const std::set< TagID > &tags)
Set the active matrix tags.
The base class for Kokkos nodal boundary conditions.
void sync(const MemcpyKind dir)
Synchronize the active tagged vectors and matrices between host and device.
Base class for deriving any boundary condition that works at nodes.
Definition: NodalBCBase.h:26
const Array< System > _systems_device
Device copy of the Kokkos systems.
Definition: KokkosSystem.h:637
KOKKOS_FUNCTION bool hasNodalBC(dof_id_type dof) const
Check whether a local DOF index is covered by a nodal BC.
Definition: KokkosSystem.h:168
Array< unsigned int > _var_fe_types
FE type ID of each variable.
Definition: KokkosSystem.h:478
Array< Array< dof_id_type > > _ghost_comm_list
Definition: KokkosSystem.h:523
Array< Array2D< Array< Real3 > > > _qp_solutions_grad
Definition: KokkosSystem.h:458
Array< Array< dof_id_type > > _local_comm_list
List of DOFs to send and receive.
Definition: KokkosSystem.h:522
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void reinit()
Allocate the quadrature point solution vectors for active variable and tags and cache quadrature poin...
KOKKOS_FUNCTION dof_id_type getNodeGlobalDofIndex(ContiguousNodeID node, unsigned int var) const
Get the global DOF index of a variable for a node.
Definition: KokkosSystem.h:260
KOKKOS_FUNCTION void operator()(const ThreadID tid) const
Kokkos function for caching variable values on element quadrature points.
Array< unsigned int > _active_variables
List of active variable numbers.
Definition: KokkosSystem.h:491
const Parallel::Communicator & _comm
Reference of the libMesh communicator.
Definition: KokkosSystem.h:432
Array< System > & _systems_host
Host reference of the Kokkos systems.
Definition: KokkosSystem.h:633
CSR format sparsity data.
Definition: KokkosSystem.h:382
Array2D< bool > _var_subdomain_active
Whether each variable is active on subdomains.
Definition: KokkosSystem.h:482
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
const libMesh::DofMap & _dof_map
Reference of the libMesh DOF map.
Definition: KokkosSystem.h:428
const dof_id_type _num_ghost_dofs
Number of ghost DOFs.
Definition: KokkosSystem.h:444
System & kokkosSystem(unsigned int sys)
Get the writeable reference of a Kokkos system.
Definition: KokkosSystem.h:626
Array< Array< unsigned int > > _coupling
Off-diagonal coupled variable numbers of each variable.
Definition: KokkosSystem.h:486
KOKKOS_FUNCTION Real & getVectorDofValue(dof_id_type dof, TagID tag) const
Get the DOF value of a tagged vector.
Definition: KokkosSystem.h:301
KOKKOS_FUNCTION const auto & getCoupling(unsigned int var) const
Get the list of off-diagonal coupled variable numbers of a variable.
Definition: KokkosSystem.h:140
Sparsity _sparsity
Matrix sparsity pattern data.
Definition: KokkosSystem.h:529
Array< unsigned int > _max_dofs_per_elem
Maximum number of DOFs per element for each variable.
Definition: KokkosSystem.h:474
Array< Array2D< dof_id_type > > _local_elem_dof_index
Local DOF index of each variable.
Definition: KokkosSystem.h:464
KOKKOS_FUNCTION unsigned int getFETypeID(unsigned int var) const
Get the FE type ID of a variable.
Definition: KokkosSystem.h:197
The Kokkos thread object that aids in converting the one-dimensional thread index into multi-dimensio...
Definition: KokkosThread.h:29
Array< bool > _matrix_tag_active
Definition: KokkosSystem.h:505
KOKKOS_FUNCTION const Array< System > & kokkosSystems() const
Get the const reference of the Kokkos systems.
Definition: KokkosSystem.h:598
uint8_t dof_id_type
Array< Matrix > _matrices
Definition: KokkosSystem.h:451