www.mooseframework.org
AnisotropicDiffusion.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 "AnisotropicDiffusion.h"
11 
13 
14 template <>
17 {
19  p.addClassDescription("Anisotropic diffusion kernel $\\nabla \\cdot -\\widetilde{k} \\nabla u$ "
20  "with weak form given by $(\\nabla \\psi_i, \\widetilde{k} \\nabla u)$.");
21  p.addRequiredParam<RealTensorValue>("tensor_coeff",
22  "The Tensor to multiply the Diffusion operator by");
23  return p;
24 }
25 
27  : Kernel(parameters), _k(getParam<RealTensorValue>("tensor_coeff"))
28 {
29 }
30 
31 Real
33 {
34  return (_k * _grad_u[_qp]) * _grad_test[_i][_qp];
35 }
36 
37 Real
39 {
40  return (_k * _grad_phi[_j][_qp]) * _grad_test[_i][_qp];
41 }
const VariableGradient & _grad_u
Holds the solution gradient at the current quadrature points.
Definition: Kernel.h:72
This kernel implements the Laplacian operator multiplied by a 2nd order tensor giving anisotropic (di...
InputParameters validParams< AnisotropicDiffusion >()
const VariablePhiGradient & _grad_phi
gradient of the shape function
Definition: Kernel.h:66
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:165
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:168
InputParameters validParams< Kernel >()
Definition: Kernel.C:24
registerMooseObject("MooseApp", AnisotropicDiffusion)
const VariableTestGradient & _grad_test
gradient of the test function
Definition: Kernel.h:60
Definition: Kernel.h:20
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: MooseTypes.h:132
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:150
AnisotropicDiffusion(const InputParameters &parameters)