https://mooseframework.inl.gov
MassLumpedTimeDerivative.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 
12 // MOOSE includes
13 #include "Assembly.h"
14 #include "MooseVariableFE.h"
15 
16 #include "libmesh/quadrature.h"
17 
19 
22 {
24  params.addClassDescription(
25  "Lumped formulation of the time derivative $\\frac{\\partial u}{\\partial t}$. Its "
26  "corresponding weak form is $\\dot{u_i}(\\psi_i, 1)$ where $\\dot{u_i}$ denotes the time "
27  "derivative of the solution coefficient associated with node $i$.");
28  return params;
29 }
30 
32  : TimeKernel(parameters), _u_dot_nodal(_var.dofValuesDot())
33 {
34 }
35 
36 Real
38 {
39  return _test[_i][_qp] * _u_dot_nodal[_i];
40 }
41 
42 Real
44 {
45  return _test[_i][_qp] * _du_dot_du[_qp];
46 }
47 
48 void
50 {
52  for (_i = 0; _i < _test.size(); _i++)
53  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
56 }
MooseVariable & _var
This is a regular kernel so we cast to a regular MooseVariable.
Definition: Kernel.h:72
const MooseArray< Real > & _JxW
The current quadrature point weight value.
Definition: KernelBase.h:52
unsigned int number() const
Get variable number coming from libMesh.
virtual void computeJacobian() override
Compute this Kernel&#39;s contribution to the diagonal Jacobian entries.
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.
registerMooseObject("MooseApp", MassLumpedTimeDerivative)
DenseMatrix< Number > _local_ke
Holds local Jacobian entries as they are accumulated by this Kernel.
const VariableTestValue & _test
the current test function
Definition: Kernel.h:75
static InputParameters validParams()
const QBase *const & _qrule
active quadrature rule
Definition: KernelBase.h:49
unsigned int _i
current index for the test function
Definition: KernelBase.h:58
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:55
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
All time kernels should inherit from this class.
Definition: TimeKernel.h:18
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...
MassLumpedTimeDerivative(const InputParameters &parameters)
const VariableValue & _u_dot_nodal
void prepareMatrixTag(Assembly &assembly, unsigned int ivar, unsigned int jvar)
Prepare data for computing element jacobian according to the active tags.
virtual Real computeQpResidual() override
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
static InputParameters validParams()
Definition: TimeKernel.C:20
const VariableValue & _du_dot_du
Derivative of u_dot with respect to u.
Definition: TimeKernel.h:40
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:43