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 computeJacobianInternal(const Derived & kernel, ResidualDatum & datum) const;
23 
24  KOKKOS_FUNCTION Real computeQpResidual(const unsigned int i,
25  const unsigned int qp,
26  ResidualDatum & datum) const;
27  KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int i,
28  const unsigned int j,
29  const unsigned int qp,
30  ResidualDatum & datum) const;
31 
32 protected:
33  const bool _lumping;
34 };
35 
36 template <typename Derived>
37 KOKKOS_FUNCTION void
38 KokkosTimeDerivative::computeJacobianInternal(const Derived & kernel, ResidualDatum & datum) const
39 {
41 
42  Real local_ke[MAX_CACHED_DOF];
43 
44  unsigned int num_batches = datum.n_idofs() * datum.n_jdofs() / MAX_CACHED_DOF;
45 
46  if ((datum.n_idofs() * datum.n_jdofs()) % MAX_CACHED_DOF)
47  ++num_batches;
48 
49  for (unsigned int batch = 0; batch < num_batches; ++batch)
50  {
51  unsigned int ijb = batch * MAX_CACHED_DOF;
52  unsigned int ije = ::Kokkos::min(ijb + MAX_CACHED_DOF, datum.n_idofs() * datum.n_jdofs());
53 
54  for (unsigned int ij = ijb; ij < ije; ++ij)
55  local_ke[ij - ijb] = 0;
56 
57  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
58  {
59  datum.reinit();
60 
61  for (unsigned int ij = ijb; ij < ije; ++ij)
62  {
63  unsigned int i = ij % datum.n_jdofs();
64  unsigned int j = ij / datum.n_jdofs();
65 
66  local_ke[ij - ijb] += datum.JxW(qp) * kernel.computeQpJacobian(i, j, qp, datum);
67  }
68  }
69 
70  for (unsigned int ij = ijb; ij < ije; ++ij)
71  {
72  unsigned int i = ij % datum.n_jdofs();
73  unsigned int j = ij / datum.n_jdofs();
74 
76  local_ke[ij - ijb], datum.elem().id, i, _lumping ? i : j, datum.jvar());
77  }
78  }
79 }
80 
81 KOKKOS_FUNCTION inline Real
83  const unsigned int qp,
84  ResidualDatum & datum) const
85 {
86  return _test(datum, i, qp) * _u_dot(datum, qp);
87 }
88 
89 KOKKOS_FUNCTION inline Real
91  const unsigned int j,
92  const unsigned int qp,
93  ResidualDatum & datum) const
94 {
95  return _test(datum, i, qp) * _phi(datum, j, qp) * _du_dot_du;
96 }
KOKKOS_FUNCTION unsigned int n_jdofs() const
Get the number of local DOFs for the coupled variable.
Definition: KokkosDatum.h:315
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
const VariableTestValue _test
Current test function.
Definition: KokkosKernel.h:155
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
The base class for Kokkos time-derivative kernels.
const Scalar< const Real > _du_dot_du
Derivative of u_dot with respect to u.
KokkosTimeDerivative(const InputParameters &parameters)
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.
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
constexpr unsigned int MAX_CACHED_DOF
Maximum number of DOFs to cache during residual and Jacobian computation.
Definition: KokkosHeader.h:79
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:80
KOKKOS_FUNCTION Real JxW(const unsigned int qp)
Get the transformed Jacobian weight.
Definition: KokkosDatum.h:139
KOKKOS_FUNCTION void computeJacobianInternal(const Derived &kernel, ResidualDatum &datum) const
KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int i, const unsigned int j, const unsigned int qp, ResidualDatum &datum) const
const VariablePhiValue _phi
Current shape function.
Definition: KokkosKernel.h:163
auto min(const L &left, const R &right)
KOKKOS_FUNCTION unsigned int jvar() const
Get the coupled variable number.
Definition: KokkosDatum.h:330
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 unsigned int n_idofs() const
Get the number of local DOFs.
Definition: KokkosDatum.h:310
ContiguousElementID id
Contiguous element ID.
Definition: KokkosMesh.h:44