https://mooseframework.inl.gov
MaterialDerivativeStdVectorRealTestKernel.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
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(
19  "Class used for testing derivatives of a std::vector<Real> material property.");
20  params.addRequiredParam<MaterialPropertyName>(
21  "material_property", "Name of material property for which derivatives are to be tested.");
22  params.addRequiredCoupledVar("args", "List of variables the material property depends on");
23  params.addRequiredParam<unsigned int>("i", "std::vector component");
24  return params;
25 }
26 
28  const InputParameters & parameters)
30  _n_vars(_coupled_moose_vars.size()),
31  _name(getParam<MaterialPropertyName>("material_property")),
32  _p(getMaterialProperty<Real>(_name)),
33  _p_off_diag_derivatives(_n_vars),
34  _p_diag_derivative(
35  getMaterialPropertyDerivative<std::vector<Real>>(_name, VariableName(_var.name()))),
36  _component_i(getParam<unsigned int>("i"))
37 {
38  for (unsigned int m = 0; m < _n_vars; ++m)
40  _name, VariableName(_coupled_moose_vars[m]->name()));
41 }
42 
43 Real
45 {
46  // This would be normally _p[_qp] * _test[_i][_qp];
47  // But this kernel sits on a variable where another kernel with Real-valued
48  // derivatives is. Our residual then must be 0, becuase we want to contribute
49  // with just jacobians without modifying the residual
50  return 0;
51 }
52 
53 Real
55 {
56  if (_p_diag_derivative[_qp].size() > 0)
58  else
59  return 0;
60 }
61 
62 Real
64 {
65  // get the coupled variable number corresponding to jvar
66  const unsigned int cvar = mapJvarToCvar(jvar);
67  if ((*_p_off_diag_derivatives[cvar])[_qp].size() > 0)
68  return (*_p_off_diag_derivatives[cvar])[_qp][_component_i] * _phi[_j][_qp] * _test[_i][_qp];
69  else
70  return 0;
71 }
static InputParameters validParams()
const MaterialProperty< std::vector< Real > > & _p_diag_derivative
material property for the diagonal derivative of the tested property
void addRequiredParam(const std::string &name, const std::string &doc_string)
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative(const std::string &base, const std::vector< VariableName > &c)
MaterialDerivativeStdVectorRealTestKernel(const InputParameters &parameters)
virtual Real computeQpOffDiagJacobian(unsigned int jvar) override
const unsigned int _component_i
the component of the std::vector material property
const std::string name
Definition: Setup.h:20
std::vector< const MaterialProperty< std::vector< Real > > * > _p_off_diag_derivatives
material properties for the off-diagonal derivatives of the tested property
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
std::vector< MooseVariableFieldBase *> _coupled_moose_vars
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addClassDescription(const std::string &doc_string)
registerMooseObject("ThermalHydraulicsTestApp", MaterialDerivativeStdVectorRealTestKernel)
Kernel for testing derivatives of a std::vector<Real> material property.
void ErrorVector unsigned int
const unsigned int _n_vars
number of nonlinear variables