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 template <typename Derived>
42 class KernelGrad : public Kernel<Derived>
43 {
44  usingKokkosKernelMembers(Derived);
45 
46 public:
48 
53 
58 
67  KOKKOS_FUNCTION Real3 precomputeQpJacobian(const unsigned int /* j */,
68  const unsigned int /* qp */,
69  ResidualDatum & /* datum */) const
70  {
71  return Real3(0);
72  }
74 
79  KOKKOS_FUNCTION void computeResidualInternal(const Derived * kernel, ResidualDatum & datum) const;
81  KOKKOS_FUNCTION void computeJacobianInternal(const Derived * kernel, ResidualDatum & datum) const;
83 
84 protected:
89  virtual bool defaultJacobian() const override
90  {
91  return &Derived::precomputeQpJacobian == &KernelGrad::precomputeQpJacobian;
92  }
93 };
94 
95 template <typename Derived>
98 {
100  return params;
101 }
102 
103 template <typename Derived>
104 KernelGrad<Derived>::KernelGrad(const InputParameters & parameters) : Kernel<Derived>(parameters)
105 {
106 }
107 
108 template <typename Derived>
109 KOKKOS_FUNCTION void
110 KernelGrad<Derived>::computeResidualInternal(const Derived * kernel, ResidualDatum & datum) const
111 {
113  datum,
114  [&](Real * local_re, const unsigned int ib, const unsigned int ie)
115  {
116  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
117  {
118  datum.reinit();
119 
120  Real3 value = datum.JxW(qp) * kernel->precomputeQpResidual(qp, datum);
121 
122  for (unsigned int i = ib; i < ie; ++i)
123  local_re[i] += value * _grad_test(datum, i, qp);
124  }
125  });
126 }
127 
128 template <typename Derived>
129 KOKKOS_FUNCTION void
130 KernelGrad<Derived>::computeJacobianInternal(const Derived * kernel, ResidualDatum & datum) const
131 {
132  Real3 value;
133 
135  datum,
136  [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
137  {
138  for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
139  {
140  datum.reinit();
141 
142  unsigned int j_old = libMesh::invalid_uint;
143 
144  for (unsigned int ij = ijb; ij < ije; ++ij)
145  {
146  unsigned int i = ij % datum.n_jdofs();
147  unsigned int j = ij / datum.n_jdofs();
148 
149  if (j != j_old)
150  {
151  value = datum.JxW(qp) * kernel->precomputeQpJacobian(j, qp, datum);
152  j_old = j;
153  }
154 
155  local_ke[ij] += value * _grad_test(datum, i, qp);
156  }
157  }
158  });
159 }
160 
161 } // namespace Kokkos
162 } // namespace Moose
163 
164 #define usingKokkosKernelGradMembers(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()
KOKKOS_FUNCTION void computeJacobianInternal(const Derived *kernel, ResidualDatum &datum) const
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
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual bool defaultJacobian() const override
Get whether precomputeQpJacobian() was not defined in the derived class.
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 gr...
KernelGrad(const InputParameters &parameters)
Constructor.
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: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 base class for a user to derive their own Kokkos kernels where the residual is of the form...
KOKKOS_FUNCTION Real3 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...
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