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 
14 template <typename Derived>
16 {
18 
19 public:
21 
23 
24  KOKKOS_FUNCTION void computeJacobianInternal(const Derived * kernel, ResidualDatum & datum) const;
25 
26  KOKKOS_FUNCTION Real computeQpResidual(const unsigned int i,
27  const unsigned int qp,
28  ResidualDatum & datum) const;
29  KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int i,
30  const unsigned int j,
31  const unsigned int qp,
32  ResidualDatum & datum) const;
33 
34 protected:
35  const bool _lumping;
36 };
37 
38 template <typename Derived>
41 {
43  params.addParam<bool>("lumping", false, "True for mass matrix lumping, false otherwise");
44  return params;
45 }
46 
47 template <typename Derived>
49  : Moose::Kokkos::TimeKernel<Derived>(parameters),
50  _lumping(this->template getParam<bool>("lumping"))
51 {
52 }
53 
54 template <typename Derived>
55 KOKKOS_FUNCTION void
57  ResidualDatum & datum) const
58 {
60 
61  Real local_ke[MAX_CACHED_DOF];
62 
63  unsigned int num_batches = datum.n_idofs() * datum.n_jdofs() / MAX_CACHED_DOF;
64 
65  if ((datum.n_idofs() * datum.n_jdofs()) % MAX_CACHED_DOF)
66  ++num_batches;
67 
68  for (unsigned int batch = 0; batch < num_batches; ++batch)
69  {
70  unsigned int ijb = batch * MAX_CACHED_DOF;
71  unsigned int ije = ::Kokkos::min(ijb + MAX_CACHED_DOF, datum.n_idofs() * datum.n_jdofs());
72 
73  for (unsigned int ij = ijb; ij < ije; ++ij)
74  local_ke[ij - ijb] = 0;
75 
76  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
77  {
78  datum.reinit();
79 
80  for (unsigned int ij = ijb; ij < ije; ++ij)
81  {
82  unsigned int i = ij % datum.n_jdofs();
83  unsigned int j = ij / datum.n_jdofs();
84 
85  local_ke[ij - ijb] += datum.JxW(qp) * kernel->computeQpJacobian(i, j, qp, datum);
86  }
87  }
88 
89  for (unsigned int ij = ijb; ij < ije; ++ij)
90  {
91  unsigned int i = ij % datum.n_jdofs();
92  unsigned int j = ij / datum.n_jdofs();
93 
94  accumulateTaggedElementalMatrix(
95  local_ke[ij - ijb], datum.elem().id, i, _lumping ? i : j, datum.jvar());
96  }
97  }
98 }
99 
100 template <typename Derived>
101 KOKKOS_FUNCTION Real
103  const unsigned int qp,
104  ResidualDatum & datum) const
105 {
106  return _test(datum, i, qp) * _u_dot(datum, qp);
107 }
108 
109 template <typename Derived>
110 KOKKOS_FUNCTION Real
112  const unsigned int j,
113  const unsigned int qp,
114  ResidualDatum & datum) const
115 {
116  return _test(datum, i, qp) * _phi(datum, j, qp) * _du_dot_du;
117 }
118 
119 class KokkosTimeDerivativeKernel final : public KokkosTimeDerivative<KokkosTimeDerivativeKernel>
120 {
121 public:
122  static InputParameters validParams();
123 
125 };
KOKKOS_FUNCTION unsigned int n_jdofs() const
Get the number of local DOFs for the coupled variable.
Definition: KokkosDatum.h:287
usingKokkosTimeKernelMembers(Derived)
static InputParameters validParams()
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
KOKKOS_FUNCTION Real computeQpResidual(const unsigned int i, const unsigned int qp, ResidualDatum &datum) const
KokkosTimeDerivativeKernel(const InputParameters &parameters)
The base class for Kokkos time-derivative kernels.
static InputParameters validParams()
KOKKOS_FUNCTION unsigned int n_qps() const
Get the number of local quadrature points.
Definition: KokkosDatum.h:95
KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int i, const unsigned int j, const unsigned int qp, ResidualDatum &datum) const
KokkosTimeDerivative(const InputParameters &parameters)
constexpr unsigned int MAX_CACHED_DOF
Maximum number of DOFs to cache during residual and Jacobian computation.
Definition: KokkosHeader.h:79
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
All time kernels should inherit from this class.
Definition: TimeKernel.h:18
static InputParameters validParams()
KOKKOS_FUNCTION const ElementInfo & elem() const
Get the element information object.
Definition: KokkosDatum.h:80
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
KOKKOS_FUNCTION Real JxW(const unsigned int qp)
Get the transformed Jacobian weight.
Definition: KokkosDatum.h:126
KOKKOS_FUNCTION void computeJacobianInternal(const Derived *kernel, ResidualDatum &datum) const
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
auto min(const L &left, const R &right)
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
KOKKOS_FUNCTION unsigned int n_idofs() const
Get the number of local DOFs.
Definition: KokkosDatum.h:282
ContiguousElementID id
Contiguous element ID.
Definition: KokkosMesh.h:44