https://mooseframework.inl.gov
KokkosKernelGrad.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 
41 class KernelGrad : public Kernel
42 {
43 public:
45 
50 
55 
64  KOKKOS_FUNCTION Real3 computeQpJacobian(const unsigned int /* j */,
65  const unsigned int /* qp */,
66  AssemblyDatum & /* datum */) const
67  {
68  return Real3(0);
69  }
76 
80  template <typename Derived>
82  KOKKOS_FUNCTION Real3 computeQpResidualShim(const Derived & kernel,
83  const unsigned int qp,
84  AssemblyDatum & datum) const
85  {
86  return kernel.computeQpResidual(qp, datum);
87  }
88  template <typename Derived>
89  KOKKOS_FUNCTION Real3 computeQpJacobianShim(const Derived & kernel,
90  const unsigned int j,
91  const unsigned int qp,
92  AssemblyDatum & datum) const
93  {
94  return kernel.computeQpJacobian(j, qp, datum);
95  }
97 
102  template <typename Derived>
104  KOKKOS_FUNCTION void computeResidualInternal(const Derived & kernel, AssemblyDatum & datum) const;
105  template <typename Derived>
106  KOKKOS_FUNCTION void computeJacobianInternal(const Derived & kernel, AssemblyDatum & datum) const;
108 };
109 
110 template <typename Derived>
111 KOKKOS_FUNCTION void
112 KernelGrad::computeResidualInternal(const Derived & kernel, AssemblyDatum & datum) const
113 {
115  datum,
116  [&](Real * local_re, const unsigned int ib, const unsigned int ie)
117  {
118  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
119  {
120  datum.reinit();
121 
122  Real3 value = datum.JxW(qp) * kernel.computeQpResidualShim(kernel, qp, datum);
123 
124  for (unsigned int i = ib; i < ie; ++i)
125  local_re[i] += value * _grad_test(datum, i, qp);
126  }
127  });
128 }
129 
130 template <typename Derived>
131 KOKKOS_FUNCTION void
132 KernelGrad::computeJacobianInternal(const Derived & kernel, AssemblyDatum & datum) const
133 {
134  Real3 value;
135 
137  datum,
138  [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
139  {
140  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
141  {
142  datum.reinit();
143 
144  unsigned int j_old = libMesh::invalid_uint;
145 
146  for (unsigned int ij = ijb; ij < ije; ++ij)
147  {
148  unsigned int i = ij % datum.n_jdofs();
149  unsigned int j = ij / datum.n_jdofs();
150 
151  if (j != j_old)
152  {
153  value = datum.JxW(qp) * kernel.computeQpJacobianShim(kernel, j, qp, datum);
154  j_old = j;
155  }
156 
157  local_ke[ij] += value * _grad_test(datum, i, qp);
158  }
159  }
160  });
161 }
162 
163 } // namespace Kokkos
164 } // namespace Moose
KOKKOS_FUNCTION void computeJacobianInternal(AssemblyDatum &datum, function body) const
The common loop structure template for computing elemental Jacobian.
const unsigned int invalid_uint
KOKKOS_FUNCTION void computeJacobianInternal(const Derived &kernel, AssemblyDatum &datum) const
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 unsigned int n_jdofs() const
Get the number of local DOFs for the coupled variable.
Definition: KokkosDatum.h:314
KernelGrad(const InputParameters &parameters)
Constructor.
static auto defaultJacobian()
Get the function pointer of the default computeQpJacobian()
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 Real3 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...
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 gr...
KOKKOS_FUNCTION Real3 computeQpResidualShim(const Derived &kernel, const unsigned int qp, AssemblyDatum &datum) const
Shims for hook methods that can be leveraged to implement static polymorphism.
static InputParameters validParams()
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const VariableTestGradient _grad_test
Gradient of the current test function.
Definition: KokkosKernel.h:191
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 Real3 computeQpJacobianShim(const Derived &kernel, const unsigned int j, const unsigned int qp, AssemblyDatum &datum) const
The base class for a user to derive their own Kokkos kernels where the residual is of the form...
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