www.mooseframework.org
TimeDerivative.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 "TimeDerivative.h"
11 
12 // MOOSE includes
13 #include "Assembly.h"
14 #include "MooseVariableFE.h"
15 
16 #include "libmesh/quadrature.h"
17 
19 
20 template <>
23 {
25  params.addClassDescription("The time derivative operator with the weak form of $(\\psi_i, "
26  "\\frac{\\partial u_h}{\\partial t})$.");
27  params.addParam<bool>("lumping", false, "True for mass matrix lumping, false otherwise");
28 
29  return params;
30 }
31 
33  : TimeKernel(parameters), _lumping(getParam<bool>("lumping"))
34 {
35 }
36 
37 Real
39 {
40  return _test[_i][_qp] * _u_dot[_qp];
41 }
42 
43 Real
45 {
46  return _test[_i][_qp] * _phi[_j][_qp] * _du_dot_du[_qp];
47 }
48 
49 void
51 {
52  if (_lumping)
53  {
55 
57  for (_i = 0; _i < _test.size(); _i++)
58  for (_j = 0; _j < _phi.size(); _j++)
59  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
61 
63  }
64  else
66 }
const VariableValue & _u_dot
Time derivative of u.
Definition: TimeKernel.h:33
virtual Real computeQpResidual() override
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
TimeDerivative(const InputParameters &parameters)
MooseVariable & _var
This is a regular kernel so we cast to a regular MooseVariable.
Definition: Kernel.h:54
const MooseArray< Real > & _JxW
The current quadrature point weight value.
Definition: KernelBase.h:159
unsigned int number() const
Get variable number coming from libMesh.
registerMooseObject("MooseApp", TimeDerivative)
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
InputParameters validParams< TimeKernel >()
Definition: TimeKernel.C:21
DenseMatrix< Number > _local_ke
Holds residual entries as they are accumulated by this Kernel.
virtual void precalculateJacobian()
Definition: KernelBase.h:122
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
Definition: KernelBase.h:139
const VariableTestValue & _test
the current test function
Definition: Kernel.h:57
virtual Real computeQpJacobian() override
Compute this Kernel&#39;s contribution to the Jacobian at the current quadrature point.
const QBase *const & _qrule
active quadrature rule
Definition: KernelBase.h:156
unsigned int _i
current index for the test function
Definition: KernelBase.h:165
void accumulateTaggedLocalMatrix()
Local Jacobian blocks will be appended by adding the current local kernel Jacobian.
const MooseArray< Real > & _coord
The scaling factor to convert from cartesian to another coordinate system (e.g rz, spherical, etc.)
Definition: KernelBase.h:162
virtual void computeJacobian() override
Compute this Kernel&#39;s contribution to the diagonal Jacobian entries.
Definition: Kernel.C:112
unsigned int _j
current index for the shape function
Definition: KernelBase.h:168
All time kernels should inherit from this class.
Definition: TimeKernel.h:24
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...
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an option parameter and a documentation string to the InputParameters object...
InputParameters validParams< TimeDerivative >()
virtual void computeJacobian() override
Compute this Kernel&#39;s contribution to the diagonal Jacobian entries.
void prepareMatrixTag(Assembly &assembly, unsigned int ivar, unsigned int jvar)
Prepare data for computing element jacobian according to the ative tags.
const VariablePhiValue & _phi
the current shape functions
Definition: Kernel.h:63
const VariableValue & _du_dot_du
Derivative of u_dot with respect to u.
Definition: TimeKernel.h:36
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:150