https://mooseframework.inl.gov
ArrayIntegratedBC.h
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 #pragma once
11 
12 #include "IntegratedBCBase.h"
13 #include "MooseVariableInterface.h"
14 #include "MooseVariableScalar.h"
15 
19 class ArrayIntegratedBC : public IntegratedBCBase, public MooseVariableInterface<RealEigenVector>
20 {
21 public:
23 
25 
26  virtual const ArrayMooseVariable & variable() const override { return _var; }
27 
28  virtual void computeResidual() override;
29  virtual void computeJacobian() override;
33  virtual void computeOffDiagJacobian(unsigned int jvar) override;
38  void computeOffDiagJacobianScalar(unsigned int jvar) override;
39 
40 protected:
44  virtual void computeQpResidual(RealEigenVector & residual) = 0;
45 
50 
55 
61 
65  virtual void initQpResidual() {}
66 
70  virtual void initQpJacobian() {}
71 
76  virtual void initQpOffDiagJacobian(const MooseVariableFEBase &) {}
77 
79 
82 
85 
88 
91 
93  const unsigned int _count;
94 
95 private:
98 };
const MooseArray< Point > & _normals
normals at quadrature points
RealEigenVector _work_vector
Work vector for residual.
Class for stuff related to variables.
Definition: Adaptivity.h:31
const ArrayVariablePhiValue & _phi
shape function values (in QPs)
virtual void computeQpResidual(RealEigenVector &residual)=0
Method for computing the residual at quadrature points, to be filled in residual. ...
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const ArrayVariableValue & _u
the values of the unknown variable this BC is acting on
This class provides an interface for common operations on field variables of both FE and FV types wit...
const unsigned int _count
Number of components of the array variable.
static InputParameters validParams()
virtual RealEigenMatrix computeQpOffDiagJacobianScalar(const MooseVariableScalar &jvar)
This is the virtual that derived classes should override for computing a full Jacobian component...
ArrayIntegratedBC(const InputParameters &parameters)
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
virtual RealEigenMatrix computeQpOffDiagJacobian(const MooseVariableFEBase &jvar)
Method for computing an off-diagonal jacobian component at quadrature points.
void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
virtual void initQpJacobian()
Put necessary evaluations depending on qp but independent on test and shape functions here...
ArrayMooseVariable & _var
virtual void initQpResidual()
Put necessary evaluations depending on qp but independent on test functions here. ...
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
Definition: MooseTypes.h:149
virtual const ArrayMooseVariable & variable() const override
Returns the variable that this object operates on.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-ivar-residual / d-jvar...
virtual void initQpOffDiagJacobian(const MooseVariableFEBase &)
Put necessary evaluations depending on qp but independent on test and shape functions here for off-di...
const ArrayVariableTestValue & _test
test function values (in QPs)
Interface for objects that need to get values of MooseVariables.
Class for scalar variables (they are different).
Base class for deriving any boundary condition of a integrated type.
const InputParameters & parameters() const
Get the parameters of the object.
Base class for deriving any boundary condition of a integrated type.
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:146
OutputTools< RealEigenVector >::VariableTestValue ArrayVariableTestValue
Definition: MooseTypes.h:359
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
virtual RealEigenVector computeQpJacobian()
Method for computing the diagonal Jacobian at quadrature points.