www.mooseframework.org
CoupledTimeDerivative.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 
10 #include "CoupledTimeDerivative.h"
11 
13 
16 {
18  params.addClassDescription("Time derivative Kernel that acts on a coupled variable. Weak form: "
19  "$(\\psi_i, \\frac{\\partial v_h}{\\partial t})$.");
20  params.addRequiredCoupledVar("v", "Coupled variable");
21  return params;
22 }
23 
25  : Kernel(parameters), _v_dot(coupledDot("v")), _dv_dot(coupledDotDu("v")), _v_var(coupled("v"))
26 {
27 }
28 
29 Real
31 {
32  return _test[_i][_qp] * _v_dot[_qp];
33 }
34 
35 Real
37 {
38  return 0.0;
39 }
40 
41 Real
43 {
44  if (jvar == _v_var)
45  return _test[_i][_qp] * _phi[_j][_qp] * _dv_dot[_qp];
46 
47  return 0.0;
48 }
static InputParameters validParams()
Definition: Kernel.C:23
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's contribution to the Jacobian at the current quadrature point.
This calculates the time derivative for a coupled variable.
CoupledTimeDerivative(const InputParameters &parameters)
const VariableTestValue & _test
the current test function
Definition: Kernel.h:75
static InputParameters validParams()
virtual Real computeQpResidual() override
Compute this Kernel's contribution to the residual at the current quadrature point.
unsigned int _i
current index for the test function
Definition: KernelBase.h:57
const unsigned int _v_var
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:60
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real computeQpOffDiagJacobian(unsigned int jvar) override
For coupling standard variables.
const VariableValue & _v_dot
Definition: Kernel.h:15
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...
registerMooseObject("MooseApp", CoupledTimeDerivative)
const VariablePhiValue & _phi
the current shape functions
Definition: Kernel.h:81
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:42