www.mooseframework.org
VectorIntegratedBC.h
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 #pragma once
11 
12 #include "IntegratedBCBase.h"
13 #include "MooseVariableInterface.h"
14 
18 class VectorIntegratedBC : public IntegratedBCBase, public MooseVariableInterface<RealVectorValue>
19 {
20 public:
22 
24 
25  virtual const VectorMooseVariable & variable() const override { return _var; }
26 
27  virtual void computeResidual() override;
28  virtual void computeJacobian() override;
32  virtual void computeOffDiagJacobian(unsigned int jvar) override;
37  void computeOffDiagJacobianScalar(unsigned int jvar) override;
38 
39 protected:
43  virtual Real computeQpResidual() = 0;
44 
48  virtual Real computeQpJacobian() { return 0; }
49 
53  virtual Real computeQpOffDiagJacobian(unsigned int /*jvar*/) { return 0; }
54 
58  virtual Real computeQpOffDiagJacobianScalar(unsigned int /*jvar*/) { return 0; }
59 
61 
64 
65  // shape functions
66 
69 
70  // test functions
71 
74 
75  // solution variable
76 
79 };
Class for stuff related to variables.
Definition: Adaptivity.h:31
virtual Real computeQpOffDiagJacobianScalar(unsigned int)
Method for computing an off-diagonal jacobian component from a scalar var.
const VectorMooseVariable & _var
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
OutputTools< RealVectorValue >::VariableValue VectorVariableValue
Definition: MooseTypes.h:319
OutputTools< RealVectorValue >::VariableTestValue VectorVariableTestValue
Definition: MooseTypes.h:329
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual Real computeQpOffDiagJacobian(unsigned int)
Method for computing an off-diagonal jacobian component at quadrature points.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-ivar-residual / d-jvar...
virtual Real computeQpJacobian()
Method for computing the diagonal Jacobian at quadrature points.
Base class for deriving any boundary condition of a integrated type.
const MooseArray< Point > & _normals
normals at quadrature points
const VectorVariableTestValue & _test
test function values (in QPs)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real computeQpResidual()=0
Method for computing the residual at quadrature points.
OutputTools< RealVectorValue >::VariablePhiValue VectorVariablePhiValue
Definition: MooseTypes.h:324
Interface for objects that need to get values of MooseVariables.
Base class for deriving any boundary condition of a integrated type.
const InputParameters & parameters() const
Get the parameters of the object.
static InputParameters validParams()
const VectorVariableValue & _u
the values of the unknown variable this BC is acting on
VectorIntegratedBC(const InputParameters &parameters)
virtual const VectorMooseVariable & variable() const override
Returns the variable that this object operates on.
const VectorVariablePhiValue & _phi
shape function values (in QPs)