https://mooseframework.inl.gov
KokkosKernel.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 "KokkosKernelBase.h"
13 
14 namespace Moose
15 {
16 namespace Kokkos
17 {
18 
39 class Kernel : public KernelBase
40 {
41 public:
43 
48 
52  virtual void computeResidual() override;
56  virtual void computeJacobian() override;
57 
62 
71  KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int /* i */,
72  const unsigned int /* j */,
73  const unsigned int /* qp */,
74  AssemblyDatum & /* datum */) const
75  {
76  return 0;
77  }
87  KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int /* i */,
88  const unsigned int /* j */,
89  const unsigned int /* jvar */,
90  const unsigned int /* qp */,
91  AssemblyDatum & /* datum */) const
92  {
93  return 0;
94  }
99  static auto defaultJacobian() { return &Kernel::computeQpJacobian; }
106 
110  template <typename Derived>
112  KOKKOS_FUNCTION Real computeQpResidualShim(const Derived & kernel,
113  const unsigned int i,
114  const unsigned int qp,
115  AssemblyDatum & datum) const
116  {
117  return kernel.computeQpResidual(i, qp, datum);
118  }
119  template <typename Derived>
120  KOKKOS_FUNCTION Real computeQpJacobianShim(const Derived & kernel,
121  const unsigned int i,
122  const unsigned int j,
123  const unsigned int qp,
124  AssemblyDatum & datum) const
125  {
126  return kernel.computeQpJacobian(i, j, qp, datum);
127  }
128  template <typename Derived>
129  KOKKOS_FUNCTION Real computeQpOffDiagJacobianShim(const Derived & kernel,
130  const unsigned int i,
131  const unsigned int j,
132  const unsigned int jvar,
133  const unsigned int qp,
134  AssemblyDatum & datum) const
135  {
136  return kernel.computeQpOffDiagJacobian(i, j, jvar, qp, datum);
137  }
139 
143  template <typename Derived>
145  KOKKOS_FUNCTION void operator()(ResidualLoop, const ThreadID tid, const Derived & kernel) const;
146  template <typename Derived>
147  KOKKOS_FUNCTION void operator()(JacobianLoop, const ThreadID tid, const Derived & kernel) const;
148  template <typename Derived>
149  KOKKOS_FUNCTION void
150  operator()(OffDiagJacobianLoop, const ThreadID tid, const Derived & kernel) const;
152 
158 
164  template <typename Derived>
165  KOKKOS_FUNCTION void computeResidualInternal(const Derived & kernel, AssemblyDatum & datum) const;
171  template <typename Derived>
172  KOKKOS_FUNCTION void computeJacobianInternal(const Derived & kernel, AssemblyDatum & datum) const;
178  template <typename Derived>
179  KOKKOS_FUNCTION void computeOffDiagJacobianInternal(const Derived & kernel,
180  AssemblyDatum & datum) const;
182 
183 protected:
208 };
209 
210 template <typename Derived>
211 KOKKOS_FUNCTION void
212 Kernel::operator()(ResidualLoop, const ThreadID tid, const Derived & kernel) const
213 {
214  auto elem = kokkosBlockElementID(tid);
215 
216  AssemblyDatum datum(elem,
218  kokkosAssembly(),
219  kokkosSystems(),
220  _kokkos_var,
221  _kokkos_var.var());
222 
223  kernel.computeResidualInternal(kernel, datum);
224 }
225 
226 template <typename Derived>
227 KOKKOS_FUNCTION void
228 Kernel::operator()(JacobianLoop, const ThreadID tid, const Derived & kernel) const
229 {
230  auto elem = kokkosBlockElementID(tid);
231 
232  AssemblyDatum datum(elem,
234  kokkosAssembly(),
235  kokkosSystems(),
236  _kokkos_var,
237  _kokkos_var.var());
238 
239  kernel.computeJacobianInternal(kernel, datum);
240 }
241 
242 template <typename Derived>
243 KOKKOS_FUNCTION void
244 Kernel::operator()(OffDiagJacobianLoop, const ThreadID tid, const Derived & kernel) const
245 {
246  auto elem = kokkosBlockElementID(_thread(tid, 1));
247 
248  auto & sys = kokkosSystem(_kokkos_var.sys());
249  auto jvar = sys.getCoupling(_kokkos_var.var())[_thread(tid, 0)];
250 
251  if (!sys.isVariableActive(jvar, kokkosMesh().getElementInfo(elem).subdomain))
252  return;
253 
254  AssemblyDatum datum(
256 
257  kernel.computeOffDiagJacobianInternal(kernel, datum);
258 }
259 
260 template <typename Derived>
261 KOKKOS_FUNCTION void
262 Kernel::computeResidualInternal(const Derived & kernel, AssemblyDatum & datum) const
263 {
265  datum,
266  [&](Real * local_re, const unsigned int ib, const unsigned int ie)
267  {
268  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
269  {
270  datum.reinit();
271 
272  for (unsigned int i = ib; i < ie; ++i)
273  local_re[i] += datum.JxW(qp) * kernel.computeQpResidualShim(kernel, i, qp, datum);
274  }
275  });
276 }
277 
278 template <typename Derived>
279 KOKKOS_FUNCTION void
280 Kernel::computeJacobianInternal(const Derived & kernel, AssemblyDatum & datum) const
281 {
283  datum,
284  [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
285  {
286  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
287  {
288  datum.reinit();
289 
290  for (unsigned int ij = ijb; ij < ije; ++ij)
291  {
292  unsigned int i = ij % datum.n_jdofs();
293  unsigned int j = ij / datum.n_jdofs();
294 
295  local_ke[ij] += datum.JxW(qp) * kernel.computeQpJacobianShim(kernel, i, j, qp, datum);
296  }
297  }
298  });
299 }
300 
301 template <typename Derived>
302 KOKKOS_FUNCTION void
303 Kernel::computeOffDiagJacobianInternal(const Derived & kernel, AssemblyDatum & datum) const
304 {
306  datum,
307  [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
308  {
309  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
310  {
311  datum.reinit();
312 
313  for (unsigned int ij = ijb; ij < ije; ++ij)
314  {
315  unsigned int i = ij % datum.n_jdofs();
316  unsigned int j = ij / datum.n_jdofs();
317 
318  local_ke[ij] += datum.JxW(qp) * kernel.computeQpOffDiagJacobianShim(
319  kernel, i, j, datum.jvar(), qp, datum);
320  }
321  }
322  });
323 }
324 
325 } // namespace Kokkos
326 } // namespace Moose
KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int, const unsigned int, const unsigned int, const unsigned int, AssemblyDatum &) const
Compute off-diagonal Jacobian contribution on a quadrature point.
Definition: KokkosKernel.h:87
KOKKOS_FUNCTION unsigned int sys(unsigned int comp=0) const
Get the system number of a component.
KOKKOS_FUNCTION void computeJacobianInternal(AssemblyDatum &datum, function body) const
The common loop structure template for computing elemental Jacobian.
virtual void computeResidual() override
Dispatch residual calculation.
const unsigned int invalid_uint
KOKKOS_FUNCTION const Assembly & kokkosAssembly() const
Get the const reference of the Kokkos assembly.
KOKKOS_FUNCTION Real computeQpOffDiagJacobianShim(const Derived &kernel, const unsigned int i, const unsigned int j, const unsigned int jvar, const unsigned int qp, AssemblyDatum &datum) const
Definition: KokkosKernel.h:129
dof_id_type ThreadID
Definition: KokkosThread.h:18
static auto defaultJacobian()
Get the function pointer of the default computeQpJacobian()
Definition: KokkosKernel.h:99
const VariableValue _u
Current solution at quadrature points.
Definition: KokkosKernel.h:203
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
const VariableTestValue _test
Current test function.
Definition: KokkosKernel.h:187
KOKKOS_FUNCTION Real computeQpJacobianShim(const Derived &kernel, const unsigned int i, const unsigned int j, const unsigned int qp, AssemblyDatum &datum) const
Definition: KokkosKernel.h:120
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
Thread _thread
Kokkos thread object.
KOKKOS_FUNCTION unsigned int jvar() const
Get the coupled variable number.
Definition: KokkosDatum.h:329
KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int, const unsigned int, const unsigned int, AssemblyDatum &) const
Default methods to prevent compile errors even when these methods were not defined in the derived cla...
Definition: KokkosKernel.h:71
KOKKOS_FUNCTION unsigned int n_jdofs() const
Get the number of local DOFs for the coupled variable.
Definition: KokkosDatum.h:314
The Kokkos wrapper classes for MOOSE-like shape function access.
KOKKOS_FUNCTION void computeResidualInternal(const Derived &kernel, AssemblyDatum &datum) const
The parallel computation bodies that can be customized in the derived class by defining them in the d...
Definition: KokkosKernel.h:262
static auto defaultOffDiagJacobian()
Get the function pointer of the default computeQpOffDiagJacobian()
Definition: KokkosKernel.h:104
KOKKOS_FUNCTION void computeJacobianInternal(const Derived &kernel, AssemblyDatum &datum) const
Compute diagonal Jacobian.
Definition: KokkosKernel.h:280
KOKKOS_FUNCTION ContiguousElementID kokkosBlockElementID(ThreadID tid) const
Get the contiguous element ID this Kokkos thread is operating on.
static InputParameters validParams()
KOKKOS_FUNCTION void operator()(ResidualLoop, const ThreadID tid, const Derived &kernel) const
The parallel computation entry functions called by Kokkos.
Definition: KokkosKernel.h:212
KOKKOS_FUNCTION unsigned int n_qps() const
Get the number of local quadrature points.
Definition: KokkosDatum.h:100
KOKKOS_FUNCTION Real computeQpResidualShim(const Derived &kernel, const unsigned int i, const unsigned int qp, AssemblyDatum &datum) const
Shims for hook methods that can be leveraged to implement static polymorphism.
Definition: KokkosKernel.h:112
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 kernels.
const VariableGradient _grad_u
Gradient of the current solution at quadrature points.
Definition: KokkosKernel.h:207
virtual void computeJacobian() override
Dispatch diagonal and off-diagonal Jacobian calculation.
KOKKOS_FUNCTION const Mesh & kokkosMesh() const
Get the const reference of the Kokkos mesh.
Definition: KokkosMesh.h:360
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const VariableTestGradient _grad_test
Gradient of the current test function.
Definition: KokkosKernel.h:191
Kernel(const InputParameters &parameters)
Constructor.
KOKKOS_FUNCTION void computeOffDiagJacobianInternal(const Derived &kernel, AssemblyDatum &datum) const
Compute off-diagonal Jacobian.
Definition: KokkosKernel.h:303
const VariablePhiGradient _grad_phi
Gradient of the current shape function.
Definition: KokkosKernel.h:199
The Kokkos object that holds thread-private data in the parallel operations of Kokkos kernels...
Definition: KokkosDatum.h:244
KOKKOS_FUNCTION Real JxW(const unsigned int qp)
Get the transformed Jacobian weight.
Definition: KokkosDatum.h:139
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...
The Kokkos wrapper classes for MOOSE-like variable value access.
KOKKOS_FUNCTION unsigned int var(unsigned int comp=0) const
Get the variable number of a component.
const VariablePhiValue _phi
Current shape function.
Definition: KokkosKernel.h:195
The base class for a user to derive their own Kokkos kernels.
Definition: KokkosKernel.h:39
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 computeResidualInternal(AssemblyDatum &datum, function body) const
The common loop structure template for computing elemental residual.
KOKKOS_FUNCTION void reinit()
Reset the reinit flag.
Definition: KokkosDatum.h:166
KOKKOS_FUNCTION const Array< System > & kokkosSystems() const
Get the const reference of the Kokkos systems.
Definition: KokkosSystem.h:608