https://mooseframework.inl.gov
Public Member Functions | List of all members
Moose::Kokkos::VariableTestGradient Class Reference

#include <KokkosVariableValue.h>

Public Member Functions

KOKKOS_FUNCTION Real3 operator() (ResidualDatum &datum, unsigned int i, unsigned int qp) const
 Get the gradient of the current test function. More...
 

Detailed Description

Definition at line 92 of file KokkosVariableValue.h.

Member Function Documentation

◆ operator()()

KOKKOS_FUNCTION Real3 Moose::Kokkos::VariableTestGradient::operator() ( ResidualDatum datum,
unsigned int  i,
unsigned int  qp 
) const
inline

Get the gradient of the current test function.

Parameters
datumThe ResidualDatum object of the current thread
iThe element-local DOF index
qpThe local quadrature point index
Returns
The gradient of the test function

Definition at line 102 of file KokkosVariableValue.h.

103  {
104  auto & elem = datum.elem();
105  auto side = datum.side();
106  auto fe = datum.ife();
107 
108  return datum.J(qp) *
109  (side == libMesh::invalid_uint
110  ? datum.assembly().getGradPhi(elem.subdomain, elem.type, fe)(i, qp)
111  : datum.assembly().getGradPhiFace(elem.subdomain, elem.type, fe)(side)(i, qp));
112  }
const unsigned int invalid_uint
KOKKOS_FUNCTION const auto & getGradPhiFace(ContiguousSubdomainID subdomain, unsigned int elem_type, unsigned int fe_type) const
Get the gradient of face shape functions of a FE type for an element type and subdomain.
KOKKOS_FUNCTION const auto & getGradPhi(ContiguousSubdomainID subdomain, unsigned int elem_type, unsigned int fe_type) const
Get the gradient of shape functions of a FE type for an element type and subdomain.
KOKKOS_FUNCTION const Real33 & J(const unsigned int qp)
Get the inverse of Jacobian matrix | dxi/dx deta/dx dzeta/dx | | dxi/dy deta/dy dzeta/dy | | dxi/dz d...
Definition: KokkosDatum.h:115
KOKKOS_FUNCTION unsigned int side() const
Get the side index.
Definition: KokkosDatum.h:90
KOKKOS_FUNCTION unsigned int ife() const
Get the variable FE type ID.
Definition: KokkosDatum.h:312
KOKKOS_FUNCTION const Assembly & assembly() const
Get the Kokkos assembly.
Definition: KokkosDatum.h:68
KOKKOS_FUNCTION const ElementInfo & elem() const
Get the element information object.
Definition: KokkosDatum.h:80

The documentation for this class was generated from the following file: