https://mooseframework.inl.gov
FunctionDiffusion.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 "FunctionDiffusion.h"
11 #include "Function.h"
12 
14 
17 {
19  params.addClassDescription("Diffusion with a function coefficient.");
20  params.addParam<FunctionName>("function", 1.0, "Function multiplier for diffusion term.");
21  params.addCoupledVar("v",
22  "Coupled concentration variable for kernel to operate on; if this "
23  "is not specified, the kernel's nonlinear variable will be used as "
24  "usual");
25  return params;
26 }
27 
29  : Diffusion(parameters),
30  _function(getFunction("function")),
31  _grad_v(isCoupled("v") ? coupledGradient("v") : _grad_u),
32  _v_var(isCoupled("v") ? getVar("v", 0) : nullptr),
33  _grad_v_phi(isCoupled("v") ? _v_var->gradPhi() : _grad_phi)
34 {
35 }
36 
37 Real
39 {
41 }
42 
43 Real
45 {
47 }
virtual Real computeQpJacobian() override
Compute this Kernel&#39;s contribution to the Jacobian at the current quadrature point.
The Laplacian operator with a function coefficient.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
FunctionDiffusion(const InputParameters &parameters)
This kernel implements the Laplacian operator: $ u $.
Definition: Diffusion.h:18
static InputParameters validParams()
Definition: Diffusion.C:15
const VariablePhiGradient & _grad_v_phi
Gradient of the shape function.
virtual Real computeQpResidual() override
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
unsigned int _i
current index for the test function
Definition: KernelBase.h:58
static InputParameters validParams()
void addCoupledVar(const std::string &name, const std::string &doc_string)
This method adds a coupled variable name pair.
const Function & _function
Function coefficient.
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
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...
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
virtual Real value(Real t, const Point &p) const
Override this to evaluate the scalar function at point (t,x,y,z), by default this returns zero...
Definition: Function.C:44
registerMooseObject("MooseApp", FunctionDiffusion)
const VariableGradient & _grad_v
Gradient of the concentration variable for kernel to operate on.
const MooseArray< Point > & _q_point
The physical location of the element&#39;s quadrature Points, indexed by _qp.
Definition: KernelBase.h:46
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:43