-
Notifications
You must be signed in to change notification settings - Fork 1k
/
BoundaryCondition.h
191 lines (154 loc) · 4.28 KB
/
BoundaryCondition.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
#ifndef BOUNDARYCONDITION_H
#define BOUNDARYCONDITION_H
// local includes
#include "Moose.h"
#include "ValidParams.h"
#include "MooseArray.h"
#include "PDEBase.h"
#include "MaterialPropertyInterface.h"
//forward declarations
class BoundaryCondition;
class MooseSystem;
class DofData;
class FaceData;
class AuxData;
namespace libMesh
{
class Elem;
template<class T> class DenseMatrix;
}
template<>
InputParameters validParams<BoundaryCondition>();
/**
* Extrapolation of a Kernel for BC usage. Children of this class should override
* computeQpResidual() for use by computeSideResidual.
*/
class BoundaryCondition :
public PDEBase,
protected MaterialPropertyInterface
{
public:
/**
* Factory constructor, takes parameters so that all derived classes can be built using the same
* constructor.
*/
BoundaryCondition(std::string name, MooseSystem & moose_system, InputParameters parameters);
virtual ~BoundaryCondition(){}
/**
* Boundary ID the BC is active on.
*
* @return The boundary ID.
*/
unsigned int boundaryID();
/**
* Computes the residual for the current side.
*/
virtual void computeResidual();
/**
* Computes the jacobian for the current side.
*/
virtual void computeJacobian();
/**
* Computes d-ivar-residual / d-jvar... storing the result in Ke.
*/
void computeJacobianBlock(DenseMatrix<Number> & Ke, unsigned int ivar, unsigned int jvar);
/**
* For use with non-integrated / nodal (Dirichlet) BC's
*/
void computeAndStoreResidual();
/**
* Computes the integral over the current side.
*/
Real computeIntegral();
protected:
/**
* Convenience reference to the DofData object inside of MooseSystem
*/
DofData & _dof_data;
/**
* Convenience reference to the FaceData object inside of MooseSystem
*/
FaceData & _face_data;
/**
* This is the virtual that derived classes should override for computing the residual.
*/
virtual Real computeQpResidual()=0;
/**
* This is the virtual that derived classes should override for computing the Jacobian.
*/
virtual Real computeQpJacobian();
/**
* This is the virtual that derived classes should override for computing an off-diagonal jacobian component.
*/
virtual Real computeQpOffDiagJacobian(unsigned int jvar);
/**
* This is the virtual that derived classes should override for computing the volume integral of kernel.
*/
virtual Real computeQpIntegral();
/**
* Boundary ID this BC is active on.
*/
unsigned int _boundary_id;
/**
* The current side as an element
* Only valid if there is a Dirichlet BC
* on this side.
*/
Elem * _side_elem;
/**
* Normal vectors at the quadrature points.
*/
const std::vector<Point>& _normals;
/**
* Current side.
*/
unsigned int & _current_side;
/**
* Current side element.
*/
const Elem * & _current_side_elem;
/**
* Current node for nodal BC's
*/
const Node * & _current_node;
/**
* Current residual vector. Only valid for nodal BC's.
*/
NumericVector<Number> * & _current_residual;
/**
* Holds the current solution at the current quadrature point on the face.
*/
MooseArray<Real> & _u;
/**
* Holds the current solution gradient at the current quadrature point on the face.
*/
MooseArray<RealGradient> & _grad_u;
/**
* Holds the current solution second derivative at the current quadrature point on the face
*/
MooseArray<RealTensor> & _second_u;
/**
* Returns a reference (that can be stored) to a coupled variable's value.
*
* @param name The name the kernel wants to refer to the variable as.
*/
MooseArray<Real> & coupledValue(std::string var_name, int i = 0);
/**
* Returns a reference (that can be stored) to a coupled variable's gradient.
*
* @param name The name the kernel wants to refer to the variable as.
*/
MooseArray<RealGradient> & coupledGradient(std::string var_name, int i = 0);
/** Side shape function.
*/
const std::vector<std::vector<Real> > & _test;
/**
* Gradient of side shape function.
*/
const std::vector<std::vector<RealGradient> > & _grad_test;
/**
* Second derivative of side shape function.
*/
const std::vector<std::vector<RealTensor> > & _second_test;
};
#endif //BOUNDARYCONDITION_H