62 virtual void compute()
override;
85 template <
typename Derived>
88 template <
typename Derived>
103 template <
typename Derived>
111 template <
typename Derived>
135 const unsigned int comp = 0)
const;
203 const unsigned int state)
override final;
206 template <
typename Derived>
219 auxkernel.computeElementInternal(auxkernel, datum);
222 template <
typename Derived>
234 auxkernel.computeNodeInternal(auxkernel, datum);
237 template <
typename Derived>
245 for (
unsigned int i = 0; i < datum.
n_dofs(); ++i)
250 for (
unsigned int j = 0; j < datum.
n_dofs(); ++j)
251 A[j + datum.
n_dofs() * i] = 0;
254 for (
unsigned int qp = 0; qp < datum.
n_qps(); ++qp)
256 const auto value = auxkernel.computeValue(qp, datum);
260 for (
unsigned int i = 0; i < datum.
n_dofs(); ++i)
262 const auto t = datum.
JxW(qp) *
_test(datum, i, qp);
266 for (
unsigned int j = 0; j < datum.
n_dofs(); ++j)
267 A[j + datum.
n_dofs() * i] += t *
_test(datum, j, qp);
279 template <
typename Derived>
283 auto value = auxkernel.computeValue(0, datum);
288 KOKKOS_FUNCTION
inline void 291 const unsigned int comp)
const 296 auto elem = datum.
elem().
id;
298 for (
unsigned int i = 0; i < datum.
n_dofs(); ++i)
299 sys.getVectorDofValue(sys.getElemLocalDofIndex(elem, i, var_num), tag) = values[i];
302 KOKKOS_FUNCTION
inline void 305 const unsigned int comp)
const 310 auto node = datum.
node();
312 sys.getVectorDofValue(sys.getNodeLocalDofIndex(node, 0, var_num), tag) =
value;
KOKKOS_FUNCTION TagID tag() const
Get the vector tag ID.
KOKKOS_FUNCTION unsigned int sys(unsigned int comp=0) const
Get the system number of a component.
KOKKOS_FUNCTION bool isNodal() const
Get whether this auxiliary kernel is nodal.
Base class for auxiliary kernels.
KOKKOS_FUNCTION ContiguousNodeID kokkosBoundaryNodeID(ThreadID tid) const
Get the contiguous node ID this Kokkos thread is operating on.
const unsigned int invalid_uint
KOKKOS_FUNCTION const Assembly & kokkosAssembly() const
Get the const reference of the Kokkos assembly.
KOKKOS_FUNCTION ContiguousElementID kokkosBlockNodeID(ThreadID tid) const
Get the contiguous node index this Kokkos thread is operating on.
The base class for a user to derive their own Kokkos auxiliary kernels.
KOKKOS_FUNCTION void setNodeSolution(const Real value, const AssemblyDatum &datum, const unsigned int comp=0) const
Set node value to the auxiliary solution vector.
const InputParameters & parameters() const
Get the parameters of the object.
virtual void compute() override
Dispatch calculation.
const bool _nodal
Flag whether this kernel is nodal.
static InputParameters validParams()
std::unique_ptr< DispatcherBase > _element_dispatcher
Kokkos functor dispatchers.
KOKKOS_INLINE_FUNCTION void choleskySolve(Real *const A, Real *const x, Real *const b, const unsigned int n)
Perform an in-place linear solve using Cholesky decomposition Matrix and right-hand-side vector are m...
std::unique_ptr< DispatcherBase > _node_dispatcher
Variable _kokkos_var
Kokkos variable.
The Kokkos interface that holds the host reference of the Kokkos systems and copies it to device duri...
The Kokkos interface that holds the host reference of the Kokkos mesh and copies it to device during ...
const bool _bnd
true if the kernel is boundary kernel, false if it is interior kernels
KOKKOS_FUNCTION void operator()(ElementLoop, const ThreadID tid, const Derived &auxkernel) const
The parallel computation entry functions called by Kokkos.
Scalar< Real > _t
TODO: Move to TransientInterface.
virtual void subdomainSetup() override final
Gets called when the subdomain changes (i.e.
KOKKOS_FUNCTION void computeNodeInternal(const Derived &auxkernel, AssemblyDatum &datum) const
Compute a node.
KOKKOS_FUNCTION ContiguousElementID kokkosBlockElementID(ThreadID tid) const
Get the contiguous element ID this Kokkos thread is operating on.
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
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 Kokkos interface that holds the host reference of the Kokkos assembly and copies it to device dur...
Scalar< int > _t_step
The number of the time step.
Scalar< Real > _dt_old
Size of the old time step.
VariableValue uOlder() const
Retrieve the older value of the variable that this kernel operates on.
constexpr unsigned int MAX_CACHED_DOF
Maximum number of DOFs to cache during residual and Jacobian computation.
const VariableTestValue _test
Current test function.
KOKKOS_FUNCTION ContiguousNodeID node() const
Get the contiguous node ID.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
KOKKOS_FUNCTION void setElementSolution(const Real *const values, const AssemblyDatum &datum, const unsigned int comp=0) const
Set element values to the auxiliary solution vector.
The Kokkos variable object that carries the coupled variable and tag information. ...
const VariableValue _u
Current solution.
KOKKOS_FUNCTION const ElementInfo & elem() const
Get the element information object.
The Kokkos object that holds thread-private data in the parallel operations of Kokkos kernels...
VariableValue uOld() const
Retrieve the old value of the variable that this kernel operates on.
KOKKOS_FUNCTION Real JxW(const unsigned int qp)
Get the transformed Jacobian weight.
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
The Kokkos wrapper classes for MOOSE-like variable value access.
KOKKOS_FUNCTION unsigned int var(unsigned int comp=0) const
Get the variable number of a component.
KOKKOS_FUNCTION void computeElementInternal(const Derived &auxkernel, AssemblyDatum &datum) const
The parallel computation bodies that can be customized in the derived class by defining them in the d...
Scalar< Real > _dt
Time step size.
KOKKOS_FUNCTION unsigned int n_dofs() const
Get the number of local DOFs.
AuxKernel(const InputParameters ¶meters)
Constructor.
void getKokkosMaterialPropertyHook(const std::string &prop_name_in, const unsigned int state) override final
Override of the MaterialPropertyInterface function to throw on material property request for nodal ke...
KOKKOS_FUNCTION void reinit()
Reset the reinit flag.
Scalar< const Real > _t_old
Old time.
ContiguousElementID id
Contiguous element ID.
KOKKOS_FUNCTION const Array< System > & kokkosSystems() const
Get the const reference of the Kokkos systems.