www.mooseframework.org
ConservativeAdvection.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 #include "ConservativeAdvection.h"
15 
16 template <>
19 {
21  params.addClassDescription("Conservative form of $\\nabla \\cdot \\vec{v} u$ which in its weak "
22  "form is given by: $(-\\nabla \\psi_i, \\vec{v} u)$.");
23  params.addRequiredParam<RealVectorValue>("velocity", "Velocity vector");
24  return params;
25 }
26 
28  : Kernel(parameters), _velocity(getParam<RealVectorValue>("velocity"))
29 {
30 }
31 
32 Real
34 {
35  return -_u[_qp] * (_velocity * _grad_test[_i][_qp]);
36 }
37 
38 Real
40 {
41  return -_phi[_j][_qp] * (_velocity * _grad_test[_i][_qp]);
42 }
VectorValue< Real > RealVectorValue
Definition: Assembly.h:40
const VariableTestGradient & _grad_test
gradient of the test function
Definition: KernelBase.h:155
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
ConservativeAdvection(const InputParameters &parameters)
void addRequiredParam(const std::string &name, const std::string &doc_string)
This method adds a parameter and documentation string to the InputParameters object that will be extr...
virtual Real computeQpResidual()
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 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 computeQpJacobian()
Compute this Kernel&#39;s contribution to the Jacobian at the current quadrature point.
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...
InputParameters validParams< ConservativeAdvection >()
const VariableValue & _u
Holds the solution at current quadrature points.
Definition: Kernel.h:51
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:131