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 
16 {
18  p.addClassDescription("Anisotropic diffusion kernel $\\nabla \\cdot -\\widetilde{k} \\nabla u$ "
19  "with weak form given by $(\\nabla \\psi_i, \\widetilde{k} \\nabla u)$.");
20  p.addRequiredParam<RealTensorValue>("tensor_coeff",
21  "The Tensor to multiply the Diffusion operator by");
22  return p;
23 }
24 
26  : Kernel(parameters), _k(getParam<RealTensorValue>("tensor_coeff"))
27 {
28 }
29 
30 Real
32 {
33  return (_k * _grad_u[_qp]) * _grad_test[_i][_qp];
34 }
35 
36 Real
38 {
39  return (_k * _grad_phi[_j][_qp]) * _grad_test[_i][_qp];
40 }
const VariableGradient & _grad_u
Holds the solution gradient at the current quadrature points.
Definition: Kernel.h:90
static InputParameters validParams()
Definition: Kernel.C:23
static InputParameters validParams()
This kernel implements the Laplacian operator multiplied by a 2nd order tensor giving anisotropic (di...
const VariablePhiGradient & _grad_phi
gradient of the shape function
Definition: Kernel.h:84
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.
TensorValue< Real > RealTensorValue
unsigned int _i
current index for the test function
Definition: KernelBase.h:57
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:60
registerMooseObject("MooseApp", AnisotropicDiffusion)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const VariableTestGradient & _grad_test
gradient of the test function
Definition: Kernel.h:78
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:42
AnisotropicDiffusion(const InputParameters &parameters)