www.mooseframework.org
IntegratedBCBase.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 "BoundaryCondition.h"
13 #include "RandomInterface.h"
16 
17 // Forward declarations
18 class IntegratedBCBase;
19 template <typename>
23 
24 template <>
26 
31  public RandomInterface,
34 {
35 public:
37 
38  virtual ~IntegratedBCBase();
39 
40  virtual void computeResidual() = 0;
41  virtual void computeJacobian() = 0;
45  virtual void computeJacobianBlock(MooseVariableFEBase & jvar) = 0;
50  virtual void computeJacobianBlockScalar(unsigned int jvar) = 0;
51 
56  virtual void computeNonlocalJacobian() {}
57 
62  virtual void computeNonlocalOffDiagJacobian(unsigned int /* jvar */) {}
63 
64 protected:
66  const Elem * const & _current_elem;
68  const Real & _current_elem_volume;
70  const unsigned int & _current_side;
72  const Elem * const & _current_side_elem;
74  const Real & _current_side_volume;
75 
77  unsigned int _qp;
79  const QBase * const & _qrule;
87  unsigned int _i, _j;
88 
91  std::vector<MooseVariableFEBase *> _save_in;
92  std::vector<AuxVariableName> _save_in_strings;
93 
96  std::vector<MooseVariableFEBase *> _diag_save_in;
97  std::vector<AuxVariableName> _diag_save_in_strings;
98 
99  virtual Real computeQpJacobian() { return 0; }
104  virtual Real computeQpOffDiagJacobian(unsigned int /*jvar*/) { return 0; }
105 };
106 
MooseVariableFE< VectorValue< Real > > VectorMooseVariable
Interface for objects that need parallel consistent random numbers without patterns over the course o...
InputParameters validParams< IntegratedBCBase >()
std::vector< MooseVariableFEBase * > _diag_save_in
Class for stuff related to variables.
Definition: Adaptivity.h:29
std::vector< MooseVariableFEBase * > _save_in
virtual Real computeQpJacobian()
const Elem *const & _current_elem
current element
virtual void computeNonlocalJacobian()
Compute this IntegratedBCBase&#39;s contribution to the diagonal Jacobian entries corresponding to nonloc...
const Elem *const & _current_side_elem
current side element
MooseVariableFE< Real > MooseVariable
IntegratedBCBase(const InputParameters &parameters)
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
virtual void computeJacobianBlockScalar(unsigned int jvar)=0
Computes jacobian block with respect to a scalar variable.
virtual void computeJacobian()=0
unsigned int _qp
quadrature point index
virtual void computeNonlocalOffDiagJacobian(unsigned int)
Computes d-residual / d-jvar...
const MooseArray< Point > & _q_point
active quadrature points
const Real & _current_elem_volume
Volume of the current element.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseObject.h:65
std::vector< AuxVariableName > _save_in_strings
const MooseArray< Real > & _coord
coordinate transformation
const Real & _current_side_volume
Volume of the current side.
bool _has_save_in
The aux variables to save the residual contributions to.
Base class for creating new types of boundary conditions.
const QBase *const & _qrule
active quadrature rule
const unsigned int & _current_side
current side of the current element
std::vector< AuxVariableName > _diag_save_in_strings
virtual void computeResidual()=0
An interface for accessing Materials.
virtual void computeJacobianBlock(MooseVariableFEBase &jvar)=0
Computes d-ivar-residual / d-jvar...
virtual Real computeQpOffDiagJacobian(unsigned int)
This is the virtual that derived classes should override for computing an off-diagonal jacobian compo...
Intermediate base class that ties together all the interfaces for getting MooseVariableFEBases with t...
virtual ~IntegratedBCBase()
Base class for deriving any boundary condition of a integrated type.
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
const MooseArray< Real > & _JxW
transformed Jacobian weights