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"
15 
22 {
23 public:
25 
27 
28  void prepareShapes(unsigned int var_num) override final;
29 
30 protected:
32  const Elem * const & _current_elem;
36  const unsigned int & _current_side;
38  const Elem * const & _current_side_elem;
43 
45  unsigned int _qp;
47  const QBase * const & _qrule;
55  unsigned int _i, _j;
56 
59  std::vector<MooseVariableFEBase *> _save_in;
60  std::vector<AuxVariableName> _save_in_strings;
61 
64  std::vector<MooseVariableFEBase *> _diag_save_in;
65  std::vector<AuxVariableName> _diag_save_in_strings;
66 };
std::vector< MooseVariableFEBase * > _diag_save_in
std::vector< MooseVariableFEBase * > _save_in
const Elem *const & _current_elem
current element
const Elem *const & _current_side_elem
current side element
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
const BoundaryID & _current_boundary_id
The currenty boundary id.
unsigned int _qp
quadrature point index
const MooseArray< Point > & _q_point
active quadrature points
const Real & _current_elem_volume
Volume of the current element.
std::vector< AuxVariableName > _save_in_strings
boundary_id_type BoundaryID
void prepareShapes(unsigned int var_num) override final
Prepare shape functions.
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
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::vector< AuxVariableName > _diag_save_in_strings
An interface for accessing Materials.
Intermediate base class that ties together all the interfaces for getting MooseVariableFEBases with t...
static InputParameters validParams()
Base class for deriving any boundary condition of a integrated type.
const InputParameters & parameters() const
Get the parameters of the object.
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
const MooseArray< Real > & _JxW
transformed Jacobian weights