www.mooseframework.org
ConvectiveFluxBC.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 "ConvectiveFluxBC.h"
16 
17 template <>
20 {
22  params.set<Real>("rate") = 7500;
23  params.set<Real>("initial") = 500;
24  params.set<Real>("final") = 500;
25  params.set<Real>("duration") = 0.0;
26  params.addClassDescription(
27  "Determines boundary values via the initial and final values, flux, and exposure duration");
28  return params;
29 }
30 
32  : IntegratedBC(parameters),
33  _initial(getParam<Real>("initial")),
34  _final(getParam<Real>("final")),
35  _rate(getParam<Real>("rate")),
36  _duration(getParam<Real>("duration"))
37 {
38 }
39 
40 Real
42 {
43  Real value;
44 
45  if (_t < _duration)
46  value = _initial + (_final - _initial) * std::sin((0.5 / _duration) * libMesh::pi * _t);
47  else
48  value = _final;
49 
50  return -(_test[_i][_qp] * _rate * (value - _u[_qp]));
51 }
52 
53 Real
55 {
56  return -(_test[_i][_qp] * _rate * (-_phi[_j][_qp]));
57 }
const VariableTestValue & _test
test function values (in QPs)
Definition: IntegratedBC.h:105
Real _initial
Ratio of u to du/dn.
InputParameters validParams< IntegratedBC >()
Definition: IntegratedBC.C:28
ConvectiveFluxBC(const InputParameters &parameters)
Factory constructor, takes parameters so that all derived classes can be built using the same constru...
InputParameters validParams< ConvectiveFluxBC >()
T & set(const std::string &name, bool quiet_mode=false)
Returns a writable reference to the named parameters.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const VariablePhiValue & _phi
shape function values (in QPs)
Definition: IntegratedBC.h:98
Base class for deriving any boundary condition of a integrated type.
Definition: IntegratedBC.h:33
virtual const VariableValue & value()
The value of the variable this object is operating on.
unsigned int _qp
quadrature point index
Definition: IntegratedBC.h:83
unsigned int _j
Definition: IntegratedBC.h:93
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 _i
i-th, j-th index for enumerating test and shape functions
Definition: IntegratedBC.h:93
virtual Real computeQpJacobian() override
const VariableValue & _u
the values of the unknown variable this BC is acting on
Definition: IntegratedBC.h:112
virtual Real computeQpResidual() override