72 template <
typename Derived>
78 ::Kokkos::abort(
"Default computeQpJacobian() should never be called. Make sure you properly " 79 "redefined this method in your class without typos.");
93 template <
typename Derived>
101 "Default computeQpOffDiagJacobian() should never be called. Make sure you properly " 102 "redefined this method in your class without typos.");
113 template <
typename Derived>
117 return &Kernel::computeQpJacobian<Derived>;
119 template <
typename Derived>
122 return &Kernel::computeQpOffDiagJacobian<Derived>;
129 template <
typename Derived>
131 KOKKOS_FUNCTION
void operator()(ResidualLoop,
const ThreadID tid,
const Derived & kernel)
const;
132 template <
typename Derived>
133 KOKKOS_FUNCTION
void operator()(JacobianLoop,
const ThreadID tid,
const Derived & kernel)
const;
134 template <
typename Derived>
150 template <
typename Derived>
157 template <
typename Derived>
164 template <
typename Derived>
196 template <
typename Derived>
211 kernel.computeResidualInternal(kernel, datum);
214 template <
typename Derived>
229 kernel.computeJacobianInternal(kernel, datum);
232 template <
typename Derived>
241 if (!sys.isVariableActive(jvar,
kokkosMesh().getElementInfo(elem).subdomain))
249 kernel.computeOffDiagJacobianInternal(kernel, datum);
252 template <
typename Derived>
258 [&](
Real * local_re,
const unsigned int ib,
const unsigned int ie)
260 for (
unsigned int qp = 0; qp < datum.
n_qps(); ++qp)
261 for (
unsigned int i = ib; i < ie; ++i)
262 local_re[i] += datum.
JxW(qp) * kernel.template computeQpResidual<Derived>(i, qp, datum);
266 template <
typename Derived>
272 [&](
Real * local_ke,
const unsigned int ib,
const unsigned int ie,
const unsigned int j)
274 for (
unsigned int qp = 0; qp < datum.
n_qps(); ++qp)
275 for (
unsigned int i = ib; i < ie; ++i)
277 datum.
JxW(qp) * kernel.template computeQpJacobian<Derived>(i, j, qp, datum);
281 template <
typename Derived>
287 [&](
Real * local_ke,
const unsigned int ib,
const unsigned int ie,
const unsigned int j)
289 for (
unsigned int qp = 0; qp < datum.
n_qps(); ++qp)
290 for (
unsigned int i = ib; i < ie; ++i)
291 local_ke[i] += datum.
JxW(qp) * kernel.template computeQpOffDiagJacobian<Derived>(
292 i, j, datum.
jvar(), qp, datum);
KOKKOS_FUNCTION unsigned int sys(unsigned int comp=0) const
Get the system number of a component.
KOKKOS_FUNCTION void computeJacobianInternal(AssemblyDatum &datum, function body) const
The common loop structure template for computing elemental Jacobian.
virtual void computeResidual() override
Dispatch residual calculation.
KOKKOS_FUNCTION void set_local_parallel(const unsigned int local_thread_id, const unsigned int num_local_threads)
Set local parallelization option.
const unsigned int invalid_uint
KOKKOS_FUNCTION const Assembly & kokkosAssembly() const
Get the const reference of the Kokkos assembly.
KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int, const unsigned int, const unsigned int, const unsigned int, AssemblyDatum &) const
Compute off-diagonal Jacobian contribution on a quadrature point.
const VariableValue _u
Current solution at quadrature points.
static auto defaultJacobian()
Functions used to check if users have overriden the hook methods, whose calculations can be skipped w...
const InputParameters & parameters() const
Get the parameters of the object.
const VariableTestValue _test
Current test function.
Thread _thread
Kokkos thread object.
KOKKOS_FUNCTION unsigned int jvar() const
Get the coupled variable number.
The Kokkos wrapper classes for MOOSE-like shape function access.
KOKKOS_FUNCTION void computeResidualInternal(const Derived &kernel, AssemblyDatum &datum) const
The parallel computation bodies that can be customized in the derived class by defining them in the d...
KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int, 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 computeJacobianInternal(const Derived &kernel, AssemblyDatum &datum) const
Compute diagonal Jacobian.
static InputParameters validParams()
KOKKOS_FUNCTION void operator()(ResidualLoop, const ThreadID tid, const Derived &kernel) const
The parallel computation entry functions called by Kokkos.
KOKKOS_FUNCTION unsigned int n_qps() const
Get the number of local quadrature points.
KOKKOS_FUNCTION const System & kokkosSystem(unsigned int sys) const
Get the const reference of a Kokkos system.
The base class for Kokkos kernels.
const VariableGradient _grad_u
Gradient of the current solution at quadrature points.
The Kokkos wrapper classes for MOOSE-like variable value access.
virtual void computeJacobian() override
Dispatch diagonal and off-diagonal Jacobian calculation.
KOKKOS_FUNCTION ContiguousElementID kokkosBlockElementID(Moose::Kokkos::ThreadID tid) const
Get the contiguous element ID this Kokkos thread is operating on.
KOKKOS_FUNCTION const Mesh & kokkosMesh() const
Get the const reference of the Kokkos mesh.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const VariableTestGradient _grad_test
Gradient of the current test function.
Kernel(const InputParameters ¶meters)
Constructor.
KOKKOS_FUNCTION void computeOffDiagJacobianInternal(const Derived &kernel, AssemblyDatum &datum) const
Compute off-diagonal Jacobian.
const VariablePhiGradient _grad_phi
Gradient of the current shape function.
The Kokkos object that holds thread-private data in the parallel operations of Kokkos kernels...
KOKKOS_FUNCTION Real JxW(const unsigned int qp)
Get the transformed Jacobian weight.
Variable _kokkos_var
Kokkos variable.
static auto defaultOffDiagJacobian()
KOKKOS_FUNCTION unsigned int var(unsigned int comp=0) const
Get the variable number of a component.
const VariablePhiValue _phi
Current shape function.
The base class for a user to derive their own Kokkos kernels.
KOKKOS_FUNCTION thread_id_type size() const
Get the total thread pool size.
KOKKOS_FUNCTION const auto & getCoupling(unsigned int var) const
Get the list of off-diagonal coupled variable numbers of a variable.
KOKKOS_FUNCTION void computeResidualInternal(AssemblyDatum &datum, function body) const
The common loop structure template for computing elemental residual.
KOKKOS_FUNCTION const Array< System > & kokkosSystems() const
Get the const reference of the Kokkos systems.