https://mooseframework.inl.gov
VectorCoupledTimeDerivative.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  "Time derivative Kernel that acts on a coupled vector variable. Weak form: "
20  "$(\\vec{\\psi}_i, \\frac{\\partial \\vec{v_h}}{\\partial t})$.");
21  params.addRequiredCoupledVar("v", "Coupled vector variable");
22  return params;
23 }
24 
26  : VectorKernel(parameters),
27  _v_dot(coupledVectorDot("v")),
28  _dv_dot(coupledVectorDotDu("v")),
29  _v_var(coupled("v"))
30 {
31 }
32 
33 Real
35 {
36  return _test[_i][_qp] * _v_dot[_qp];
37 }
38 
39 Real
41 {
42  return 0.0;
43 }
44 
45 Real
47 {
48  if (jvar == _v_var)
49  return _test[_i][_qp] * _phi[_j][_qp] * _dv_dot[_qp];
50 
51  return 0.0;
52 }
const VectorVariableTestValue & _test
the current test function
Definition: VectorKernel.h:65
VectorCoupledTimeDerivative(const InputParameters &parameters)
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
static InputParameters validParams()
virtual Real computeQpOffDiagJacobian(unsigned int jvar) override
This is the virtual that derived classes should override for computing an off-diagonal Jacobian compo...
const VectorVariableValue & _v_dot
virtual Real computeQpResidual() override
Compute this Kernel's contribution to the residual at the current quadrature point.
virtual Real computeQpJacobian() override
Compute this Kernel's contribution to the Jacobian at the current quadrature point.
unsigned int _i
current index for the test function
Definition: KernelBase.h:58
registerMooseObject("MooseApp", VectorCoupledTimeDerivative)
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
This method adds a coupled variable name pair.
unsigned int _j
current index for the shape function
Definition: KernelBase.h:61
const VectorVariablePhiValue & _phi
the current shape functions
Definition: VectorKernel.h:71
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static InputParameters validParams()
Definition: VectorKernel.C:21
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...
This calculates the time derivative for a coupled vector variable.
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:43