Reference documentation for deal.II version 8.4.2
elasticity.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2015 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__integrators_elasticity_h
17 #define dealii__integrators_elasticity_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/exceptions.h>
22 #include <deal.II/base/quadrature.h>
23 #include <deal.II/lac/full_matrix.h>
24 #include <deal.II/fe/mapping.h>
25 #include <deal.II/fe/fe_values.h>
26 #include <deal.II/meshworker/dof_info.h>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 namespace LocalIntegrators
31 {
39  namespace Elasticity
40  {
47  template <int dim>
48  inline void cell_matrix (
50  const FEValuesBase<dim> &fe,
51  const double factor = 1.)
52  {
53  const unsigned int n_dofs = fe.dofs_per_cell;
54 
55  AssertDimension(fe.get_fe().n_components(), dim);
56  AssertDimension(M.m(), n_dofs);
57  AssertDimension(M.n(), n_dofs);
58 
59  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
60  {
61  const double dx = factor * fe.JxW(k);
62  for (unsigned int i=0; i<n_dofs; ++i)
63  for (unsigned int j=0; j<n_dofs; ++j)
64  for (unsigned int d1=0; d1<dim; ++d1)
65  for (unsigned int d2=0; d2<dim; ++d2)
66  M(i,j) += dx * .25 *
67  (fe.shape_grad_component(j,k,d1)[d2] + fe.shape_grad_component(j,k,d2)[d1]) *
68  (fe.shape_grad_component(i,k,d1)[d2] + fe.shape_grad_component(i,k,d2)[d1]);
69  }
70  }
71 
72 
78  template <int dim, typename number>
79  inline void
81  Vector<number> &result,
82  const FEValuesBase<dim> &fe,
83  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &input,
84  double factor = 1.)
85  {
86  const unsigned int nq = fe.n_quadrature_points;
87  const unsigned int n_dofs = fe.dofs_per_cell;
88  AssertDimension(fe.get_fe().n_components(), dim);
89 
91  Assert(result.size() == n_dofs, ExcDimensionMismatch(result.size(), n_dofs));
92 
93  for (unsigned int k=0; k<nq; ++k)
94  {
95  const double dx = factor * fe.JxW(k);
96  for (unsigned int i=0; i<n_dofs; ++i)
97  for (unsigned int d1=0; d1<dim; ++d1)
98  for (unsigned int d2=0; d2<dim; ++d2)
99  {
100  result(i) += dx * .25 *
101  (input[d1][k][d2] + input[d2][k][d1]) *
102  (fe.shape_grad_component(i,k,d1)[d2] + fe.shape_grad_component(i,k,d2)[d1]);
103  }
104  }
105  }
106 
107 
114  template <int dim>
115  inline void nitsche_matrix (
117  const FEValuesBase<dim> &fe,
118  double penalty,
119  double factor = 1.)
120  {
121  const unsigned int n_dofs = fe.dofs_per_cell;
122 
123  AssertDimension(fe.get_fe().n_components(), dim);
124  AssertDimension(M.m(), n_dofs);
125  AssertDimension(M.n(), n_dofs);
126 
127  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
128  {
129  const double dx = factor * fe.JxW(k);
130  const Tensor<1,dim> n = fe.normal_vector(k);
131  for (unsigned int i=0; i<n_dofs; ++i)
132  for (unsigned int j=0; j<n_dofs; ++j)
133  for (unsigned int d1=0; d1<dim; ++d1)
134  {
135  const double u = fe.shape_value_component(j,k,d1);
136  const double v = fe.shape_value_component(i,k,d1);
137  M(i,j) += dx * 2. * penalty * u * v;
138  for (unsigned int d2=0; d2<dim; ++d2)
139  {
140  // v . nabla u n
141  M(i,j) -= .5*dx* fe.shape_grad_component(j,k,d1)[d2] *n[d2]* v;
142  // v (nabla u)^T n
143  M(i,j) -= .5*dx* fe.shape_grad_component(j,k,d2)[d1] *n[d2]* v;
144  // u nabla v n
145  M(i,j) -= .5*dx* fe.shape_grad_component(i,k,d1)[d2] *n[d2]* u;
146  // u (nabla v)^T n
147  M(i,j) -= .5*dx* fe.shape_grad_component(i,k,d2)[d1] *n[d2]* u;
148  }
149  }
150  }
151  }
152 
169  template <int dim, typename number>
171  Vector<number> &result,
172  const FEValuesBase<dim> &fe,
173  const VectorSlice<const std::vector<std::vector<double> > > &input,
174  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput,
175  const VectorSlice<const std::vector<std::vector<double> > > &data,
176  double penalty,
177  double factor = 1.)
178  {
179  const unsigned int n_dofs = fe.dofs_per_cell;
183 
184  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
185  {
186  const double dx = factor * fe.JxW(k);
187  const Tensor<1,dim> n = fe.normal_vector(k);
188  for (unsigned int i=0; i<n_dofs; ++i)
189  for (unsigned int d1=0; d1<dim; ++d1)
190  {
191  const double u= input[d1][k];
192  const double v= fe.shape_value_component(i,k,d1);
193  const double g= data[d1][k];
194  result(i) += dx * 2.*penalty * (u-g) * v;
195 
196  for (unsigned int d2=0; d2<dim; ++d2)
197  {
198  // v . nabla u n
199  result(i) -= .5*dx* v * Dinput[d1][k][d2] * n[d2];
200  // v . (nabla u)^T n
201  result(i) -= .5*dx* v * Dinput[d2][k][d1] * n[d2];
202  // u nabla v n
203  result(i) -= .5*dx * (u-g) * fe.shape_grad_component(i,k,d1)[d2] * n[d2];
204  // u (nabla v)^T n
205  result(i) -= .5*dx * (u-g) * fe.shape_grad_component(i,k,d2)[d1] * n[d2];
206  }
207  }
208  }
209  }
210 
226  template <int dim, typename number>
228  Vector<number> &result,
229  const FEValuesBase<dim> &fe,
230  const VectorSlice<const std::vector<std::vector<double> > > &input,
231  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput,
232  double penalty,
233  double factor = 1.)
234  {
235  const unsigned int n_dofs = fe.dofs_per_cell;
238 
239  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
240  {
241  const double dx = factor * fe.JxW(k);
242  const Tensor<1,dim> n = fe.normal_vector(k);
243  for (unsigned int i=0; i<n_dofs; ++i)
244  for (unsigned int d1=0; d1<dim; ++d1)
245  {
246  const double u= input[d1][k];
247  const double v= fe.shape_value_component(i,k,d1);
248  result(i) += dx * 2.*penalty * u * v;
249 
250  for (unsigned int d2=0; d2<dim; ++d2)
251  {
252  // v . nabla u n
253  result(i) -= .5*dx* v * Dinput[d1][k][d2] * n[d2];
254  // v . (nabla u)^T n
255  result(i) -= .5*dx* v * Dinput[d2][k][d1] * n[d2];
256  // u nabla v n
257  result(i) -= .5*dx * u * fe.shape_grad_component(i,k,d1)[d2] * n[d2];
258  // u (nabla v)^T n
259  result(i) -= .5*dx * u * fe.shape_grad_component(i,k,d2)[d1] * n[d2];
260  }
261  }
262  }
263  }
264 
268  template <int dim>
269  inline void ip_matrix (
270  FullMatrix<double> &M11,
271  FullMatrix<double> &M12,
272  FullMatrix<double> &M21,
273  FullMatrix<double> &M22,
274  const FEValuesBase<dim> &fe1,
275  const FEValuesBase<dim> &fe2,
276  const double pen,
277  const double int_factor = 1.,
278  const double ext_factor = -1.)
279  {
280  const unsigned int n_dofs = fe1.dofs_per_cell;
281 
282  AssertDimension(fe1.get_fe().n_components(), dim);
283  AssertDimension(fe2.get_fe().n_components(), dim);
284  AssertDimension(M11.m(), n_dofs);
285  AssertDimension(M11.n(), n_dofs);
286  AssertDimension(M12.m(), n_dofs);
287  AssertDimension(M12.n(), n_dofs);
288  AssertDimension(M21.m(), n_dofs);
289  AssertDimension(M21.n(), n_dofs);
290  AssertDimension(M22.m(), n_dofs);
291  AssertDimension(M22.n(), n_dofs);
292 
293  const double nu1 = int_factor;
294  const double nu2 = (ext_factor < 0) ? int_factor : ext_factor;
295  const double penalty = .5 * pen * (nu1 + nu2);
296 
297  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
298  {
299  const double dx = fe1.JxW(k);
300  const Tensor<1,dim> n = fe1.normal_vector(k);
301  for (unsigned int i=0; i<n_dofs; ++i)
302  for (unsigned int j=0; j<n_dofs; ++j)
303  for (unsigned int d1=0; d1<dim; ++d1)
304  {
305  const double u1 = fe1.shape_value_component(j,k,d1);
306  const double u2 = fe2.shape_value_component(j,k,d1);
307  const double v1 = fe1.shape_value_component(i,k,d1);
308  const double v2 = fe2.shape_value_component(i,k,d1);
309 
310  M11(i,j) += dx * penalty * u1*v1;
311  M12(i,j) -= dx * penalty * u2*v1;
312  M21(i,j) -= dx * penalty * u1*v2;
313  M22(i,j) += dx * penalty * u2*v2;
314 
315  for (unsigned int d2=0; d2<dim; ++d2)
316  {
317  // v . nabla u n
318  M11(i,j) -= .25 * dx * nu1 * fe1.shape_grad_component(j,k,d1)[d2] * n[d2] * v1;
319  M12(i,j) -= .25 * dx * nu2 * fe2.shape_grad_component(j,k,d1)[d2] * n[d2] * v1;
320  M21(i,j) += .25 * dx * nu1 * fe1.shape_grad_component(j,k,d1)[d2] * n[d2] * v2;
321  M22(i,j) += .25 * dx * nu2 * fe2.shape_grad_component(j,k,d1)[d2] * n[d2] * v2;
322  // v (nabla u)^T n
323  M11(i,j) -= .25 * dx * nu1 * fe1.shape_grad_component(j,k,d2)[d1] * n[d2] * v1;
324  M12(i,j) -= .25 * dx * nu2 * fe2.shape_grad_component(j,k,d2)[d1] * n[d2] * v1;
325  M21(i,j) += .25 * dx * nu1 * fe1.shape_grad_component(j,k,d2)[d1] * n[d2] * v2;
326  M22(i,j) += .25 * dx * nu2 * fe2.shape_grad_component(j,k,d2)[d1] * n[d2] * v2;
327  // u nabla v n
328  M11(i,j) -= .25 * dx * nu1 * fe1.shape_grad_component(i,k,d1)[d2] * n[d2] * u1;
329  M12(i,j) += .25 * dx * nu1 * fe1.shape_grad_component(i,k,d1)[d2] * n[d2] * u2;
330  M21(i,j) -= .25 * dx * nu2 * fe2.shape_grad_component(i,k,d1)[d2] * n[d2] * u1;
331  M22(i,j) += .25 * dx * nu2 * fe2.shape_grad_component(i,k,d1)[d2] * n[d2] * u2;
332  // u (nabla v)^T n
333  M11(i,j) -= .25 * dx * nu1 * fe1.shape_grad_component(i,k,d2)[d1] * n[d2] * u1;
334  M12(i,j) += .25 * dx * nu1 * fe1.shape_grad_component(i,k,d2)[d1] * n[d2] * u2;
335  M21(i,j) -= .25 * dx * nu2 * fe2.shape_grad_component(i,k,d2)[d1] * n[d2] * u1;
336  M22(i,j) += .25 * dx * nu2 * fe2.shape_grad_component(i,k,d2)[d1] * n[d2] * u2;
337  }
338  }
339  }
340  }
347  template<int dim, typename number>
348  void
350  Vector<number> &result1,
351  Vector<number> &result2,
352  const FEValuesBase<dim> &fe1,
353  const FEValuesBase<dim> &fe2,
354  const VectorSlice<const std::vector<std::vector<double> > > &input1,
355  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput1,
356  const VectorSlice<const std::vector<std::vector<double> > > &input2,
357  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput2,
358  double pen,
359  double int_factor = 1.,
360  double ext_factor = -1.)
361  {
362  const unsigned int n1 = fe1.dofs_per_cell;
363 
364  AssertDimension(fe1.get_fe().n_components(), dim);
365  AssertDimension(fe2.get_fe().n_components(), dim);
370 
371  const double nu1 = int_factor;
372  const double nu2 = (ext_factor < 0) ? int_factor : ext_factor;
373  const double penalty = .5 * pen * (nu1 + nu2);
374 
375 
376  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
377  {
378  const double dx = fe1.JxW(k);
379  const Tensor<1,dim> n = fe1.normal_vector(k);
380 
381  for (unsigned int i=0; i<n1; ++i)
382  for (unsigned int d1=0; d1<dim; ++d1)
383  {
384  const double v1 = fe1.shape_value_component(i,k,d1);
385  const double v2 = fe2.shape_value_component(i,k,d1);
386  const double u1 = input1[d1][k];
387  const double u2 = input2[d1][k];
388 
389  result1(i) += dx * penalty * u1*v1;
390  result1(i) -= dx * penalty * u2*v1;
391  result2(i) -= dx * penalty * u1*v2;
392  result2(i) += dx * penalty * u2*v2;
393 
394  for (unsigned int d2=0; d2<dim; ++d2)
395  {
396  // v . nabla u n
397  result1(i) -= .25*dx* (nu1*Dinput1[d1][k][d2]+nu2*Dinput2[d1][k][d2]) * n[d2] * v1;
398  result2(i) += .25*dx* (nu1*Dinput1[d1][k][d2]+nu2*Dinput2[d1][k][d2]) * n[d2] * v2;
399  // v . (nabla u)^T n
400  result1(i) -= .25*dx* (nu1*Dinput1[d2][k][d1]+nu2*Dinput2[d2][k][d1]) * n[d2] * v1;
401  result2(i) += .25*dx* (nu1*Dinput1[d2][k][d1]+nu2*Dinput2[d2][k][d1]) * n[d2] * v2;
402  // u nabla v n
403  result1(i) -= .25*dx* nu1*fe1.shape_grad_component(i,k,d1)[d2] * n[d2] * (u1-u2);
404  result2(i) -= .25*dx* nu2*fe2.shape_grad_component(i,k,d1)[d2] * n[d2] * (u1-u2);
405  // u (nabla v)^T n
406  result1(i) -= .25*dx* nu1*fe1.shape_grad_component(i,k,d2)[d1] * n[d2] * (u1-u2);
407  result2(i) -= .25*dx* nu2*fe2.shape_grad_component(i,k,d2)[d1] * n[d2] * (u1-u2);
408  }
409  }
410  }
411  }
412  }
413 }
414 
415 DEAL_II_NAMESPACE_CLOSE
416 
417 #endif
size_type m() const
void nitsche_residual_homogeneous(Vector< number > &result, const FEValuesBase< dim > &fe, const VectorSlice< const std::vector< std::vector< double > > > &input, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim > > > > &Dinput, double penalty, double factor=1.)
Definition: elasticity.h:227
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1052
const unsigned int dofs_per_cell
Definition: fe_values.h:1464
#define AssertVectorVectorDimension(vec, dim1, dim2)
Definition: exceptions.h:1063
const FiniteElement< dim, spacedim > & get_fe() const
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
Library of integrals over cells and faces.
Definition: advection.h:30
void cell_residual(Vector< number > &result, const FEValuesBase< dim > &fe, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim > > > > &input, double factor=1.)
Definition: elasticity.h:80
size_type n() const
#define Assert(cond, exc)
Definition: exceptions.h:294
std::size_t size() const
void ip_matrix(FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const double pen, const double int_factor=1., const double ext_factor=-1.)
Definition: elasticity.h:269
const unsigned int n_quadrature_points
Definition: fe_values.h:1457
void nitsche_residual(Vector< number > &result, const FEValuesBase< dim > &fe, const VectorSlice< const std::vector< std::vector< double > > > &input, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim > > > > &Dinput, const VectorSlice< const std::vector< std::vector< double > > > &data, double penalty, double factor=1.)
Definition: elasticity.h:170
void nitsche_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
Definition: elasticity.h:115
double JxW(const unsigned int quadrature_point) const
void ip_residual(Vector< number > &result1, Vector< number > &result2, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const VectorSlice< const std::vector< std::vector< double > > > &input1, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim > > > > &Dinput1, const VectorSlice< const std::vector< std::vector< double > > > &input2, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim > > > > &Dinput2, double pen, double int_factor=1., double ext_factor=-1.)
Definition: elasticity.h:349
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
Definition: elasticity.h:48
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const
Tensor< 1, spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const