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::Kokkos
15 {
16 
39 class Kernel : public KernelBase
40 {
41 public:
43 
48 
52  virtual void computeResidual() override;
56  virtual void computeJacobian() override;
57 
62 
72  template <typename Derived>
73  KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int /* i */,
74  const unsigned int /* j */,
75  const unsigned int /* qp */,
76  AssemblyDatum & /* datum */) const
77  {
78  ::Kokkos::abort("Default computeQpJacobian() should never be called. Make sure you properly "
79  "redefined this method in your class without typos.");
80 
81  return 0;
82  }
93  template <typename Derived>
94  KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int /* i */,
95  const unsigned int /* j */,
96  const unsigned int /* jvar */,
97  const unsigned int /* qp */,
98  AssemblyDatum & /* datum */) const
99  {
100  ::Kokkos::abort(
101  "Default computeQpOffDiagJacobian() should never be called. Make sure you properly "
102  "redefined this method in your class without typos.");
103 
104  return 0;
105  }
107 
113  template <typename Derived>
115  static auto defaultJacobian()
116  {
117  return &Kernel::computeQpJacobian<Derived>;
118  }
119  template <typename Derived>
121  {
122  return &Kernel::computeQpOffDiagJacobian<Derived>;
123  }
125 
129  template <typename Derived>
131  KOKKOS_FUNCTION void operator()(ResidualLoop, const ThreadID tid, const Derived & kernel) const;
132  template <typename Derived>
133  KOKKOS_FUNCTION void operator()(JacobianLoop, const ThreadID tid, const Derived & kernel) const;
134  template <typename Derived>
135  KOKKOS_FUNCTION void
136  operator()(OffDiagJacobianLoop, const ThreadID tid, const Derived & kernel) const;
138 
144 
150  template <typename Derived>
151  KOKKOS_FUNCTION void computeResidualInternal(const Derived & kernel, AssemblyDatum & datum) const;
157  template <typename Derived>
158  KOKKOS_FUNCTION void computeJacobianInternal(const Derived & kernel, AssemblyDatum & datum) const;
164  template <typename Derived>
165  KOKKOS_FUNCTION void computeOffDiagJacobianInternal(const Derived & kernel,
166  AssemblyDatum & datum) const;
168 
169 protected:
194 };
195 
196 template <typename Derived>
197 KOKKOS_FUNCTION void
198 Kernel::operator()(ResidualLoop, const ThreadID tid, const Derived & kernel) const
199 {
200  auto elem = kokkosBlockElementID(_thread(tid, 1));
201 
202  AssemblyDatum datum(elem,
204  kokkosAssembly(),
205  kokkosSystems(),
206  _kokkos_var,
207  _kokkos_var.var());
208 
209  datum.set_local_parallel(_thread(tid, 0), _thread.size(0));
210 
211  kernel.computeResidualInternal(kernel, datum);
212 }
213 
214 template <typename Derived>
215 KOKKOS_FUNCTION void
216 Kernel::operator()(JacobianLoop, const ThreadID tid, const Derived & kernel) const
217 {
218  auto elem = kokkosBlockElementID(_thread(tid, 1));
219 
220  AssemblyDatum datum(elem,
222  kokkosAssembly(),
223  kokkosSystems(),
224  _kokkos_var,
225  _kokkos_var.var());
226 
227  datum.set_local_parallel(_thread(tid, 0), _thread.size(0));
228 
229  kernel.computeJacobianInternal(kernel, datum);
230 }
231 
232 template <typename Derived>
233 KOKKOS_FUNCTION void
234 Kernel::operator()(OffDiagJacobianLoop, const ThreadID tid, const Derived & kernel) const
235 {
236  auto elem = kokkosBlockElementID(_thread(tid, 2));
237 
238  auto & sys = kokkosSystem(_kokkos_var.sys());
239  auto jvar = sys.getCoupling(_kokkos_var.var())[_thread(tid, 1)];
240 
241  if (!sys.isVariableActive(jvar, kokkosMesh().getElementInfo(elem).subdomain))
242  return;
243 
244  AssemblyDatum datum(
246 
247  datum.set_local_parallel(_thread(tid, 0), _thread.size(0));
248 
249  kernel.computeOffDiagJacobianInternal(kernel, datum);
250 }
251 
252 template <typename Derived>
253 KOKKOS_FUNCTION void
254 Kernel::computeResidualInternal(const Derived & kernel, AssemblyDatum & datum) const
255 {
257  datum,
258  [&](Real * local_re, const unsigned int ib, const unsigned int ie)
259  {
260  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
261  for (unsigned int i = ib; i < ie; ++i)
262  local_re[i] += datum.JxW(qp) * kernel.template computeQpResidual<Derived>(i, qp, datum);
263  });
264 }
265 
266 template <typename Derived>
267 KOKKOS_FUNCTION void
268 Kernel::computeJacobianInternal(const Derived & kernel, AssemblyDatum & datum) const
269 {
271  datum,
272  [&](Real * local_ke, const unsigned int ib, const unsigned int ie, const unsigned int j)
273  {
274  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
275  for (unsigned int i = ib; i < ie; ++i)
276  local_ke[i] +=
277  datum.JxW(qp) * kernel.template computeQpJacobian<Derived>(i, j, qp, datum);
278  });
279 }
280 
281 template <typename Derived>
282 KOKKOS_FUNCTION void
283 Kernel::computeOffDiagJacobianInternal(const Derived & kernel, AssemblyDatum & datum) const
284 {
286  datum,
287  [&](Real * local_ke, const unsigned int ib, const unsigned int ie, const unsigned int j)
288  {
289  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
290  for (unsigned int i = ib; i < ie; ++i)
291  local_ke[i] += datum.JxW(qp) * kernel.template computeQpOffDiagJacobian<Derived>(
292  i, j, datum.jvar(), qp, datum);
293  });
294 }
295 
296 } // namespace Moose::Kokkos
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.
KOKKOS_FUNCTION void set_local_parallel(const unsigned int local_thread_id, const unsigned int num_local_threads)
Set local parallelization option.
Definition: KokkosDatum.h:187
const unsigned int invalid_uint
KOKKOS_FUNCTION const Assembly & kokkosAssembly() const
Get the const reference of the Kokkos assembly.
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:94
const VariableValue _u
Current solution at quadrature points.
Definition: KokkosKernel.h:189
static auto defaultJacobian()
Functions used to check if users have overriden the hook methods, whose calculations can be skipped w...
Definition: KokkosKernel.h:115
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
const VariableTestValue _test
Current test function.
Definition: KokkosKernel.h:173
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:456
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:254
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:73
KOKKOS_FUNCTION void computeJacobianInternal(const Derived &kernel, AssemblyDatum &datum) const
Compute diagonal Jacobian.
Definition: KokkosKernel.h:268
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:198
KOKKOS_FUNCTION unsigned int n_qps() const
Get the number of local quadrature points.
Definition: KokkosDatum.h:118
KOKKOS_FUNCTION const System & kokkosSystem(unsigned int sys) const
Get the const reference of a Kokkos system.
Definition: KokkosSystem.h:792
The base class for Kokkos kernels.
const VariableGradient _grad_u
Gradient of the current solution at quadrature points.
Definition: KokkosKernel.h:193
The Kokkos wrapper classes for MOOSE-like variable value access.
virtual void computeJacobian() override
Dispatch diagonal and off-diagonal Jacobian calculation.
KOKKOS_FUNCTION ContiguousElementID kokkosBlockElementID(Moose::Kokkos::ThreadID tid) const
Get the contiguous element ID this Kokkos thread is operating on.
KOKKOS_FUNCTION const Mesh & kokkosMesh() const
Get the const reference of the Kokkos mesh.
Definition: KokkosMesh.h:452
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const VariableTestGradient _grad_test
Gradient of the current test function.
Definition: KokkosKernel.h:177
Kernel(const InputParameters &parameters)
Constructor.
KOKKOS_FUNCTION void computeOffDiagJacobianInternal(const Derived &kernel, AssemblyDatum &datum) const
Compute off-diagonal Jacobian.
Definition: KokkosKernel.h:283
const VariablePhiGradient _grad_phi
Gradient of the current shape function.
Definition: KokkosKernel.h:185
The Kokkos object that holds thread-private data in the parallel operations of Kokkos kernels...
Definition: KokkosDatum.h:364
KOKKOS_FUNCTION Real JxW(const unsigned int qp)
Get the transformed Jacobian weight.
Definition: KokkosDatum.h:311
Variable _kokkos_var
Kokkos variable.
static auto defaultOffDiagJacobian()
Definition: KokkosKernel.h:120
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:181
The base class for a user to derive their own Kokkos kernels.
Definition: KokkosKernel.h:39
KOKKOS_FUNCTION thread_id_type size() const
Get the total thread pool size.
Definition: KokkosThread.h:50
KOKKOS_FUNCTION const auto & getCoupling(unsigned int var) const
Get the list of off-diagonal coupled variable numbers of a variable.
Definition: KokkosSystem.h:146
KOKKOS_FUNCTION void computeResidualInternal(AssemblyDatum &datum, function body) const
The common loop structure template for computing elemental residual.
KOKKOS_FUNCTION const Array< System > & kokkosSystems() const
Get the const reference of the Kokkos systems.
Definition: KokkosSystem.h:775
dof_id_type ThreadID
Definition: KokkosThread.h:22