https://mooseframework.inl.gov
ArrayVacuumBC.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 "ArrayVacuumBC.h"
11 
13 
16 {
18  params.addParam<RealEigenVector>("alpha", "Ratio between directional gradient and solution");
19  params.addClassDescription("Imposes the Robin boundary condition $\\partial_n "
20  "\\vec{u}=-\\frac{\\vec{\\alpha}}{2}\\vec{u}$.");
21  return params;
22 }
23 
25  : ArrayIntegratedBC(parameters),
26  _alpha(isParamValid("alpha") ? getParam<RealEigenVector>("alpha")
27  : RealEigenVector::Ones(_count))
28 {
29  if (_alpha.size() != _count)
30  paramError(
31  "alpha", "Number of values must equal number of variable components (", _count, ").");
32 
33  _alpha /= 2;
34 }
35 
36 void
38 {
39  residual = _alpha.cwiseProduct(_u[_qp]) * _test[_i][_qp];
40 }
41 
44 {
45  return _test[_i][_qp] * _phi[_j][_qp] * _alpha;
46 }
const ArrayVariablePhiValue & _phi
shape function values (in QPs)
void paramError(const std::string &param, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
Definition: MooseBase.h:435
virtual RealEigenVector computeQpJacobian() override
Method for computing the diagonal Jacobian at quadrature points.
Definition: ArrayVacuumBC.C:43
ArrayVacuumBC(const InputParameters &parameters)
Definition: ArrayVacuumBC.C:24
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
const ArrayVariableValue & _u
the values of the unknown variable this BC is acting on
const unsigned int _count
Number of components of the array variable.
static InputParameters validParams()
unsigned int _qp
quadrature point index
virtual void computeQpResidual(RealEigenVector &residual) override
Method for computing the residual at quadrature points, to be filled in residual. ...
Definition: ArrayVacuumBC.C:37
registerMooseObject("MooseApp", ArrayVacuumBC)
const ArrayVariableTestValue & _test
test function values (in QPs)
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...
Base class for deriving any boundary condition of a integrated type.
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:146
static InputParameters validParams()
Definition: ArrayVacuumBC.C:15
RealEigenVector _alpha
Ratio of u to du/dn.
Definition: ArrayVacuumBC.h:26