https://mooseframework.inl.gov
LowerDIntegratedBC.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 "IntegratedBC.h"
13 
15 {
16 public:
18 
20 
21  virtual const MooseVariable & variable() const override { return _var; }
22  const MooseVariable & lowerDVariable() const { return _lowerd_var; }
23 
24  virtual void computeResidual() override;
25  virtual void computeJacobian() override;
26  virtual void computeOffDiagJacobian(unsigned int jvar) override;
27 
28 protected:
32  virtual Real computeLowerDQpResidual() = 0;
33 
38 
44 
49  const unsigned int jvar_num);
50 
55  const MooseVariableFEBase &)
56  {
57  return 0;
58  }
59 
63  virtual void initLowerDQpResidual() {}
64 
69 
75  const MooseVariableFEBase &)
76  {
77  }
78 
87 };
const MooseVariable & lowerDVariable() const
virtual void initLowerDQpResidual()
Put necessary evaluations depending on qp but independent on test functions here. ...
virtual Real computeLowerDQpResidual()=0
Method for computing the Lower part of residual at quadrature points.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void initLowerDQpJacobian(Moose::ConstraintJacobianType)
Put necessary evaluations depending on qp but independent on test and shape functions here...
This class provides an interface for common operations on field variables of both FE and FV types wit...
const VariableTestValue & _test_lambda
test functions
virtual const MooseVariable & variable() const override
Returns the variable that this object operates on.
LowerDIntegratedBC(const InputParameters &parameters)
virtual void computeLowerDJacobian(Moose::ConstraintJacobianType type)
Method for computing the LowerLower, PrimaryLower and LowerPrimary parts of Jacobian.
virtual Real computeLowerDQpJacobian(Moose::ConstraintJacobianType)=0
Method for computing the LowerLower, PrimaryLower and LowerPrimary parts of Jacobian at quadrature po...
virtual void computeJacobian() override
Compute this object's contribution to the diagonal Jacobian entries.
OutputTools< Real >::VariableTestValue VariableTestValue
Definition: MooseTypes.h:324
virtual Real computeLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType, const MooseVariableFEBase &)
Method for computing an off-diagonal jacobian component at quadrature points.
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:51
const VariableValue & _lambda
Holds the current solution at the current quadrature point on the face.
Base class for deriving any boundary condition of a integrated type.
Definition: IntegratedBC.h:18
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-ivar-residual / d-jvar...
static InputParameters validParams()
OutputTools< Real >::VariableValue VariableValue
Definition: MooseTypes.h:314
MooseVariable & _var
Definition: IntegratedBC.h:82
void computeLowerDOffDiagJacobian(Moose::ConstraintJacobianType type, const unsigned int jvar_num)
Method for computing an off-diagonal jacobian component.
forward declarations
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const VariablePhiValue & _phi_lambda
Shape functions.
ConstraintJacobianType
Definition: MooseTypes.h:796
const MooseVariable & _lowerd_var
Variable this kernel operates on.
const InputParameters & parameters() const
Get the parameters of the object.
virtual void initLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType, const MooseVariableFEBase &)
Put necessary evaluations depending on qp but independent on test and shape functions here for off-di...
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.