https://mooseframework.inl.gov
KokkosNodalBC.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 "KokkosNodalBCBase.h"
13 
14 namespace Moose
15 {
16 namespace Kokkos
17 {
18 
36 class NodalBC : public NodalBCBase
37 {
38 public:
40 
45 
49  virtual void computeResidual() override;
53  virtual void computeJacobian() override;
54 
59 
66  KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int /* qp */,
67  AssemblyDatum & /* datum */) const
68  {
69  return 1;
70  }
78  KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int /* jvar */,
79  const unsigned int /* qp */,
80  AssemblyDatum & /* datum */) const
81  {
82  return 0;
83  }
88  static auto defaultJacobian() { return &NodalBC::computeQpJacobian; }
95 
99  template <typename Derived>
101  KOKKOS_FUNCTION Real computeQpResidualShim(const Derived & bc,
102  const unsigned int qp,
103  AssemblyDatum & datum) const
104  {
105  return bc.computeQpResidual(qp, datum);
106  }
107  template <typename Derived>
108  KOKKOS_FUNCTION Real computeQpJacobianShim(const Derived & bc,
109  const unsigned int qp,
110  AssemblyDatum & datum) const
111  {
112  return bc.computeQpJacobian(qp, datum);
113  }
114  template <typename Derived>
115  KOKKOS_FUNCTION Real computeQpOffDiagJacobianShim(const Derived & bc,
116  const unsigned int jvar,
117  const unsigned int qp,
118  AssemblyDatum & datum) const
119  {
120  return bc.computeQpOffDiagJacobian(jvar, qp, datum);
121  }
123 
127  template <typename Derived>
129  KOKKOS_FUNCTION void operator()(ResidualLoop, const ThreadID tid, const Derived & bc) const;
130  template <typename Derived>
131  KOKKOS_FUNCTION void operator()(JacobianLoop, const ThreadID tid, const Derived & bc) const;
132  template <typename Derived>
133  KOKKOS_FUNCTION void
134  operator()(OffDiagJacobianLoop, const ThreadID tid, const Derived & bc) const;
136 
137 protected:
142 };
143 
144 template <typename Derived>
145 KOKKOS_FUNCTION void
146 NodalBC::operator()(ResidualLoop, const ThreadID tid, const Derived & bc) const
147 {
148  auto node = kokkosBoundaryNodeID(tid);
149  auto & sys = kokkosSystem(_kokkos_var.sys());
150 
151  if (!sys.isNodalDefined(node, _kokkos_var.var()))
152  return;
153 
155 
156  Real local_re = bc.computeQpResidualShim(bc, 0, datum);
157 
158  accumulateTaggedNodalResidual(false, local_re, node);
159 }
160 
161 template <typename Derived>
162 KOKKOS_FUNCTION void
163 NodalBC::operator()(JacobianLoop, const ThreadID tid, const Derived & bc) const
164 {
165  auto node = kokkosBoundaryNodeID(tid);
166  auto & sys = kokkosSystem(_kokkos_var.sys());
167 
168  if (!sys.isNodalDefined(node, _kokkos_var.var()))
169  return;
170 
172 
173  Real local_ke = bc.computeQpJacobianShim(bc, 0, datum);
174 
175  // This initializes the row to zero except the diagonal
176  accumulateTaggedNodalMatrix(false, local_ke, node, _kokkos_var.var());
177 }
178 
179 template <typename Derived>
180 KOKKOS_FUNCTION void
181 NodalBC::operator()(OffDiagJacobianLoop, const ThreadID tid, const Derived & bc) const
182 {
183  auto node = kokkosBoundaryNodeID(_thread(tid, 1));
184  auto & sys = kokkosSystem(_kokkos_var.sys());
185  auto jvar = sys.getCoupling(_kokkos_var.var())[_thread(tid, 0)];
186 
187  if (!sys.isNodalDefined(node, _kokkos_var.var()))
188  return;
189 
190  AssemblyDatum datum(node, kokkosAssembly(), kokkosSystems(), _kokkos_var, jvar);
191 
192  Real local_ke = bc.computeQpOffDiagJacobianShim(bc, jvar, 0, datum);
193 
194  accumulateTaggedNodalMatrix(true, local_ke, node, jvar);
195 }
196 
197 } // namespace Kokkos
198 } // namespace Moose
KOKKOS_FUNCTION Real computeQpResidualShim(const Derived &bc, const unsigned int qp, AssemblyDatum &datum) const
Shims for hook methods that can be leveraged to implement static polymorphism.
KOKKOS_FUNCTION unsigned int sys(unsigned int comp=0) const
Get the system number of a component.
KOKKOS_FUNCTION void operator()(ResidualLoop, const ThreadID tid, const Derived &bc) const
The parallel computation entry functions called by Kokkos.
KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int, const unsigned int, AssemblyDatum &) const
Compute off-diagonal Jacobian contribution on a node.
Definition: KokkosNodalBC.h:78
KOKKOS_FUNCTION ContiguousNodeID kokkosBoundaryNodeID(ThreadID tid) const
Get the contiguous node ID this Kokkos thread is operating on.
NodalBC(const InputParameters &parameters)
Constructor.
KOKKOS_FUNCTION const Assembly & kokkosAssembly() const
Get the const reference of the Kokkos assembly.
dof_id_type ThreadID
Definition: KokkosThread.h:18
static InputParameters validParams()
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
KOKKOS_FUNCTION void accumulateTaggedNodalResidual(const bool add, const Real local_re, const ContiguousNodeID node, const unsigned int comp=0) const
Accumulate or set local nodal residual contribution to tagged vectors.
KOKKOS_FUNCTION Real computeQpOffDiagJacobianShim(const Derived &bc, const unsigned int jvar, const unsigned int qp, AssemblyDatum &datum) const
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
Thread _thread
Kokkos thread object.
static auto defaultJacobian()
Get the function pointer of the default computeQpJacobian()
Definition: KokkosNodalBC.h:88
KOKKOS_FUNCTION Real computeQpJacobianShim(const Derived &bc, const unsigned int qp, AssemblyDatum &datum) const
virtual void computeResidual() override
Dispatch residual calculation.
KOKKOS_FUNCTION const System & kokkosSystem(unsigned int sys) const
Get the const reference of a Kokkos system.
Definition: KokkosSystem.h:625
The base class for Kokkos nodal boundary conditions.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const VariableValue _u
Current solution at nodes.
The Kokkos object that holds thread-private data in the parallel operations of Kokkos kernels...
Definition: KokkosDatum.h:244
Variable _kokkos_var
Kokkos variable.
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
virtual void computeJacobian() override
Dispatch diagonal and off-diagonal Jacobian calculation.
The Kokkos wrapper classes for MOOSE-like variable value access.
static auto defaultOffDiagJacobian()
Get the function pointer of the default computeQpOffDiagJacobian()
Definition: KokkosNodalBC.h:93
KOKKOS_FUNCTION unsigned int var(unsigned int comp=0) const
Get the variable number of a component.
KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int, AssemblyDatum &) const
Default methods to prevent compile errors even when these methods were not defined in the derived cla...
Definition: KokkosNodalBC.h:66
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
KOKKOS_FUNCTION void accumulateTaggedNodalMatrix(const bool add, const Real local_ke, const ContiguousNodeID node, const unsigned int jvar, const unsigned int comp=0) const
Accumulate or set local nodal Jacobian contribution to tagged matrices.
The base class for a user to derive their own Kokkos nodal boundary conditions.
Definition: KokkosNodalBC.h:36
KOKKOS_FUNCTION const Array< System > & kokkosSystems() const
Get the const reference of the Kokkos systems.
Definition: KokkosSystem.h:608