https://mooseframework.inl.gov
KokkosTimeDerivative.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
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 "KokkosTimeKernel.h"
13 
15 {
16 public:
18 
20 
21  template <typename Derived>
22  KOKKOS_FUNCTION void computeResidualInternal(const Derived & kernel, AssemblyDatum & datum) const;
23  template <typename Derived>
24  KOKKOS_FUNCTION void computeJacobianInternal(const Derived & kernel, AssemblyDatum & datum) const;
25 
26  template <typename Derived>
27  KOKKOS_FUNCTION Real computeQpResidual(const unsigned int qp, AssemblyDatum & datum) const;
28  template <typename Derived>
29  KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int j,
30  const unsigned int qp,
31  AssemblyDatum & datum) const;
32  template <typename Derived>
33  KOKKOS_FUNCTION Real computeQpJacobianDummy(const unsigned int,
34  const unsigned int,
35  AssemblyDatum &) const
36  {
37  return 0;
38  }
39 
40  // The computeQpJacobian() of this class has the same signature with that of KernelValue (without
41  // test function index), but this class itself derives from TimeKernel whose base class is Kernel
42  // (with test function index). Therefore, their function pointers cannot be compared. Instead, we
43  // provide the function pointer of a dummy function to the dispatcher registry to make it think
44  // that non-default computeQpJacobian() was implemented.
45  template <typename Derived>
46  static auto defaultJacobian()
47  {
48  return &KokkosTimeDerivative::computeQpJacobianDummy<Derived>;
49  }
50 
51 protected:
52  const bool _lumping;
53 };
54 
55 template <typename Derived>
56 KOKKOS_FUNCTION void
57 KokkosTimeDerivative::computeResidualInternal(const Derived & kernel, AssemblyDatum & datum) const
58 {
59  ResidualObject::computeResidualInternal(
60  datum,
61  [&](Real * local_re, const unsigned int ib, const unsigned int ie)
62  {
63  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
64  {
65  Real value = datum.JxW(qp) * kernel.template computeQpResidual<Derived>(qp, datum);
66 
67  for (unsigned int i = ib; i < ie; ++i)
68  local_re[i] += value * _test(datum, i, qp);
69  }
70  });
71 }
72 
73 template <typename Derived>
74 KOKKOS_FUNCTION void
75 KokkosTimeDerivative::computeJacobianInternal(const Derived & kernel, AssemblyDatum & datum) const
76 {
78 
79  Real local_ke[MAX_CACHED_DOF];
80 
81  for (unsigned int j = datum.local_thread_id(); j < datum.n_jdofs();
82  j += datum.num_local_threads())
83  {
84  unsigned int num_batches = datum.n_idofs() / MAX_CACHED_DOF;
85 
86  if (datum.n_idofs() % MAX_CACHED_DOF)
87  ++num_batches;
88 
89  for (unsigned int batch = 0; batch < num_batches; ++batch)
90  {
91  unsigned int ib = batch * MAX_CACHED_DOF;
92  unsigned int ie = ::Kokkos::min(ib + MAX_CACHED_DOF, datum.n_idofs());
93 
94  for (unsigned int i = ib; i < ie; ++i)
95  local_ke[i - ib] = 0;
96 
97  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
98  {
99  Real value = datum.JxW(qp) * kernel.template computeQpJacobian<Derived>(j, qp, datum);
100 
101  for (unsigned int i = ib; i < ie; ++i)
102  local_ke[i - ib] += value * _test(datum, i, qp);
103  }
104 
105  for (unsigned int i = ib; i < ie; ++i)
107  local_ke[i - ib], datum.elem().id, i, _lumping ? i : j, datum.jvar());
108  }
109  }
110 }
111 
112 template <typename Derived>
113 KOKKOS_FUNCTION Real
114 KokkosTimeDerivative::computeQpResidual(const unsigned int qp, AssemblyDatum & datum) const
115 {
116  return _u_dot(datum, qp);
117 }
118 
119 template <typename Derived>
120 KOKKOS_FUNCTION Real
122  const unsigned int qp,
123  AssemblyDatum & datum) const
124 {
125  return _phi(datum, j, qp) * _du_dot_du;
126 }
KOKKOS_FUNCTION Real computeQpResidual(const unsigned int qp, AssemblyDatum &datum) const
KOKKOS_FUNCTION unsigned int num_local_threads() const
Get the number of local threads.
Definition: KokkosDatum.h:202
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...
KOKKOS_FUNCTION unsigned int jvar() const
Get the coupled variable number.
Definition: KokkosDatum.h:456
KOKKOS_FUNCTION unsigned int n_jdofs() const
Get the number of local DOFs for the coupled variable.
Definition: KokkosDatum.h:436
The base class for Kokkos time-derivative kernels.
KOKKOS_FUNCTION Real computeQpJacobianDummy(const unsigned int, const unsigned int, AssemblyDatum &) const
const Scalar< const Real > _du_dot_du
Derivative of u_dot with respect to u.
KokkosTimeDerivative(const InputParameters &parameters)
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
KOKKOS_FUNCTION void accumulateTaggedElementalMatrix(const Real local_ke, const ContiguousElementID elem, const unsigned int i, const unsigned int j, const unsigned int jvar, const unsigned int comp=0) const
Accumulate local elemental Jacobian contribution to tagged matrices.
KOKKOS_FUNCTION unsigned int n_idofs() const
Get the number of local DOFs.
Definition: KokkosDatum.h:431
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:118
KOKKOS_FUNCTION unsigned int local_thread_id() const
Get the current local thread ID.
Definition: KokkosDatum.h:197
KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int j, const unsigned int qp, AssemblyDatum &datum) const
KOKKOS_FUNCTION void computeJacobianInternal(const Derived &kernel, AssemblyDatum &datum) const
constexpr unsigned int MAX_CACHED_DOF
Maximum number of DOFs to cache during residual and Jacobian computation.
Definition: KokkosHeader.h:77
static InputParameters validParams()
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
KOKKOS_FUNCTION const ElementInfo & elem() const
Get the element information object.
Definition: KokkosDatum.h:84
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
const VariablePhiValue _phi
Current shape function.
Definition: KokkosKernel.h:181
auto min(const L &left, const R &right)
ContiguousElementID id
Contiguous element ID.
Definition: KokkosMesh.h:42
KOKKOS_FUNCTION void computeResidualInternal(const Derived &kernel, AssemblyDatum &datum) const