www.mooseframework.org
AnisotropicDiffusion.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 
15 #include "AnisotropicDiffusion.h"
16 
17 template <>
20 {
22  p.addClassDescription("Anisotropic diffusion kernel $\\nabla \\cdot -\\widetilde{k} \\nabla u$ "
23  "with weak form given by $(\\nabla \\psi_i, \\widetilde{k} \\nabla u)$.");
24  p.addRequiredParam<RealTensorValue>("tensor_coeff",
25  "The Tensor to multiply the Diffusion operator by");
26  return p;
27 }
28 
30  : Kernel(parameters), _k(getParam<RealTensorValue>("tensor_coeff"))
31 {
32 }
33 
34 Real
36 {
37  return (_k * _grad_u[_qp]) * _grad_test[_i][_qp];
38 }
39 
40 Real
42 {
43  return (_k * _grad_phi[_j][_qp]) * _grad_test[_i][_qp];
44 }
const VariableGradient & _grad_u
Holds the solution gradient at the current quadrature points.
Definition: Kernel.h:54
const VariablePhiGradient & _grad_phi
gradient of the shape function
Definition: KernelBase.h:161
const VariableTestGradient & _grad_test
gradient of the test function
Definition: KernelBase.h:155
InputParameters validParams< AnisotropicDiffusion >()
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
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 computeQpJacobian() override
Compute this Kernel&#39;s contribution to the Jacobian at the current quadrature point.
unsigned int _i
current index for the test function
Definition: KernelBase.h:146
virtual Real computeQpResidual() override
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
unsigned int _j
current index for the shape function
Definition: KernelBase.h:149
InputParameters validParams< Kernel >()
Definition: Kernel.C:30
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...
TensorValue< Real > RealTensorValue
Definition: Assembly.h:45
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:131
AnisotropicDiffusion(const InputParameters &parameters)