https://mooseframework.inl.gov
Diffusion.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 "Diffusion.h"
11 
12 registerMooseObject("MooseApp", Diffusion);
13 
16 {
18  params.addClassDescription("The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak "
19  "form of $(\\nabla \\phi_i, \\nabla u_h)$.");
20  return params;
21 }
22 
23 Diffusion::Diffusion(const InputParameters & parameters) : Kernel(parameters) {}
24 
25 Real
27 {
28  return _grad_u[_qp] * _grad_test[_i][_qp];
29 }
30 
31 Real
33 {
34  return _grad_phi[_j][_qp] * _grad_test[_i][_qp];
35 }
const VariableGradient & _grad_u
Holds the solution gradient at the current quadrature points.
Definition: Kernel.h:90
static InputParameters validParams()
Definition: Kernel.C:24
const VariablePhiGradient & _grad_phi
gradient of the shape function
Definition: Kernel.h:84
virtual Real computeQpJacobian() override
Compute this Kernel's contribution to the Jacobian at the current quadrature point.
Definition: Diffusion.C:32
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
registerMooseObject("MooseApp", Diffusion)
This kernel implements the Laplacian operator: $ u $.
Definition: Diffusion.h:18
static InputParameters validParams()
Definition: Diffusion.C:15
virtual Real computeQpResidual() override
Compute this Kernel's contribution to the residual at the current quadrature point.
Definition: Diffusion.C:26
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
const VariableTestGradient & _grad_test
gradient of the test function
Definition: Kernel.h:78
Diffusion(const InputParameters &parameters)
Definition: Diffusion.C:23
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...
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:43