www.mooseframework.org
ConvectiveFluxBC.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 "ConvectiveFluxBC.h"
11 
13 
16 {
18  params.set<Real>("rate") = 7500;
19  params.set<Real>("initial") = 500;
20  params.set<Real>("final") = 500;
21  params.set<Real>("duration") = 0.0;
22  params.addClassDescription(
23  "Determines boundary values via the initial and final values, flux, and exposure duration");
24  return params;
25 }
26 
28  : IntegratedBC(parameters),
29  _initial(getParam<Real>("initial")),
30  _final(getParam<Real>("final")),
31  _rate(getParam<Real>("rate")),
32  _duration(getParam<Real>("duration"))
33 {
34 }
35 
36 Real
38 {
39  Real value;
40 
41  if (_t < _duration)
43  else
44  value = _final;
45 
46  return -(_test[_i][_qp] * _rate * (value - _u[_qp]));
47 }
48 
49 Real
51 {
52  return -(_test[_i][_qp] * _rate * (-_phi[_j][_qp]));
53 }
const VariableTestValue & _test
test function values (in QPs)
Definition: IntegratedBC.h:97
Real _initial
Ratio of u to du/dn.
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sin(_arg) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tan
ConvectiveFluxBC(const InputParameters &parameters)
T & set(const std::string &name, bool quiet_mode=false)
Returns a writable reference to the named parameters.
static InputParameters validParams()
Definition: IntegratedBC.C:22
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
unsigned int _i
i-th, j-th index for enumerating test and shape functions
registerMooseObject("MooseApp", ConvectiveFluxBC)
const VariablePhiValue & _phi
shape function values (in QPs)
Definition: IntegratedBC.h:90
unsigned int _qp
quadrature point index
Base class for deriving any boundary condition of a integrated type.
Definition: IntegratedBC.h:18
virtual const OutputTools< Real >::VariableValue & value()
The value of the variable this object is operating on.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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 computeQpJacobian() override
Method for computing the diagonal Jacobian at quadrature points.
static InputParameters validParams()
Factory constructor, takes parameters so that all derived classes can be built using the same constru...
const VariableValue & _u
the values of the unknown variable this BC is acting on
Definition: IntegratedBC.h:104
virtual Real computeQpResidual() override
Method for computing the residual at quadrature points.
const Real pi