www.mooseframework.org
MaterialDerivativeTestKernel.C
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 
11 
13 
16 {
18  params.addClassDescription("Class used for testing derivatives of a scalar material property.");
19  return params;
20 }
21 
24 {
25 }
26 
27 Real
29 {
30  return _p[_qp] * _test[_i][_qp];
31 }
32 
33 Real
35 {
36  return _p_diag_derivative[_qp] * _phi[_j][_qp] * _test[_i][_qp];
37 }
38 
39 Real
41 {
42  // get the coupled variable number corresponding to jvar
43  const unsigned int cvar = mapJvarToCvar(jvar);
44  return (*_p_off_diag_derivatives[cvar])[_qp] * _phi[_j][_qp] * _test[_i][_qp];
45 }
virtual Real computeQpResidual() override
Compute this Kernel's contribution to the residual at the current quadrature point.
This kernel is used for testing derivatives of a material property.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual Real computeQpJacobian() override
Compute this Kernel's contribution to the Jacobian at the current quadrature point.
This kernel is used for testing derivatives of a material property.
const MaterialProperty< Real > & _p
material property for which to test derivatives
const MaterialProperty< Real > & _p_diag_derivative
material property for the diagonal derivative of the tested property
MaterialDerivativeTestKernel(const InputParameters &parameters)
const VariableTestValue & _test
the current test function
Definition: Kernel.h:75
virtual Real computeQpOffDiagJacobian(unsigned int jvar) override
For coupling standard variables.
unsigned int _i
current index for the test function
Definition: KernelBase.h:57
unsigned int mapJvarToCvar(unsigned int jvar)
Return index into the _coupled_moose_vars array for a given jvar.
unsigned int _j
current index for the shape function
Definition: KernelBase.h:60
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
const VariablePhiValue & _phi
the current shape functions
Definition: Kernel.h:81
std::vector< const MaterialProperty< Real > *> _p_off_diag_derivatives
material properties for the off-diagonal derivatives of the tested property
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:42
registerMooseObject("MooseApp", MaterialDerivativeTestKernel)