https://mooseframework.inl.gov
VectorDiffusion.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 
10 #include "VectorDiffusion.h"
11 
13 
16 {
18  params.addClassDescription(
19  "The Laplacian operator ($-\\nabla \\cdot \\nabla \\vec{u}$), with the weak "
20  "form of $(\\nabla \\vec{\\phi_i}, \\nabla \\vec{u_h})$.");
21  return params;
22 }
23 
25 
26 Real
28 {
29  return _grad_u[_qp].contract(_grad_test[_i][_qp]);
30 }
31 
32 Real
34 {
35  return _grad_phi[_j][_qp].contract(_grad_test[_i][_qp]);
36 }
This kernel implements the Laplacian operator: $ {u} {}$.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const VectorVariablePhiGradient & _grad_phi
gradient of the shape function
Definition: VectorKernel.h:74
const VectorVariableTestGradient & _grad_test
gradient of the test function
Definition: VectorKernel.h:68
const VectorVariableGradient & _grad_u
Holds the solution gradient at current quadrature points.
Definition: VectorKernel.h:80
VectorDiffusion(const InputParameters &parameters)
registerMooseObject("MooseApp", VectorDiffusion)
unsigned int _i
current index for the test function
Definition: KernelBase.h:58
unsigned int _j
current index for the shape function
Definition: KernelBase.h:61
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real computeQpJacobian() override
Compute this Kernel's contribution to the Jacobian at the current quadrature point.
static InputParameters validParams()
Definition: VectorKernel.C:21
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...
virtual Real computeQpResidual() override
Compute this Kernel's contribution to the residual at the current quadrature point.
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:43
static InputParameters validParams()