www.mooseframework.org
CoupledTimeDerivative.C
Go to the documentation of this file.
1 /****************************************************************/
2 /* DO NOT MODIFY THIS HEADER */
3 /* MOOSE - Multiphysics Object Oriented Simulation Environment */
4 /* */
5 /* (c) 2010 Battelle Energy Alliance, LLC */
6 /* ALL RIGHTS RESERVED */
7 /* */
8 /* Prepared by Battelle Energy Alliance, LLC */
9 /* Under Contract No. DE-AC07-05ID14517 */
10 /* With the U. S. Department of Energy */
11 /* */
12 /* See COPYRIGHT for full restrictions */
13 /****************************************************************/
14 
15 #include "CoupledTimeDerivative.h"
16 
17 template <>
20 {
22  params.addClassDescription("Time derivative Kernel that acts on a coupled variable. Weak form: "
23  "$(\\psi_i, \\frac{\\partial v_h}{\\partial t})$.");
24  params.addRequiredCoupledVar("v", "Coupled variable");
25  return params;
26 }
27 
29  : Kernel(parameters), _v_dot(coupledDot("v")), _dv_dot(coupledDotDu("v")), _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 VariableTestValue & _test
the current test function
Definition: KernelBase.h:152
const VariableValue & _dv_dot
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual Real computeQpJacobian() override
Compute this Kernel&#39;s contribution to the Jacobian at the current quadrature point.
CoupledTimeDerivative(const InputParameters &parameters)
InputParameters validParams< CoupledTimeDerivative >()
virtual Real computeQpResidual() override
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
unsigned int _i
current index for the test function
Definition: KernelBase.h:146
const unsigned int _v_var
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
This method adds a coupled variable name pair.
const VariablePhiValue & _phi
the current shape functions
Definition: KernelBase.h:158
unsigned int _j
current index for the shape function
Definition: KernelBase.h:149
InputParameters validParams< Kernel >()
Definition: Kernel.C:30
virtual Real computeQpOffDiagJacobian(unsigned int jvar) override
This is the virtual that derived classes should override for computing an off-diagonal Jacobian compo...
const VariableValue & _v_dot
Definition: Kernel.h:25
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...
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:131