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 class KernelValue : public Kernel
41 {
42 public:
44 
49 
54 
63  KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int /* j */,
64  const unsigned int /* qp */,
65  AssemblyDatum & /* datum */) const
66  {
67  return 0;
68  }
75 
79  template <typename Derived>
81  KOKKOS_FUNCTION Real computeQpResidualShim(const Derived & kernel,
82  const unsigned int qp,
83  AssemblyDatum & datum) const
84  {
85  return kernel.computeQpResidual(qp, datum);
86  }
87  template <typename Derived>
88  KOKKOS_FUNCTION Real computeQpJacobianShim(const Derived & kernel,
89  const unsigned int j,
90  const unsigned int qp,
91  AssemblyDatum & datum) const
92  {
93  return kernel.computeQpJacobian(j, qp, datum);
94  }
96 
101  template <typename Derived>
103  KOKKOS_FUNCTION void computeResidualInternal(const Derived & kernel, AssemblyDatum & datum) const;
104  template <typename Derived>
105  KOKKOS_FUNCTION void computeJacobianInternal(const Derived & kernel, AssemblyDatum & datum) const;
107 };
108 
109 template <typename Derived>
110 KOKKOS_FUNCTION void
111 KernelValue::computeResidualInternal(const Derived & kernel, AssemblyDatum & datum) const
112 {
114  datum,
115  [&](Real * local_re, const unsigned int ib, const unsigned int ie)
116  {
117  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
118  {
119  datum.reinit();
120 
121  Real value = datum.JxW(qp) * kernel.computeQpResidualShim(kernel, qp, datum);
122 
123  for (unsigned int i = ib; i < ie; ++i)
124  local_re[i] += value * _test(datum, i, qp);
125  }
126  });
127 }
128 
129 template <typename Derived>
130 KOKKOS_FUNCTION void
131 KernelValue::computeJacobianInternal(const Derived & kernel, AssemblyDatum & datum) const
132 {
133  Real value = 0;
134 
136  datum,
137  [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
138  {
139  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
140  {
141  datum.reinit();
142 
143  unsigned int j_old = libMesh::invalid_uint;
144 
145  for (unsigned int ij = ijb; ij < ije; ++ij)
146  {
147  unsigned int i = ij % datum.n_jdofs();
148  unsigned int j = ij / datum.n_jdofs();
149 
150  if (j != j_old)
151  {
152  value = datum.JxW(qp) * kernel.computeQpJacobianShim(kernel, j, qp, datum);
153  j_old = j;
154  }
155 
156  local_ke[ij] += value * _test(datum, i, qp);
157  }
158  }
159  });
160 }
161 
162 } // namespace Kokkos
163 } // namespace Moose
static auto defaultJacobian()
Get the function pointer of the default computeQpJacobian()
KOKKOS_FUNCTION void computeJacobianInternal(AssemblyDatum &datum, function body) const
The common loop structure template for computing elemental Jacobian.
KOKKOS_FUNCTION void computeResidualInternal(const Derived &kernel, AssemblyDatum &datum) const
The parallel computation bodies that hide the base class methods to optimize for factoring out the te...
const unsigned int invalid_uint
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
const VariableTestValue _test
Current test function.
Definition: KokkosKernel.h:187
KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int, const unsigned int, AssemblyDatum &) const
Default methods to prevent compile errors even when these methods were not defined in the derived cla...
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...
KOKKOS_FUNCTION unsigned int n_jdofs() const
Get the number of local DOFs for the coupled variable.
Definition: KokkosDatum.h:314
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
KOKKOS_FUNCTION unsigned int n_qps() const
Get the number of local quadrature points.
Definition: KokkosDatum.h:100
KOKKOS_FUNCTION Real computeQpResidualShim(const Derived &kernel, const unsigned int qp, AssemblyDatum &datum) const
Shims for hook methods that can be leveraged to implement static polymorphism.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
KOKKOS_FUNCTION Real computeQpJacobianShim(const Derived &kernel, const unsigned int j, const unsigned int qp, AssemblyDatum &datum) const
KOKKOS_FUNCTION void computeJacobianInternal(const Derived &kernel, AssemblyDatum &datum) const
KernelValue(const InputParameters &parameters)
Constructor.
The Kokkos object that holds thread-private data in the parallel operations of Kokkos kernels...
Definition: KokkosDatum.h:244
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 base class for a user to derive their own Kokkos kernels.
Definition: KokkosKernel.h:39
KOKKOS_FUNCTION void computeResidualInternal(AssemblyDatum &datum, function body) const
The common loop structure template for computing elemental residual.
KOKKOS_FUNCTION void reinit()
Reset the reinit flag.
Definition: KokkosDatum.h:166
static InputParameters validParams()