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 
52 template <typename Derived>
53 class Kernel : public KernelBase
54 {
55 public:
57 
62 
66  virtual void computeResidual() override;
70  virtual void computeJacobian() override;
71 
76 
85  KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int /* i */,
86  const unsigned int /* j */,
87  const unsigned int /* qp */,
88  ResidualDatum & /* datum */) const
89  {
90  return 0;
91  }
101  KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int /* i */,
102  const unsigned int /* j */,
103  const unsigned int /* jvar */,
104  const unsigned int /* qp */,
105  ResidualDatum & /* datum */) const
106  {
107  return 0;
108  }
110 
114  KOKKOS_FUNCTION void operator()(ResidualLoop, const ThreadID tid) const;
116  KOKKOS_FUNCTION void operator()(JacobianLoop, const ThreadID tid) const;
117  KOKKOS_FUNCTION void operator()(OffDiagJacobianLoop, const ThreadID tid) const;
119 
125 
131  KOKKOS_FUNCTION void computeResidualInternal(const Derived * kernel, ResidualDatum & datum) const;
137  KOKKOS_FUNCTION void computeJacobianInternal(const Derived * kernel, ResidualDatum & datum) const;
143  KOKKOS_FUNCTION void computeOffDiagJacobianInternal(const Derived * kernel,
144  ResidualDatum & datum) const;
146 
147 protected:
172 
177  virtual bool defaultJacobian() const
178  {
179  return &Derived::computeQpJacobian == &Kernel::computeQpJacobian;
180  }
185  virtual bool defaultOffDiagJacobian() const
186  {
187  return &Derived::computeQpOffDiagJacobian == &Kernel::computeQpOffDiagJacobian;
188  }
189 };
190 
191 template <typename Derived>
194 {
196  return params;
197 }
198 
199 template <typename Derived>
202  _test(),
203  _grad_test(),
204  _phi(),
205  _grad_phi(),
206  _u(_var),
207  _grad_u(_var)
208 {
210 }
211 
212 template <typename Derived>
213 void
215 {
216  ::Kokkos::RangePolicy<ResidualLoop, ExecSpace, ::Kokkos::IndexType<ThreadID>> policy(
217  0, numKokkosBlockElements());
218  ::Kokkos::parallel_for(policy, *static_cast<Derived *>(this));
219  ::Kokkos::fence();
220 }
221 
222 template <typename Derived>
223 void
225 {
226  if (!defaultJacobian())
227  {
228  ::Kokkos::RangePolicy<JacobianLoop, ExecSpace, ::Kokkos::IndexType<ThreadID>> policy(
229  0, numKokkosBlockElements());
230  ::Kokkos::parallel_for(policy, *static_cast<Derived *>(this));
231  ::Kokkos::fence();
232  }
233 
234  if (!defaultOffDiagJacobian())
235  {
236  auto & sys = kokkosSystem(_kokkos_var.sys());
237 
238  _thread.resize({sys.getCoupling(_kokkos_var.var()).size(), numKokkosBlockElements()});
239 
240  ::Kokkos::RangePolicy<OffDiagJacobianLoop, ExecSpace, ::Kokkos::IndexType<ThreadID>> policy(
241  0, _thread.size());
242  ::Kokkos::parallel_for(policy, *static_cast<Derived *>(this));
243  ::Kokkos::fence();
244  }
245 }
246 
247 template <typename Derived>
248 KOKKOS_FUNCTION void
250 {
251  auto kernel = static_cast<const Derived *>(this);
252  auto elem = kokkosBlockElementID(tid);
253 
254  ResidualDatum datum(elem, kokkosAssembly(), kokkosSystems(), _kokkos_var, _kokkos_var.var());
255 
256  kernel->computeResidualInternal(kernel, datum);
257 }
258 
259 template <typename Derived>
260 KOKKOS_FUNCTION void
262 {
263  auto kernel = static_cast<const Derived *>(this);
264  auto elem = kokkosBlockElementID(tid);
265 
266  ResidualDatum datum(elem, kokkosAssembly(), kokkosSystems(), _kokkos_var, _kokkos_var.var());
267 
268  kernel->computeJacobianInternal(kernel, datum);
269 }
270 
271 template <typename Derived>
272 KOKKOS_FUNCTION void
274 {
275  auto kernel = static_cast<const Derived *>(this);
276  auto elem = kokkosBlockElementID(_thread(tid, 1));
277 
278  auto & sys = kokkosSystem(_kokkos_var.sys());
279  auto jvar = sys.getCoupling(_kokkos_var.var())[_thread(tid, 0)];
280 
281  if (!sys.isVariableActive(jvar, kokkosMesh().getElementInfo(elem).subdomain))
282  return;
283 
284  ResidualDatum datum(elem, kokkosAssembly(), kokkosSystems(), _kokkos_var, jvar);
285 
286  kernel->computeOffDiagJacobianInternal(kernel, datum);
287 }
288 
289 template <typename Derived>
290 KOKKOS_FUNCTION void
291 Kernel<Derived>::computeResidualInternal(const Derived * kernel, ResidualDatum & datum) const
292 {
293  ResidualObject::computeResidualInternal(
294  datum,
295  [&](Real * local_re, const unsigned int ib, const unsigned int ie)
296  {
297  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
298  {
299  datum.reinit();
300 
301  for (unsigned int i = ib; i < ie; ++i)
302  local_re[i] += datum.JxW(qp) * kernel->computeQpResidual(i, qp, datum);
303  }
304  });
305 }
306 
307 template <typename Derived>
308 KOKKOS_FUNCTION void
309 Kernel<Derived>::computeJacobianInternal(const Derived * kernel, ResidualDatum & datum) const
310 {
311  ResidualObject::computeJacobianInternal(
312  datum,
313  [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
314  {
315  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
316  {
317  datum.reinit();
318 
319  for (unsigned int ij = ijb; ij < ije; ++ij)
320  {
321  unsigned int i = ij % datum.n_jdofs();
322  unsigned int j = ij / datum.n_jdofs();
323 
324  local_ke[ij] += datum.JxW(qp) * kernel->computeQpJacobian(i, j, qp, datum);
325  }
326  }
327  });
328 }
329 
330 template <typename Derived>
331 KOKKOS_FUNCTION void
333 {
334  ResidualObject::computeJacobianInternal(
335  datum,
336  [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
337  {
338  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
339  {
340  datum.reinit();
341 
342  for (unsigned int ij = ijb; ij < ije; ++ij)
343  {
344  unsigned int i = ij % datum.n_jdofs();
345  unsigned int j = ij / datum.n_jdofs();
346 
347  local_ke[ij] +=
348  datum.JxW(qp) * kernel->computeQpOffDiagJacobian(i, j, datum.jvar(), qp, datum);
349  }
350  }
351  });
352 }
353 
354 } // namespace Kokkos
355 } // namespace Moose
356 
357 #define usingKokkosKernelMembers(T) \
358  usingKokkosKernelBaseMembers; \
359  \
360 protected: \
361  using Moose::Kokkos::Kernel<T>::_test; \
362  using Moose::Kokkos::Kernel<T>::_grad_test; \
363  using Moose::Kokkos::Kernel<T>::_phi; \
364  using Moose::Kokkos::Kernel<T>::_grad_phi; \
365  using Moose::Kokkos::Kernel<T>::_u; \
366  using Moose::Kokkos::Kernel<T>::_grad_u; \
367  \
368 public: \
369  using Moose::Kokkos::Kernel<T>::operator()
VarFieldType
Definition: MooseTypes.h:722
KOKKOS_FUNCTION unsigned int n_jdofs() const
Get the number of local DOFs for the coupled variable.
Definition: KokkosDatum.h:287
virtual void computeResidual() override
Dispatch residual calculation.
Definition: KokkosKernel.h:214
KOKKOS_FUNCTION void computeOffDiagJacobianInternal(const Derived *kernel, ResidualDatum &datum) const
Compute off-diagonal Jacobian.
Definition: KokkosKernel.h:332
virtual bool defaultJacobian() const
Get whether computeQpJacobian() was not defined in the derived class.
Definition: KokkosKernel.h:177
virtual void computeResidual() override
Compute this Kernel&#39;s contribution to the residual.
Definition: Kernel.C:92
static InputParameters validParams()
Definition: KokkosKernel.h:193
MooseVariable & _var
This is a regular kernel so we cast to a regular MooseVariable.
Definition: Kernel.h:72
const VariablePhiValue _phi
Current shape function.
Definition: KokkosKernel.h:159
dof_id_type ThreadID
Definition: KokkosThread.h:18
KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int, const unsigned int, const unsigned int, ResidualDatum &) const
Default methods to prevent compile errors even when these methods were not defined in the derived cla...
Definition: KokkosKernel.h:85
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
static InputParameters validParams()
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
KOKKOS_FUNCTION void operator()(ResidualLoop, const ThreadID tid) const
The parallel computation entry functions called by Kokkos.
Definition: KokkosKernel.h:249
const VariableValue _u
Current solution at quadrature points.
Definition: KokkosKernel.h:167
const VariableGradient _grad_u
Gradient of the current solution at quadrature points.
Definition: KokkosKernel.h:171
The Kokkos wrapper classes for MOOSE-like shape function access.
virtual void computeJacobian() override
Dispatch diagonal and off-diagonal Jacobian calculation.
Definition: KokkosKernel.h:224
KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int, const unsigned int, const unsigned int, const unsigned int, ResidualDatum &) const
Compute off-diagonal Jacobian contribution on a quadrature point.
Definition: KokkosKernel.h:101
Kernel(const InputParameters &parameters)
Constructor.
Definition: KokkosKernel.h:200
KOKKOS_FUNCTION void computeResidualInternal(const Derived *kernel, ResidualDatum &datum) const
The parallel computation bodies that can be customized in the derived class by defining them in the d...
Definition: KokkosKernel.h:291
KOKKOS_FUNCTION void computeJacobianInternal(const Derived *kernel, ResidualDatum &datum) const
Compute diagonal Jacobian.
Definition: KokkosKernel.h:309
KOKKOS_FUNCTION unsigned int n_qps() const
Get the number of local quadrature points.
Definition: KokkosDatum.h:95
const VariablePhiGradient _grad_phi
Gradient of the current shape function.
Definition: KokkosKernel.h:163
The base class for Kokkos kernels.
virtual void computeJacobian() override
Compute this Kernel&#39;s contribution to the diagonal Jacobian entries.
Definition: Kernel.C:115
const VariableTestGradient _grad_test
Gradient of the current test function.
Definition: KokkosKernel.h:155
void addMooseVariableDependency(MooseVariableFieldBase *var)
Call this function to add the passed in MooseVariableFieldBase as a variable that this object depends...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
KOKKOS_FUNCTION Real JxW(const unsigned int qp)
Get the transformed Jacobian weight.
Definition: KokkosDatum.h:126
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
virtual bool defaultOffDiagJacobian() const
Get whether computeQpOffDiagJacobian() was not defined in the derived class.
Definition: KokkosKernel.h:185
The Kokkos wrapper classes for MOOSE-like variable value access.
The base class for a user to derive their own Kokkos kernels.
Definition: KokkosKernel.h:53
KOKKOS_FUNCTION unsigned int jvar() const
Get the coupled variable number.
Definition: KokkosDatum.h:302
The Kokkos object that holds thread-private data in the parallel operations of Kokkos residual object...
Definition: KokkosDatum.h:222
KOKKOS_FUNCTION void reinit()
Reset the reinit flag.
Definition: KokkosDatum.h:147
const VariableTestValue _test
Current test function.
Definition: KokkosKernel.h:151