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 
14 template <>
17 {
19  params.set<Real>("rate") = 7500;
20  params.set<Real>("initial") = 500;
21  params.set<Real>("final") = 500;
22  params.set<Real>("duration") = 0.0;
23  params.addClassDescription(
24  "Determines boundary values via the initial and final values, flux, and exposure duration");
25  return params;
26 }
27 
29  : IntegratedBC(parameters),
30  _initial(getParam<Real>("initial")),
31  _final(getParam<Real>("final")),
32  _rate(getParam<Real>("rate")),
33  _duration(getParam<Real>("duration"))
34 {
35 }
36 
37 Real
39 {
40  Real value;
41 
42  if (_t < _duration)
43  value = _initial + (_final - _initial) * std::sin((0.5 / _duration) * libMesh::pi * _t);
44  else
45  value = _final;
46 
47  return -(_test[_i][_qp] * _rate * (value - _u[_qp]));
48 }
49 
50 Real
52 {
53  return -(_test[_i][_qp] * _rate * (-_phi[_j][_qp]));
54 }
const VariableTestValue & _test
test function values (in QPs)
Definition: IntegratedBC.h:68
Real _initial
Ratio of u to du/dn.
InputParameters validParams< IntegratedBC >()
Definition: IntegratedBC.C:23
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...
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:61
unsigned int _qp
quadrature point index
Base class for deriving any boundary condition of a integrated type.
Definition: IntegratedBC.h:24
virtual const OutputTools< Real >::VariableValue & value()
The value of the variable this object is operating on.
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
const VariableValue & _u
the values of the unknown variable this BC is acting on
Definition: IntegratedBC.h:75
virtual Real computeQpResidual() override
method for computing the residual at quadrature points