https://mooseframework.inl.gov
IntegratedBC.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 
19 {
20 public:
22 
24 
25  virtual const MooseVariable & 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  virtual void computeResidualAndJacobian() override;
39 
40 protected:
44  virtual Real computeQpResidual() = 0;
45 
49  virtual Real computeQpJacobian() { return 0; }
50 
54  virtual Real computeQpOffDiagJacobian(unsigned int /*jvar*/) { return 0; }
55 
59  virtual Real computeQpOffDiagJacobianScalar(unsigned int jvar)
60  {
61  // Backwards compatibility
62  return computeQpOffDiagJacobian(jvar);
63  }
64 
68  virtual void precalculateQpResidual() {}
69 
74  virtual void precalculateQpJacobian() {}
75 
80  virtual void precalculateQpOffDiagJacobian(const MooseVariableFEBase & /*jvar*/) {}
81 
83 
86 
87  // shape functions
88 
93 
94  // test functions
95 
100 
101  // unknown
102 
104  const VariableValue & _u;
107 };
const VariableTestValue & _test
test function values (in QPs)
Definition: IntegratedBC.h:97
OutputTools< Real >::VariableGradient VariableGradient
Definition: MooseTypes.h:315
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-ivar-residual / d-jvar...
Definition: IntegratedBC.C:139
virtual Real computeQpJacobian()
Method for computing the diagonal Jacobian at quadrature points.
Definition: IntegratedBC.h:49
const MooseArray< Point > & _normals
normals at quadrature points
Definition: IntegratedBC.h:85
Class for stuff related to variables.
Definition: Adaptivity.h:31
virtual Real computeQpOffDiagJacobian(unsigned int)
Method for computing an off-diagonal jacobian component at quadrature points.
Definition: IntegratedBC.h:54
IntegratedBC(const InputParameters &parameters)
Definition: IntegratedBC.C:29
virtual const MooseVariable & variable() const override
Returns the variable that this object operates on.
Definition: IntegratedBC.h:25
const VariableGradient & _grad_u
the gradient of the unknown variable this BC is acting on
Definition: IntegratedBC.h:106
static InputParameters validParams()
Definition: IntegratedBC.C:23
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
This class provides an interface for common operations on field variables of both FE and FV types wit...
const VariablePhiValue & _phi
shape function values (in QPs)
Definition: IntegratedBC.h:90
const VariablePhiGradient & _grad_phi
gradients of shape functions (in QPs)
Definition: IntegratedBC.h:92
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
Definition: IntegratedBC.C:109
virtual Real computeQpResidual()=0
Method for computing the residual at quadrature points.
OutputTools< Real >::VariableTestValue VariableTestValue
Definition: MooseTypes.h:324
virtual void precalculateQpResidual()
Insertion point for evaluations that depend on qp but are independent of the test functions...
Definition: IntegratedBC.h:68
void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
Definition: IntegratedBC.C:171
virtual Real computeQpOffDiagJacobianScalar(unsigned int jvar)
Method for computing an off-diagonal jacobian component from a scalar var.
Definition: IntegratedBC.h:59
Base class for deriving any boundary condition of a integrated type.
Definition: IntegratedBC.h:18
virtual void computeResidualAndJacobian() override
Compute this object&#39;s contribution to the residual and Jacobian simultaneously.
Definition: IntegratedBC.C:185
MooseVariable & _var
Definition: IntegratedBC.h:82
OutputTools< Real >::VariablePhiGradient VariablePhiGradient
Definition: MooseTypes.h:320
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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.
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
Definition: IntegratedBC.C:85
const VariableTestGradient & _grad_test
gradients of test functions (in QPs)
Definition: IntegratedBC.h:99
OutputTools< Real >::VariableTestGradient VariableTestGradient
Definition: MooseTypes.h:325
const VariableValue & _u
the values of the unknown variable this BC is acting on
Definition: IntegratedBC.h:104
virtual void precalculateQpJacobian()
Insertion point for evaluations that depend on qp but are independent of the test and shape functions...
Definition: IntegratedBC.h:74
virtual void precalculateQpOffDiagJacobian(const MooseVariableFEBase &)
Insertion point for evaluations that depend on qp but are independent of the test and shape functions...
Definition: IntegratedBC.h:80