Go to the documentation of this file.
13 #include "libmesh/quadrature.h"
32 const unsigned int n_test =
_test.size();
36 for (
_i = 0;
_i < n_test;
_i++)
44 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
46 var->sys().solution().add_vector(
_local_re, var->dofIndices());
57 const unsigned int n_test =
_test.size();
62 for (
_i = 0;
_i < n_test;
_i++)
70 const unsigned int rows = ke.m();
71 DenseVector<Number> diag(rows);
72 for (
unsigned int i = 0; i < rows; i++)
75 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
77 var->sys().solution().add_vector(diag, var->dofIndices());
84 size_t jvar_num = jvar.
number();
95 for (
_j = 0;
_j < phi_size;
_j++)
111 return RealGradient(0.0);
unsigned int _i
current index for the test function
virtual const OutputTools< Real >::VariableValue & value()
The value of the variable this object is operating on.
MooseVariableFEBase & getVariable(THREAD_ID tid, const std::string &var_name)
Gets a reference to a variable of with specified name.
std::vector< MooseVariableFEBase * > _diag_save_in
const VariablePhiValue & _phi
the current shape functions
Assembly & _assembly
Reference to this Kernel's assembly object.
virtual RealGradient precomputeQpJacobian()
Called before forming the jacobian for an element.
DenseVector< Number > _local_re
Holds residual entries as they are accumulated by this Kernel.
unsigned int _j
current index for the shape function
DenseVector< Number > & residualBlock(unsigned int var_num, TagID tag_id=0)
Get local residual block for a variable and a tag.
static InputParameters validParams()
virtual Real computeQpOffDiagJacobian(unsigned int)
This is the virtual that derived classes should override for computing an off-diagonal Jacobian compo...
unsigned int _qp
The current quadrature point index.
KernelGrad(const InputParameters ¶meters)
Factory constructor initializes all internal references needed for residual computation.
virtual void computeOffDiagJacobian(MooseVariableFEBase &jvar) override
Computes d-residual / d-jvar... storing the result in Ke.
The KernelGrad class is responsible for calculating the residuals in form:
std::vector< MooseVariableFEBase * > _save_in
virtual void computeResidual() override
Compute this Kernel's contribution to the residual.
THREAD_ID _tid
The thread ID for this kernel.
const MooseArray< Real > & _JxW
The current quadrature point weight value.
const VariableTestValue & _test
the current test function
defineLegacyParams(KernelGrad)
const VariableTestGradient & _grad_test
gradient of the test function
SystemBase & _sys
Reference to the EquationSystem object.
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
virtual Real computeQpResidual() override
Compute this Kernel's contribution to the residual at the current quadrature point.
virtual void computeJacobian() override
Compute this Kernel's contribution to the diagonal Jacobian entries.
DenseMatrix< Number > _local_ke
Holds residual entries as they are accumulated by this Kernel.
const QBase *const & _qrule
active quadrature rule
static InputParameters validParams()
MooseVariable & _var
This is a regular kernel so we cast to a regular MooseVariable.
DenseMatrix< Number > & jacobianBlock(unsigned int ivar, unsigned int jvar, TagID tag=0)
Get local Jacobian block for a pair of variables and a tag.
virtual RealGradient precomputeQpResidual()=0
Called before forming the residual for an element.
bool _has_save_in
The aux variables to save the residual contributions to.
const MooseArray< Real > & _coord
The scaling factor to convert from cartesian to another coordinate system (e.g rz,...
unsigned int number() const
Get variable number coming from libMesh.