https://mooseframework.inl.gov
KokkosTimeKernel.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 "KokkosKernel.h"
13 
14 namespace Moose
15 {
16 namespace Kokkos
17 {
18 
22 class TimeKernel : public Kernel
23 {
24 public:
26 
31 
39  KOKKOS_FUNCTION void computeResidualAdditional(const unsigned int /* ib */,
40  const unsigned int /* ie */,
41  ResidualDatum & /* datum */,
42  Real * /* local_re */) const
43  {
44  }
45 
50  template <typename Derived>
51  KOKKOS_FUNCTION void computeResidualInternal(const Derived & kernel, ResidualDatum & datum) const;
52 
53 protected:
62 };
63 
64 template <typename Derived>
65 KOKKOS_FUNCTION void
66 TimeKernel::computeResidualInternal(const Derived & kernel, ResidualDatum & datum) const
67 {
69  datum,
70  [&](Real * local_re, const unsigned int ib, const unsigned int ie)
71  {
72  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
73  {
74  datum.reinit();
75 
76  for (unsigned int i = ib; i < ie; ++i)
77  local_re[i] += datum.JxW(qp) * kernel.computeQpResidual(i, qp, datum);
78  }
79 
80  kernel.computeResidualAdditional(ib, ie, datum, local_re);
81  });
82 }
83 
84 } // namespace Kokkos
85 } // namespace Moose
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...
The base class for Kokkos time-derivative kernels.
const Scalar< const Real > _du_dot_du
Derivative of u_dot with respect to u.
TimeKernel(const InputParameters &parameters)
Constructor.
KOKKOS_FUNCTION void computeResidualAdditional(const unsigned int, const unsigned int, ResidualDatum &, Real *) const
Hook for additional computation for residual after the standard calls.
const VariableValue _u_dot
Time derivative of the current solution at quadrature points.
KOKKOS_FUNCTION unsigned int n_qps() const
Get the number of local quadrature points.
Definition: KokkosDatum.h:100
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
KOKKOS_FUNCTION void computeResidualInternal(ResidualDatum &datum, function body) const
The common loop structure template for computing elemental residual.
KOKKOS_FUNCTION Real JxW(const unsigned int qp)
Get the transformed Jacobian weight.
Definition: KokkosDatum.h:139
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.
The base class for a user to derive their own Kokkos kernels.
Definition: KokkosKernel.h:40
The Kokkos object that holds thread-private data in the parallel operations of Kokkos residual object...
Definition: KokkosDatum.h:245
KOKKOS_FUNCTION void reinit()
Reset the reinit flag.
Definition: KokkosDatum.h:166
KOKKOS_FUNCTION void computeResidualInternal(const Derived &kernel, ResidualDatum &datum) const
The parallel computation body that hides the base class method to allow additional computation for re...