libMesh
assembly.C
Go to the documentation of this file.
1 // local includes
2 #include "assembly.h"
3 #include "rb_classes.h"
4 
5 // libMesh includes
6 #include "libmesh/sparse_matrix.h"
7 #include "libmesh/numeric_vector.h"
8 #include "libmesh/dense_matrix.h"
9 #include "libmesh/dense_submatrix.h"
10 #include "libmesh/dense_vector.h"
11 #include "libmesh/dense_subvector.h"
12 #include "libmesh/fe.h"
13 #include "libmesh/fe_interface.h"
14 #include "libmesh/fe_base.h"
15 #include "libmesh/elem_assembly.h"
16 #include "libmesh/quadrature_gauss.h"
17 #include "libmesh/boundary_info.h"
18 #include "libmesh/node.h"
19 
20 #ifdef LIBMESH_ENABLE_DIRICHLET
21 
22 // Bring in bits from the libMesh namespace.
23 // Just the bits we're using, since this is a header.
27 using libMesh::Number;
28 using libMesh::Point;
29 using libMesh::RBTheta;
30 using libMesh::Real;
32 
33 // Kronecker delta function
34 inline Real kronecker_delta(unsigned int i,
35  unsigned int j)
36 {
37  return i == j ? 1. : 0.;
38 }
39 
40 Real elasticity_tensor(unsigned int i,
41  unsigned int j,
42  unsigned int k,
43  unsigned int l)
44 {
45  // Define the Poisson ratio and Young's modulus
46  const Real nu = 0.3;
47  const Real E = 1.;
48 
49  // Define the Lame constants (lambda_1 and lambda_2) based on nu and E
50  const Real lambda_1 = E * nu / ((1. + nu) * (1. - 2.*nu));
51  const Real lambda_2 = 0.5 * E / (1. + nu);
52 
53  return lambda_1 * kronecker_delta(i, j) * kronecker_delta(k, l)
54  + lambda_2 * (kronecker_delta(i, k) * kronecker_delta(j, l) + kronecker_delta(i, l) * kronecker_delta(j, k));
55 }
56 
58 {
59  const unsigned int n_components = rb_sys.n_vars();
60 
61  // make sure we have three components
62  libmesh_assert_equal_to (n_components, 3);
63 
64  const unsigned int u_var = rb_sys.u_var;
65  const unsigned int v_var = rb_sys.v_var;
66  const unsigned int w_var = rb_sys.w_var;
67 
68  FEBase * elem_fe = nullptr;
69  c.get_element_fe(u_var, elem_fe);
70 
71  const std::vector<Real> & JxW = elem_fe->get_JxW();
72 
73  // The velocity shape function gradients at interior
74  // quadrature points.
75  const std::vector<std::vector<RealGradient>> & dphi = elem_fe->get_dphi();
76 
77  // Now we will build the affine operator
78  unsigned int n_qpoints = c.get_element_qrule().n_points();
79 
80  std::vector<unsigned int> n_var_dofs(n_components);
81  n_var_dofs[u_var] = c.get_dof_indices(u_var).size();
82  n_var_dofs[v_var] = c.get_dof_indices(v_var).size();
83  n_var_dofs[w_var] = c.get_dof_indices(w_var).size();
84 
85  for (unsigned int C_i = 0; C_i < n_components; C_i++)
86  {
87  unsigned int C_j = 0;
88  for (unsigned int C_k = 0; C_k < n_components; C_k++)
89  for (unsigned int C_l = 1; C_l < n_components; C_l++)
90  {
91  Real C_ijkl = elasticity_tensor(C_i, C_j, C_k, C_l);
92  for (unsigned int qp=0; qp<n_qpoints; qp++)
93  for (unsigned int i=0; i<n_var_dofs[C_i]; i++)
94  for (unsigned int j=0; j<n_var_dofs[C_k]; j++)
95  (c.get_elem_jacobian(C_i,C_k))(i,j) +=
96  JxW[qp]*(C_ijkl * dphi[i][qp](C_j)*dphi[j][qp](C_l));
97  }
98  }
99 
100  for (unsigned int C_i = 0; C_i < n_components; C_i++)
101  for (unsigned int C_j = 1; C_j < n_components; C_j++)
102  for (unsigned int C_k = 0; C_k < n_components; C_k++)
103  {
104  unsigned int C_l = 0;
105 
106  Real C_ijkl = elasticity_tensor(C_i, C_j, C_k, C_l);
107  for (unsigned int qp=0; qp<n_qpoints; qp++)
108  for (unsigned int i=0; i<n_var_dofs[C_i]; i++)
109  for (unsigned int j=0; j<n_var_dofs[C_k]; j++)
110  (c.get_elem_jacobian(C_i,C_k))(i,j) +=
111  JxW[qp]*(C_ijkl * dphi[i][qp](C_j)*dphi[j][qp](C_l));
112  }
113 
114 }
115 
117 {
118  const unsigned int n_components = rb_sys.n_vars();
119 
120  // make sure we have three components
121  libmesh_assert_equal_to (n_components, 3);
122 
123  const unsigned int u_var = rb_sys.u_var;
124  const unsigned int v_var = rb_sys.v_var;
125  const unsigned int w_var = rb_sys.w_var;
126 
127  FEBase * elem_fe = nullptr;
128  c.get_element_fe(u_var, elem_fe);
129 
130  const std::vector<Real> & JxW = elem_fe->get_JxW();
131 
132  // The velocity shape function gradients at interior
133  // quadrature points.
134  const std::vector<std::vector<RealGradient>> & dphi = elem_fe->get_dphi();
135 
136  // Now we will build the affine operator
137  unsigned int n_qpoints = c.get_element_qrule().n_points();
138 
139  std::vector<unsigned int> n_var_dofs(n_components);
140  n_var_dofs[u_var] = c.get_dof_indices(u_var).size();
141  n_var_dofs[v_var] = c.get_dof_indices(v_var).size();
142  n_var_dofs[w_var] = c.get_dof_indices(w_var).size();
143 
144  for (unsigned int C_i = 0; C_i < n_components; C_i++)
145  for (unsigned int C_j = 1; C_j < n_components; C_j++)
146  for (unsigned int C_k = 0; C_k < n_components; C_k++)
147  for (unsigned int C_l = 1; C_l < n_components; C_l++)
148  {
149  Real C_ijkl = elasticity_tensor(C_i,C_j,C_k,C_l);
150  for (unsigned int qp=0; qp<n_qpoints; qp++)
151  for (unsigned int i=0; i<n_var_dofs[C_i]; i++)
152  for (unsigned int j=0; j<n_var_dofs[C_k]; j++)
153  (c.get_elem_jacobian(C_i,C_k))(i,j) +=
154  JxW[qp]*(C_ijkl * dphi[i][qp](C_j)*dphi[j][qp](C_l));
155  }
156 }
157 
159 {
160  const unsigned int n_components = rb_sys.n_vars();
161 
162  // make sure we have three components
163  libmesh_assert_equal_to (n_components, 3);
164 
165  const unsigned int u_var = rb_sys.u_var;
166  const unsigned int v_var = rb_sys.v_var;
167  const unsigned int w_var = rb_sys.w_var;
168 
169  FEBase * elem_fe = nullptr;
170  c.get_element_fe(u_var, elem_fe);
171 
172  const std::vector<Real> & JxW = elem_fe->get_JxW();
173 
174  // The velocity shape function gradients at interior
175  // quadrature points.
176  const std::vector<std::vector<RealGradient>> & dphi = elem_fe->get_dphi();
177 
178  // Now we will build the affine operator
179  unsigned int n_qpoints = c.get_element_qrule().n_points();
180 
181  std::vector<unsigned int> n_var_dofs(n_components);
182  n_var_dofs[u_var] = c.get_dof_indices(u_var).size();
183  n_var_dofs[v_var] = c.get_dof_indices(v_var).size();
184  n_var_dofs[w_var] = c.get_dof_indices(w_var).size();
185 
186  for (unsigned int C_i = 0; C_i < n_components; C_i++)
187  {
188  unsigned int C_j = 0;
189 
190  for (unsigned int C_k = 0; C_k < n_components; C_k++)
191  {
192  unsigned int C_l = 0;
193 
194  Real C_ijkl = elasticity_tensor(C_i,C_j,C_k,C_l);
195  for (unsigned int qp=0; qp<n_qpoints; qp++)
196  for (unsigned int i=0; i<n_var_dofs[C_i]; i++)
197  for (unsigned int j=0; j<n_var_dofs[C_k]; j++)
198  (c.get_elem_jacobian(C_i,C_k))(i,j) +=
199  JxW[qp]*(C_ijkl * dphi[i][qp](C_j)*dphi[j][qp](C_l));
200  }
201  }
202 }
203 
205 {
207  (&c.get_elem(), c.side, BOUNDARY_ID_MAX_X))
208  {
209  const unsigned int u_var = 0;
210 
211  FEBase * side_fe = nullptr;
212  c.get_side_fe(u_var, side_fe);
213 
214  const std::vector<Real> & JxW_side = side_fe->get_JxW();
215 
216  const std::vector<std::vector<Real>> & phi_side = side_fe->get_phi();
217 
218  // The number of local degrees of freedom in each variable
219  const unsigned int n_u_dofs = c.get_dof_indices(u_var).size();
220 
221  // Now we will build the affine operator
222  unsigned int n_qpoints = c.get_side_qrule().n_points();
223  DenseSubVector<Number> & Fu = c.get_elem_residual(u_var);
224 
225  for (unsigned int qp=0; qp < n_qpoints; qp++)
226  for (unsigned int i=0; i < n_u_dofs; i++)
227  Fu(i) += JxW_side[qp] * (1. * phi_side[i][qp]);
228  }
229 }
230 
232 {
234  (&c.get_elem(), c.side, BOUNDARY_ID_MAX_X))
235  {
236  const unsigned int u_var = 0;
237  const unsigned int v_var = 1;
238 
239  FEBase * side_fe = nullptr;
240  c.get_side_fe(u_var, side_fe);
241 
242  const std::vector<Real> & JxW_side = side_fe->get_JxW();
243 
244  const std::vector<std::vector<Real>> & phi_side = side_fe->get_phi();
245 
246  // The number of local degrees of freedom in each variable
247  const unsigned int n_v_dofs = c.get_dof_indices(u_var).size();
248 
249  // Now we will build the affine operator
250  unsigned int n_qpoints = c.get_side_qrule().n_points();
251  DenseSubVector<Number> & Fv = c.get_elem_residual(v_var);
252 
253  for (unsigned int qp=0; qp < n_qpoints; qp++)
254  for (unsigned int i=0; i < n_v_dofs; i++)
255  Fv(i) += JxW_side[qp] * (1. * phi_side[i][qp]);
256  }
257 }
258 
260 {
262  (&c.get_elem(), c.side, BOUNDARY_ID_MAX_X))
263  {
264  const unsigned int u_var = 0;
265  const unsigned int w_var = 2;
266 
267  FEBase * side_fe = nullptr;
268  c.get_side_fe(u_var, side_fe);
269 
270  const std::vector<Real> & JxW_side = side_fe->get_JxW();
271 
272  const std::vector<std::vector<Real>> & phi_side = side_fe->get_phi();
273 
274  // The number of local degrees of freedom in each variable
275  const unsigned int n_w_dofs = c.get_dof_indices(w_var).size();
276 
277  // Now we will build the affine operator
278  unsigned int n_qpoints = c.get_side_qrule().n_points();
279  DenseSubVector<Number> & Fw = c.get_elem_residual(w_var);
280 
281  for (unsigned int qp=0; qp < n_qpoints; qp++)
282  for (unsigned int i=0; i < n_w_dofs; i++)
283  Fw(i) += JxW_side[qp] * (1. * phi_side[i][qp]);
284  }
285 }
286 
287 void AssemblyPointLoadX::get_nodal_rhs_values(std::map<numeric_index_type, Number> & values,
288  const System & sys,
289  const Node & node)
290 {
291  // First clear the values map
292  values.clear();
293 
294  if (sys.get_mesh().get_boundary_info().has_boundary_id
295  (&node, NODE_BOUNDARY_ID))
296  {
297  numeric_index_type dof_index =
298  node.dof_number(sys.number(), sys.variable_number("u"), 0);
299  values[dof_index] = 1.;
300  }
301 }
302 
303 void AssemblyPointLoadY::get_nodal_rhs_values(std::map<numeric_index_type, Number> & values,
304  const System & sys,
305  const Node & node)
306 {
307  // First clear the values map
308  values.clear();
309 
310  if (sys.get_mesh().get_boundary_info().has_boundary_id
311  (&node, NODE_BOUNDARY_ID))
312  {
313  numeric_index_type dof_index =
314  node.dof_number(sys.number(), sys.variable_number("v"), 0);
315  values[dof_index] = 1.;
316  }
317 }
318 
319 void AssemblyPointLoadZ::get_nodal_rhs_values(std::map<numeric_index_type, Number> & values,
320  const System & sys,
321  const Node & node)
322 {
323  // First clear the values map
324  values.clear();
325 
326  if (sys.get_mesh().get_boundary_info().has_boundary_id
327  (&node, NODE_BOUNDARY_ID))
328  {
329  numeric_index_type dof_index =
330  node.dof_number(sys.number(), sys.variable_number("w"), 0);
331  values[dof_index] = 1.;
332  }
333 }
334 
336 {
337  const unsigned int u_var = rb_sys.u_var;
338  const unsigned int v_var = rb_sys.v_var;
339  const unsigned int w_var = rb_sys.w_var;
340 
341  FEBase * elem_fe = nullptr;
342  c.get_element_fe(u_var, elem_fe);
343 
344  const std::vector<Real> & JxW = elem_fe->get_JxW();
345 
346  // The velocity shape function gradients at interior
347  // quadrature points.
348  const std::vector<std::vector<RealGradient>>& dphi = elem_fe->get_dphi();
349 
350  // The number of local degrees of freedom in each variable
351  const unsigned int n_u_dofs = c.get_dof_indices(u_var).size();
352  const unsigned int n_v_dofs = c.get_dof_indices(v_var).size();
353  const unsigned int n_w_dofs = c.get_dof_indices(w_var).size();
354 
355  // Now we will build the affine operator
356  unsigned int n_qpoints = c.get_element_qrule().n_points();
357 
358  DenseSubMatrix<Number> & Kuu = c.get_elem_jacobian(u_var, u_var);
359  DenseSubMatrix<Number> & Kvv = c.get_elem_jacobian(v_var, v_var);
360  DenseSubMatrix<Number> & Kww = c.get_elem_jacobian(w_var, w_var);
361 
362  for (unsigned int qp=0; qp<n_qpoints; qp++)
363  {
364  for (unsigned int i=0; i<n_u_dofs; i++)
365  for (unsigned int j=0; j<n_u_dofs; j++)
366  Kuu(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
367 
368  for (unsigned int i=0; i<n_v_dofs; i++)
369  for (unsigned int j=0; j<n_v_dofs; j++)
370  Kvv(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
371 
372  for (unsigned int i=0; i<n_w_dofs; i++)
373  for (unsigned int j=0; j<n_w_dofs; j++)
374  Kww(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
375  }
376 }
377 
378 #endif // LIBMESH_ENABLE_DIRICHLET
elasticity_tensor
Real elasticity_tensor(unsigned int i, unsigned int j, unsigned int k, unsigned int l)
Definition: assembly.C:40
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::DiffContext::get_elem_residual
const DenseVector< Number > & get_elem_residual() const
Const accessor for element residual.
Definition: diff_context.h:249
libMesh::FEMContext::get_element_qrule
const QBase & get_element_qrule() const
Accessor for element interior quadrature rule for the dimension of the current _elem.
Definition: fem_context.h:790
libMesh::System::n_vars
unsigned int n_vars() const
Definition: system.h:2155
libMesh::FEMContext::side
unsigned char side
Current side for side_* to examine.
Definition: fem_context.h:1010
ElasticityRBConstruction::v_var
unsigned int v_var
Definition: rb_classes.h:123
libMesh::MeshBase::get_boundary_info
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:132
libMesh::FEInterface
This class provides an encapsulated access to all static public member functions of finite element cl...
Definition: fe_interface.h:65
libMesh::DiffContext::get_dof_indices
const std::vector< dof_id_type > & get_dof_indices() const
Accessor for element dof indices.
Definition: diff_context.h:367
AssemblyF2::boundary_assembly
virtual void boundary_assembly(FEMContext &c)
Perform the element boundary assembly.
Definition: assembly.C:259
InnerProductAssembly::interior_assembly
virtual void interior_assembly(FEMContext &c)
Perform the element interior assembly.
Definition: assembly.C:335
AssemblyPointLoadY::get_nodal_rhs_values
virtual void get_nodal_rhs_values(std::map< numeric_index_type, Number > &values, const System &sys, const Node &node)
Definition: assembly.C:303
libMesh::QBase::n_points
unsigned int n_points() const
Definition: quadrature.h:126
libMesh::RealGradient
RealVectorValue RealGradient
Definition: hp_coarsentest.h:49
AssemblyPointLoadX::get_nodal_rhs_values
virtual void get_nodal_rhs_values(std::map< numeric_index_type, Number > &values, const System &sys, const Node &node)
Definition: assembly.C:287
kronecker_delta
Real kronecker_delta(unsigned int i, unsigned int j)
Definition: assembly.C:34
libMesh::FEMContext::get_elem
const Elem & get_elem() const
Accessor for current Elem object.
Definition: fem_context.h:896
AssemblyF0::boundary_assembly
virtual void boundary_assembly(FEMContext &c)
Perform the element boundary assembly.
Definition: assembly.C:204
AssemblyF1::boundary_assembly
virtual void boundary_assembly(FEMContext &c)
Perform the element boundary assembly.
Definition: assembly.C:231
libMesh::FEBase
FEGenericBase< Real > FEBase
Definition: exact_error_estimator.h:39
AssemblyA0::interior_assembly
virtual void interior_assembly(FEMContext &c)
Perform the element interior assembly.
Definition: assembly.C:57
libMesh::System::get_mesh
const MeshBase & get_mesh() const
Definition: system.h:2083
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
ElasticityRBConstruction::w_var
unsigned int w_var
Definition: rb_classes.h:124
libMesh::numeric_index_type
dof_id_type numeric_index_type
Definition: id_types.h:99
ElasticityAssembly::rb_sys
ElasticityRBConstruction & rb_sys
The ElasticityRBConstruction object that will use this assembly.
Definition: assembly.h:56
AssemblyA2::interior_assembly
virtual void interior_assembly(FEMContext &c)
Perform the element interior assembly.
Definition: assembly.C:158
libMesh::FEMContext::get_element_fe
void get_element_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for interior finite element object for variable var for the largest dimension in the mesh.
Definition: fem_context.h:275
libMesh::ElemAssembly
ElemAssembly provides a per-element (interior and boundary) assembly functionality.
Definition: elem_assembly.h:38
AssemblyA1::interior_assembly
virtual void interior_assembly(FEMContext &c)
Perform the element interior assembly.
Definition: assembly.C:116
libMesh::DiffContext::get_elem_jacobian
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
Definition: diff_context.h:283
libMesh::RBTheta
This class is part of the rbOOmit framework.
Definition: rb_theta.h:46
libMesh::FEMContext::get_side_qrule
const QBase & get_side_qrule() const
Accessor for element side quadrature rule for the dimension of the current _elem.
Definition: fem_context.h:797
libMesh::FEMContext::get_side_fe
void get_side_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for edge/face (2D/3D) finite element object for variable var for the largest dimension in th...
Definition: fem_context.h:312
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
AssemblyPointLoadZ::get_nodal_rhs_values
virtual void get_nodal_rhs_values(std::map< numeric_index_type, Number > &values, const System &sys, const Node &node)
Definition: assembly.C:319
libMesh::BoundaryInfo::has_boundary_id
bool has_boundary_id(const Node *const node, const boundary_id_type id) const
Definition: boundary_info.C:972
ElasticityRBConstruction::u_var
unsigned int u_var
Variable numbers.
Definition: rb_classes.h:122
libMesh::FEMContext
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62