www.mooseframework.org
FunctionNeumannBC.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 "FunctionNeumannBC.h"
16 #include "Function.h"
17 
18 template <>
21 {
23  params.addRequiredParam<FunctionName>("function", "The function.");
24  params.addClassDescription("Imposes the integrated boundary condition "
25  "$\\frac{\\partial u}{\\partial n}=h(t,\\vec{x})$, "
26  "where $h$ is a (possibly) time and space-dependent MOOSE Function.");
27  return params;
28 }
29 
31  : IntegratedBC(parameters), _func(getFunction("function"))
32 {
33 }
34 
35 Real
37 {
38  return -_test[_i][_qp] * _func.value(_t, _q_point[_qp]);
39 }
const VariableTestValue & _test
test function values (in QPs)
Definition: IntegratedBC.h:105
virtual Real value(Real t, const Point &p)
Override this to evaluate the scalar function at point (t,x,y,z), by default this returns zero...
Definition: Function.C:43
InputParameters validParams< IntegratedBC >()
Definition: IntegratedBC.C:28
InputParameters validParams< FunctionNeumannBC >()
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...
Base class for deriving any boundary condition of a integrated type.
Definition: IntegratedBC.h:33
Function & _func
The function being used for setting the value.
FunctionNeumannBC(const InputParameters &parameters)
unsigned int _qp
quadrature point index
Definition: IntegratedBC.h:83
const MooseArray< Point > & _q_point
active quadrature points
Definition: IntegratedBC.h:87
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...
virtual Real computeQpResidual() override
unsigned int _i
i-th, j-th index for enumerating test and shape functions
Definition: IntegratedBC.h:93