https://mooseframework.inl.gov
KokkosKernelValue.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 
40 template <typename Derived>
41 class KernelValue : public Kernel<Derived>
42 {
43  usingKokkosKernelMembers(Derived);
44 
45 public:
47 
52 
57 
66  KOKKOS_FUNCTION Real precomputeQpJacobian(const unsigned int /* j */,
67  const unsigned int /* qp */,
68  ResidualDatum & /* datum */) const
69  {
70  return 0;
71  }
73 
78  KOKKOS_FUNCTION void computeResidualInternal(const Derived * kernel, ResidualDatum & datum) const;
80  KOKKOS_FUNCTION void computeJacobianInternal(const Derived * kernel, ResidualDatum & datum) const;
82 
83 protected:
88  virtual bool defaultJacobian() const override
89  {
90  return &Derived::precomputeQpJacobian == &KernelValue::precomputeQpJacobian;
91  }
92 };
93 
94 template <typename Derived>
97 {
99  return params;
100 }
101 
102 template <typename Derived>
103 KernelValue<Derived>::KernelValue(const InputParameters & parameters) : Kernel<Derived>(parameters)
104 {
105 }
106 
107 template <typename Derived>
108 KOKKOS_FUNCTION void
109 KernelValue<Derived>::computeResidualInternal(const Derived * kernel, ResidualDatum & datum) const
110 {
112  datum,
113  [&](Real * local_re, const unsigned int ib, const unsigned int ie)
114  {
115  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
116  {
117  datum.reinit();
118 
119  Real value = datum.JxW(qp) * kernel->precomputeQpResidual(qp, datum);
120 
121  for (unsigned int i = ib; i < ie; ++i)
122  local_re[i] += value * _test(datum, i, qp);
123  }
124  });
125 }
126 
127 template <typename Derived>
128 KOKKOS_FUNCTION void
129 KernelValue<Derived>::computeJacobianInternal(const Derived * kernel, ResidualDatum & datum) const
130 {
131  Real value = 0;
132 
134  datum,
135  [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
136  {
137  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
138  {
139  datum.reinit();
140 
141  unsigned int j_old = libMesh::invalid_uint;
142 
143  for (unsigned int ij = ijb; ij < ije; ++ij)
144  {
145  unsigned int i = ij % datum.n_jdofs();
146  unsigned int j = ij / datum.n_jdofs();
147 
148  if (j != j_old)
149  {
150  value = datum.JxW(qp) * kernel->precomputeQpJacobian(j, qp, datum);
151  j_old = j;
152  }
153 
154  local_ke[ij] += value * _test(datum, i, qp);
155  }
156  }
157  });
158 }
159 
160 } // namespace Kokkos
161 } // namespace Moose
162 
163 #define usingKokkosKernelValueMembers(T) usingKokkosKernelMembers(T)
KOKKOS_FUNCTION unsigned int n_jdofs() const
Get the number of local DOFs for the coupled variable.
Definition: KokkosDatum.h:287
static InputParameters validParams()
Definition: KokkosKernel.h:193
const unsigned int invalid_uint
KOKKOS_FUNCTION void computeJacobianInternal(ResidualDatum &datum, function body) const
The common loop structure template for computing elemental Jacobian.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
KOKKOS_FUNCTION void computeJacobianInternal(const Derived *kernel, ResidualDatum &datum) const
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
The base class for a user to derive their own Kokkos kernels where the residual is of the form...
KernelValue(const InputParameters &parameters)
Constructor.
KOKKOS_FUNCTION Real precomputeQpJacobian(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...
static InputParameters validParams()
virtual bool defaultJacobian() const override
Get whether precomputeQpJacobian() was not defined in the derived class.
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
KOKKOS_FUNCTION void computeResidualInternal(const Derived *kernel, ResidualDatum &datum) const
The parallel computation bodies that hide the base class methods to optimize for factoring out the te...
KOKKOS_FUNCTION unsigned int n_qps() const
Get the number of local quadrature points.
Definition: KokkosDatum.h:95
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:126
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
The base class for a user to derive their own Kokkos kernels.
Definition: KokkosKernel.h:53
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