Line data Source code
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 "AnisotropicDiffusion.h" 11 : 12 : registerMooseObject("MooseApp", AnisotropicDiffusion); 13 : 14 : InputParameters 15 14344 : AnisotropicDiffusion::validParams() 16 : { 17 14344 : InputParameters p = Kernel::validParams(); 18 14344 : 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 14344 : p.addRequiredParam<RealTensorValue>("tensor_coeff", 21 : "The Tensor to multiply the Diffusion operator by"); 22 14344 : return p; 23 0 : } 24 : 25 39 : AnisotropicDiffusion::AnisotropicDiffusion(const InputParameters & parameters) 26 39 : : Kernel(parameters), _k(getParam<RealTensorValue>("tensor_coeff")) 27 : { 28 39 : } 29 : 30 : Real 31 705024 : AnisotropicDiffusion::computeQpResidual() 32 : { 33 705024 : return (_k * _grad_u[_qp]) * _grad_test[_i][_qp]; 34 : } 35 : 36 : Real 37 155648 : AnisotropicDiffusion::computeQpJacobian() 38 : { 39 155648 : return (_k * _grad_phi[_j][_qp]) * _grad_test[_i][_qp]; 40 : }