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);
57  void sync(const std::vector<TagID> & tags, const MemcpyKind dir);
58  void sync(const TagID tag, const MemcpyKind dir);
60 
64  void reinit();
69  void setActiveVariables(const std::set<MooseVariableFieldBase *> & vars);
74  void setActiveSolutionTags(const std::set<TagID> & tags);
79  void setActiveResidualTags(const std::set<TagID> & tags);
84  void setActiveMatrixTags(const std::set<TagID> & tags);
97  {
98  _active_residual_tags.destroy();
99  _residual_tag_active = false;
100  }
105  {
106  _active_matrix_tags.destroy();
107  _matrix_tag_active = false;
108  }
113  auto & getSystem() { return _system; }
115  const auto & getSystem() const { return _system; }
117 
121  const auto & getDofMap() const { return _system.dofMap(); }
126  const auto & getComm() const { return _comm; }
131  const auto & getLocalCommList() const { return _local_comm_list; }
136  const auto & getGhostCommList() const { return _ghost_comm_list; }
141  const auto & getSparsity() const { return _sparsity; }
147  KOKKOS_FUNCTION const auto & getCoupling(unsigned int var) const { return _coupling[var]; }
154  KOKKOS_FUNCTION bool isVariableActive(unsigned int var, ContiguousSubdomainID subdomain) const
155  {
156  return _var_subdomain_active(var, subdomain);
157  }
163  KOKKOS_FUNCTION bool isResidualTagActive(TagID tag) const { return _residual_tag_active[tag]; }
169  KOKKOS_FUNCTION bool isMatrixTagActive(TagID tag) const { return _matrix_tag_active[tag]; }
175  KOKKOS_FUNCTION bool hasNodalBC(dof_id_type dof) const
176  {
177  return _nbc_dof.isAlloc() && _nbc_dof[dof];
178  }
185  KOKKOS_FUNCTION bool hasNodalBCResidualTag(dof_id_type dof, TagID tag) const
186  {
187  return _nbc_residual_tag_dof[tag].isAlloc() && _nbc_residual_tag_dof[tag][dof];
188  }
195  KOKKOS_FUNCTION bool hasNodalBCMatrixTag(dof_id_type dof, TagID tag) const
196  {
197  return _nbc_matrix_tag_dof[tag].isAlloc() && _nbc_matrix_tag_dof[tag][dof];
198  }
204  KOKKOS_FUNCTION unsigned int getFETypeID(unsigned int var) const { return _var_fe_types[var]; }
209  KOKKOS_FUNCTION dof_id_type getNumLocalDofs() const { return _num_local_dofs; }
214  KOKKOS_FUNCTION dof_id_type getNumGhostDofs() const { return _num_ghost_dofs; }
221  KOKKOS_FUNCTION unsigned int getMaxDofsPerElem(unsigned int var) const
222  {
223  return _max_dofs_per_elem[var];
224  }
233  unsigned int i,
234  unsigned int var) const
235  {
236  return _local_elem_dof_index[var](elem, i);
237  }
246  unsigned int i,
247  unsigned int var) const
248  {
249  return _local_node_dof_index[var][node] + i;
250  }
259  unsigned int i,
260  unsigned int var) const
261  {
263  }
270  KOKKOS_FUNCTION dof_id_type getNodeGlobalDofIndex(ContiguousNodeID node, unsigned int var) const
271  {
273  }
279  KOKKOS_FUNCTION dof_id_type localToGlobalDofIndex(dof_id_type dof) const
280  {
281  return _local_to_global_dof_index[dof];
282  }
289  KOKKOS_FUNCTION bool isNodalDefined(ContiguousNodeID node, unsigned int var) const
290  {
292  }
298  KOKKOS_FUNCTION auto & getVector(TagID tag) const { return _vectors[tag]; }
304  KOKKOS_FUNCTION auto & getMatrix(TagID tag) const { return _matrices[tag]; }
311  KOKKOS_FUNCTION Real & getVectorDofValue(dof_id_type dof, TagID tag) const
312  {
313  return _vectors[tag][dof];
314  }
323  KOKKOS_FUNCTION Real &
324  getVectorQpValue(ElementInfo info, dof_id_type qp, unsigned int var, TagID tag) const
325  {
326  return _qp_solutions[tag](info.subdomain, var)[qp];
327  }
336  KOKKOS_FUNCTION Real3 &
337  getVectorQpGrad(ElementInfo info, dof_id_type qp, unsigned int var, TagID tag) const
338  {
339  return _qp_solutions_grad[tag](info.subdomain, var)[qp];
340  }
350  KOKKOS_FUNCTION Real getVectorQpValueFace(const ElementInfo info,
351  const unsigned int side,
352  const unsigned int qp,
353  const unsigned int var,
354  const TagID tag) const;
365  KOKKOS_FUNCTION Real3 getVectorQpGradFace(const ElementInfo info,
366  const unsigned int side,
367  const Real33 jacobian,
368  const unsigned int qp,
369  const unsigned int var,
370  const TagID tag) const;
378  KOKKOS_FUNCTION Real & getMatrixValue(dof_id_type row, dof_id_type col, TagID tag) const
379  {
380  return _matrices[tag](row, col);
381  }
382 
386  KOKKOS_FUNCTION void operator()(const ThreadID tid) const;
387 #endif
388 
392  struct Sparsity
393  {
397  };
398 
399 private:
403  void setupVariables();
407  void setupDofs();
411  void setupSparsity();
415  void checkNodalBCs();
421  void getNodalBCDofs(const NodalBCBase * nbc, Array<bool> & dofs);
422 
434  const MooseMesh & _mesh;
446  const unsigned int _num_vars;
455 
463 
470 
477 
497 
510 
517 
528 
535 
540 };
541 
542 #ifdef MOOSE_KOKKOS_SCOPE
543 KOKKOS_FUNCTION inline Real
545  ElementInfo info, unsigned int side, unsigned int qp, unsigned int var, TagID tag) const
546 {
547  auto fe = _var_fe_types[var];
548  auto n_dofs = kokkosAssembly().getNumDofs(info.type, fe);
549  auto & phi = kokkosAssembly().getPhiFace(info.subdomain, info.type, fe)(side);
550 
551  Real value = 0;
552 
553  for (unsigned int i = 0; i < n_dofs; ++i)
554  value += getVectorDofValue(getElemLocalDofIndex(info.id, i, var), tag) * phi(i, qp);
555 
556  return value;
557 }
558 KOKKOS_FUNCTION inline Real3
560  unsigned int side,
561  Real33 jacobian,
562  unsigned int qp,
563  unsigned int var,
564  TagID tag) const
565 {
566  auto fe = _var_fe_types[var];
567  auto n_dofs = kokkosAssembly().getNumDofs(info.type, fe);
568  auto & grad_phi = kokkosAssembly().getGradPhiFace(info.subdomain, info.type, fe)(side);
569 
570  Real3 grad = 0;
571 
572  for (unsigned int i = 0; i < n_dofs; ++i)
573  grad += getVectorDofValue(getElemLocalDofIndex(info.id, i, var), tag) *
574  (jacobian * grad_phi(i, qp));
575 
576  return grad;
577 }
578 #endif
579 
587 {
588 public:
593  SystemHolder(Array<System> & systems) : _systems_host(systems), _systems_device(systems) {}
597  SystemHolder(const SystemHolder & holder)
599  {
600  }
601 
602 #ifdef MOOSE_KOKKOS_SCOPE
603 
608  KOKKOS_FUNCTION const Array<System> & kokkosSystems() const
609  {
610  KOKKOS_IF_ON_HOST(return _systems_host;)
611 
612  return _systems_device;
613  }
625  KOKKOS_FUNCTION const System & kokkosSystem(unsigned int sys) const
626  {
627  KOKKOS_IF_ON_HOST(return _systems_host[sys];)
628 
629  return _systems_device[sys];
630  }
636  System & kokkosSystem(unsigned int sys) { return _systems_host[sys]; }
637 #endif
638 
639 private:
648 };
649 
650 } // namespace Kokkos
651 } // namespace Moose
SystemHolder(Array< System > &systems)
Constructor.
Definition: KokkosSystem.h:593
KOKKOS_FUNCTION auto & getVector(TagID tag) const
Get a tagged Kokkos vector.
Definition: KokkosSystem.h:298
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:232
KOKKOS_FUNCTION dof_id_type getNodeLocalDofIndex(ContiguousNodeID node, unsigned int i, unsigned int var) const
Get the local DOF index of a variable for a node.
Definition: KokkosSystem.h:245
Array< Array2D< Array< Real > > > _qp_solutions
Cached elemental quadrature values and gradients.
Definition: KokkosSystem.h:467
const auto & getLocalCommList() const
Get the list of local DOF indices to communicate.
Definition: KokkosSystem.h:131
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:378
void clearActiveResidualTags()
Clear the cached active residual tags.
Definition: KokkosSystem.h:96
Array< Array< bool > > _nbc_matrix_tag_dof
Definition: KokkosSystem.h:526
Array< dof_id_type > _local_to_global_dof_index
Map from local DOF index to global DOF index.
Definition: KokkosSystem.h:480
Array< bool > _residual_tag_active
Flag whether each tag is active.
Definition: KokkosSystem.h:514
Array< bool > _nbc_dof
Flag whether each DOF is covered by a nodal BC.
Definition: KokkosSystem.h:520
void setActiveSolutionTags(const std::set< TagID > &tags)
Set the active solution tags.
Array< Array< bool > > _nbc_residual_tag_dof
Flag whether each DOF is covered by a nodal BC for extra residual tags.
Definition: KokkosSystem.h:525
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:426
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:136
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:221
KOKKOS_FUNCTION dof_id_type getNumGhostDofs() const
Get the number of ghost DOFs.
Definition: KokkosSystem.h:214
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:337
KOKKOS_FUNCTION dof_id_type getNumLocalDofs() const
Get the number of local DOFs.
Definition: KokkosSystem.h:209
char ** vars
const auto & getComm() const
Get the libMesh communicator.
Definition: KokkosSystem.h:126
const MooseMesh & _mesh
Reference of the MOOSE mesh.
Definition: KokkosSystem.h:434
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.
Array< Array< dof_id_type > > _local_node_dof_index
Definition: KokkosSystem.h:475
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:586
const auto & getSparsity() const
Get the sparisty pattern data.
Definition: KokkosSystem.h:141
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:544
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:508
const auto & getSystem() const
Definition: KokkosSystem.h:115
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:430
void clearActiveMatrixTags()
Clear the cached active matrix tags.
Definition: KokkosSystem.h:104
void setupVariables()
Setup variable data.
Array< Vector > _vectors
Kokkos vectors and matrices on device.
Definition: KokkosSystem.h:460
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:279
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:154
static constexpr dof_id_type invalid_id
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:163
const auto & getDofMap() const
Get the libMesh DOF map.
Definition: KokkosSystem.h:121
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:446
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:559
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:304
KOKKOS_FUNCTION bool isNodalDefined(ContiguousNodeID node, unsigned int var) const
Get whether a variable is defined on a node.
Definition: KokkosSystem.h:289
void clearActiveSolutionTags()
Clear the cached active solution tags.
Definition: KokkosSystem.h:92
void clearActiveVariables()
Clear the cached active variables.
Definition: KokkosSystem.h:88
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:625
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:618
void setupDofs()
Setup DOF data.
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:185
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:324
Array< TagID > _active_residual_tags
Definition: KokkosSystem.h:507
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:195
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:258
SystemHolder(const SystemHolder &holder)
Copy constructor.
Definition: KokkosSystem.h:597
const dof_id_type _num_local_dofs
Number of local DOFs.
Definition: KokkosSystem.h:450
KOKKOS_FUNCTION bool isMatrixTagActive(TagID tag) const
Check whether a matrix tag is active.
Definition: KokkosSystem.h:169
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:647
KOKKOS_FUNCTION bool hasNodalBC(dof_id_type dof) const
Check whether a local DOF index is covered by a nodal BC.
Definition: KokkosSystem.h:175
Array< unsigned int > _var_fe_types
FE type ID of each variable.
Definition: KokkosSystem.h:488
Array< Array< dof_id_type > > _ghost_comm_list
Definition: KokkosSystem.h:533
Array< Array2D< Array< Real3 > > > _qp_solutions_grad
Definition: KokkosSystem.h:468
Array< Array< dof_id_type > > _local_comm_list
List of DOFs to send and receive.
Definition: KokkosSystem.h:532
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void reinit()
Allocate the quadrature point vectors for active variable and tags and cache quadrature point values...
auto & getSystem()
Get the MOOSE system.
Definition: KokkosSystem.h:114
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:270
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:501
const Parallel::Communicator & _comm
Reference of the libMesh communicator.
Definition: KokkosSystem.h:442
Array< System > & _systems_host
Host reference of the Kokkos systems.
Definition: KokkosSystem.h:643
CSR format sparsity data.
Definition: KokkosSystem.h:392
Array2D< bool > _var_subdomain_active
Whether each variable is active on subdomains.
Definition: KokkosSystem.h:492
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:438
const dof_id_type _num_ghost_dofs
Number of ghost DOFs.
Definition: KokkosSystem.h:454
System & kokkosSystem(unsigned int sys)
Get the writeable reference of a Kokkos system.
Definition: KokkosSystem.h:636
Array< Array< unsigned int > > _coupling
Off-diagonal coupled variable numbers of each variable.
Definition: KokkosSystem.h:496
KOKKOS_FUNCTION Real & getVectorDofValue(dof_id_type dof, TagID tag) const
Get the DOF value of a tagged vector.
Definition: KokkosSystem.h:311
KOKKOS_FUNCTION const auto & getCoupling(unsigned int var) const
Get the list of off-diagonal coupled variable numbers of a variable.
Definition: KokkosSystem.h:147
Array< TagID > _active_solution_tags
List of active tags.
Definition: KokkosSystem.h:506
Sparsity _sparsity
Matrix sparsity pattern data.
Definition: KokkosSystem.h:539
Array< unsigned int > _max_dofs_per_elem
Maximum number of DOFs per element for each variable.
Definition: KokkosSystem.h:484
Array< Array2D< dof_id_type > > _local_elem_dof_index
Local DOF index of each variable.
Definition: KokkosSystem.h:474
KOKKOS_FUNCTION unsigned int getFETypeID(unsigned int var) const
Get the FE type ID of a variable.
Definition: KokkosSystem.h:204
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:515
KOKKOS_FUNCTION const Array< System > & kokkosSystems() const
Get the const reference of the Kokkos systems.
Definition: KokkosSystem.h:608
uint8_t dof_id_type
Array< Matrix > _matrices
Definition: KokkosSystem.h:461