Reference documentation for deal.II version 8.5.1
elasticity.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2016 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 
159  template <int dim>
162  const FEValuesBase<dim> &fe,
163  double penalty,
164  double factor = 1.)
165  {
166  const unsigned int n_dofs = fe.dofs_per_cell;
167 
168  AssertDimension(fe.get_fe().n_components(), dim);
169  AssertDimension(M.m(), n_dofs);
170  AssertDimension(M.n(), n_dofs);
171 
172  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
173  {
174  const double dx = factor * fe.JxW(k);
175  const Tensor<1,dim> n = fe.normal_vector(k);
176  for (unsigned int i=0; i<n_dofs; ++i)
177  for (unsigned int j=0; j<n_dofs; ++j)
178  {
179  double udotn = 0.;
180  double vdotn = 0.;
181  double ngradun = 0.;
182  double ngradvn = 0.;
183 
184  for (unsigned int d=0; d<dim; ++d)
185  {
186  udotn += n[d]*fe.shape_value_component(j,k,d);
187  vdotn += n[d]*fe.shape_value_component(i,k,d);
188  ngradun += n*fe.shape_grad_component(j,k,d)*n[d];
189  ngradvn += n*fe.shape_grad_component(i,k,d)*n[d];
190  }
191  for (unsigned int d1=0; d1<dim; ++d1)
192  {
193  const double u = fe.shape_value_component(j,k,d1) - udotn * n[d1];
194  const double v = fe.shape_value_component(i,k,d1) - vdotn * n[d1];
195  M(i,j) += dx * 2. * penalty * u * v;
196  // Correct the gradients below and subtract normal component
197  M(i,j) += dx * (ngradun * v + ngradvn * u);
198  for (unsigned int d2=0; d2<dim; ++d2)
199  {
200  // v . nabla u n
201  M(i,j) -= .5*dx* fe.shape_grad_component(j,k,d1)[d2] *n[d2]* v;
202  // v (nabla u)^T n
203  M(i,j) -= .5*dx* fe.shape_grad_component(j,k,d2)[d1] *n[d2]* v;
204  // u nabla v n
205  M(i,j) -= .5*dx* fe.shape_grad_component(i,k,d1)[d2] *n[d2]* u;
206  // u (nabla v)^T n
207  M(i,j) -= .5*dx* fe.shape_grad_component(i,k,d2)[d1] *n[d2]* u;
208  }
209  }
210  }
211  }
212  }
213 
230  template <int dim, typename number>
232  Vector<number> &result,
233  const FEValuesBase<dim> &fe,
234  const VectorSlice<const std::vector<std::vector<double> > > &input,
235  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput,
236  const VectorSlice<const std::vector<std::vector<double> > > &data,
237  double penalty,
238  double factor = 1.)
239  {
240  const unsigned int n_dofs = fe.dofs_per_cell;
244 
245  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
246  {
247  const double dx = factor * fe.JxW(k);
248  const Tensor<1,dim> n = fe.normal_vector(k);
249  for (unsigned int i=0; i<n_dofs; ++i)
250  for (unsigned int d1=0; d1<dim; ++d1)
251  {
252  const double u= input[d1][k];
253  const double v= fe.shape_value_component(i,k,d1);
254  const double g= data[d1][k];
255  result(i) += dx * 2.*penalty * (u-g) * v;
256 
257  for (unsigned int d2=0; d2<dim; ++d2)
258  {
259  // v . nabla u n
260  result(i) -= .5*dx* v * Dinput[d1][k][d2] * n[d2];
261  // v . (nabla u)^T n
262  result(i) -= .5*dx* v * Dinput[d2][k][d1] * n[d2];
263  // u nabla v n
264  result(i) -= .5*dx * (u-g) * fe.shape_grad_component(i,k,d1)[d2] * n[d2];
265  // u (nabla v)^T n
266  result(i) -= .5*dx * (u-g) * fe.shape_grad_component(i,k,d2)[d1] * n[d2];
267  }
268  }
269  }
270  }
271 
278  template <int dim, typename number>
280  Vector<number> &result,
281  const FEValuesBase<dim> &fe,
282  const VectorSlice<const std::vector<std::vector<double> > > &input,
283  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput,
284  const VectorSlice<const std::vector<std::vector<double> > > &data,
285  double penalty,
286  double factor = 1.)
287  {
288  const unsigned int n_dofs = fe.dofs_per_cell;
292 
293  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
294  {
295  const double dx = factor * fe.JxW(k);
296  const Tensor<1,dim> n = fe.normal_vector(k);
297  for (unsigned int i=0; i<n_dofs; ++i)
298  {
299  double udotn = 0.;
300  double gdotn = 0.;
301  double vdotn = 0.;
302  double ngradun = 0.;
303  double ngradvn = 0.;
304 
305  for (unsigned int d=0; d<dim; ++d)
306  {
307  udotn += n[d]*input[d][k];
308  gdotn += n[d]*data[d][k];
309  vdotn += n[d]*fe.shape_value_component(i,k,d);
310  ngradun += n*Dinput[d][k]*n[d];
311  ngradvn += n*fe.shape_grad_component(i,k,d)*n[d];
312  }
313  for (unsigned int d1=0; d1<dim; ++d1)
314  {
315  const double u= input[d1][k] - udotn*n[d1];
316  const double v= fe.shape_value_component(i,k,d1) - vdotn*n[d1];
317  const double g= data[d1][k] - gdotn*n[d1];
318  result(i) += dx * 2.*penalty * (u-g) * v;
319  // Correct the gradients below and subtract normal component
320  result(i) += dx * (ngradun * v + ngradvn * (u-g));
321  for (unsigned int d2=0; d2<dim; ++d2)
322  {
323 
324  // v . nabla u n
325  result(i) -= .5*dx* Dinput[d1][k][d2] *n[d2]* v;
326  // v (nabla u)^T n
327  result(i) -= .5*dx* Dinput[d2][k][d1] *n[d2]* v;
328  // u nabla v n
329  result(i) -= .5*dx* (u-g) * fe.shape_grad_component(i,k,d1)[d2] *n[d2];
330  // u (nabla v)^T n
331  result(i) -= .5*dx* (u-g) * fe.shape_grad_component(i,k,d2)[d1] *n[d2];
332  }
333  }
334  }
335  }
336  }
337 
353  template <int dim, typename number>
355  Vector<number> &result,
356  const FEValuesBase<dim> &fe,
357  const VectorSlice<const std::vector<std::vector<double> > > &input,
358  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput,
359  double penalty,
360  double factor = 1.)
361  {
362  const unsigned int n_dofs = fe.dofs_per_cell;
365 
366  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
367  {
368  const double dx = factor * fe.JxW(k);
369  const Tensor<1,dim> n = fe.normal_vector(k);
370  for (unsigned int i=0; i<n_dofs; ++i)
371  for (unsigned int d1=0; d1<dim; ++d1)
372  {
373  const double u= input[d1][k];
374  const double v= fe.shape_value_component(i,k,d1);
375  result(i) += dx * 2.*penalty * u * v;
376 
377  for (unsigned int d2=0; d2<dim; ++d2)
378  {
379  // v . nabla u n
380  result(i) -= .5*dx* v * Dinput[d1][k][d2] * n[d2];
381  // v . (nabla u)^T n
382  result(i) -= .5*dx* v * Dinput[d2][k][d1] * n[d2];
383  // u nabla v n
384  result(i) -= .5*dx * u * fe.shape_grad_component(i,k,d1)[d2] * n[d2];
385  // u (nabla v)^T n
386  result(i) -= .5*dx * u * fe.shape_grad_component(i,k,d2)[d1] * n[d2];
387  }
388  }
389  }
390  }
391 
395  template <int dim>
396  inline void ip_matrix (
397  FullMatrix<double> &M11,
398  FullMatrix<double> &M12,
399  FullMatrix<double> &M21,
400  FullMatrix<double> &M22,
401  const FEValuesBase<dim> &fe1,
402  const FEValuesBase<dim> &fe2,
403  const double pen,
404  const double int_factor = 1.,
405  const double ext_factor = -1.)
406  {
407  const unsigned int n_dofs = fe1.dofs_per_cell;
408 
409  AssertDimension(fe1.get_fe().n_components(), dim);
410  AssertDimension(fe2.get_fe().n_components(), dim);
411  AssertDimension(M11.m(), n_dofs);
412  AssertDimension(M11.n(), n_dofs);
413  AssertDimension(M12.m(), n_dofs);
414  AssertDimension(M12.n(), n_dofs);
415  AssertDimension(M21.m(), n_dofs);
416  AssertDimension(M21.n(), n_dofs);
417  AssertDimension(M22.m(), n_dofs);
418  AssertDimension(M22.n(), n_dofs);
419 
420  const double nu1 = int_factor;
421  const double nu2 = (ext_factor < 0) ? int_factor : ext_factor;
422  const double penalty = .5 * pen * (nu1 + nu2);
423 
424  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
425  {
426  const double dx = fe1.JxW(k);
427  const Tensor<1,dim> n = fe1.normal_vector(k);
428  for (unsigned int i=0; i<n_dofs; ++i)
429  for (unsigned int j=0; j<n_dofs; ++j)
430  for (unsigned int d1=0; d1<dim; ++d1)
431  {
432  const double u1 = fe1.shape_value_component(j,k,d1);
433  const double u2 = fe2.shape_value_component(j,k,d1);
434  const double v1 = fe1.shape_value_component(i,k,d1);
435  const double v2 = fe2.shape_value_component(i,k,d1);
436 
437  M11(i,j) += dx * penalty * u1*v1;
438  M12(i,j) -= dx * penalty * u2*v1;
439  M21(i,j) -= dx * penalty * u1*v2;
440  M22(i,j) += dx * penalty * u2*v2;
441 
442  for (unsigned int d2=0; d2<dim; ++d2)
443  {
444  // v . nabla u n
445  M11(i,j) -= .25 * dx * nu1 * fe1.shape_grad_component(j,k,d1)[d2] * n[d2] * v1;
446  M12(i,j) -= .25 * dx * nu2 * fe2.shape_grad_component(j,k,d1)[d2] * n[d2] * v1;
447  M21(i,j) += .25 * dx * nu1 * fe1.shape_grad_component(j,k,d1)[d2] * n[d2] * v2;
448  M22(i,j) += .25 * dx * nu2 * fe2.shape_grad_component(j,k,d1)[d2] * n[d2] * v2;
449  // v (nabla u)^T n
450  M11(i,j) -= .25 * dx * nu1 * fe1.shape_grad_component(j,k,d2)[d1] * n[d2] * v1;
451  M12(i,j) -= .25 * dx * nu2 * fe2.shape_grad_component(j,k,d2)[d1] * n[d2] * v1;
452  M21(i,j) += .25 * dx * nu1 * fe1.shape_grad_component(j,k,d2)[d1] * n[d2] * v2;
453  M22(i,j) += .25 * dx * nu2 * fe2.shape_grad_component(j,k,d2)[d1] * n[d2] * v2;
454  // u nabla v n
455  M11(i,j) -= .25 * dx * nu1 * fe1.shape_grad_component(i,k,d1)[d2] * n[d2] * u1;
456  M12(i,j) += .25 * dx * nu1 * fe1.shape_grad_component(i,k,d1)[d2] * n[d2] * u2;
457  M21(i,j) -= .25 * dx * nu2 * fe2.shape_grad_component(i,k,d1)[d2] * n[d2] * u1;
458  M22(i,j) += .25 * dx * nu2 * fe2.shape_grad_component(i,k,d1)[d2] * n[d2] * u2;
459  // u (nabla v)^T n
460  M11(i,j) -= .25 * dx * nu1 * fe1.shape_grad_component(i,k,d2)[d1] * n[d2] * u1;
461  M12(i,j) += .25 * dx * nu1 * fe1.shape_grad_component(i,k,d2)[d1] * n[d2] * u2;
462  M21(i,j) -= .25 * dx * nu2 * fe2.shape_grad_component(i,k,d2)[d1] * n[d2] * u1;
463  M22(i,j) += .25 * dx * nu2 * fe2.shape_grad_component(i,k,d2)[d1] * n[d2] * u2;
464  }
465  }
466  }
467  }
474  template<int dim, typename number>
475  void
477  Vector<number> &result1,
478  Vector<number> &result2,
479  const FEValuesBase<dim> &fe1,
480  const FEValuesBase<dim> &fe2,
481  const VectorSlice<const std::vector<std::vector<double> > > &input1,
482  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput1,
483  const VectorSlice<const std::vector<std::vector<double> > > &input2,
484  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput2,
485  double pen,
486  double int_factor = 1.,
487  double ext_factor = -1.)
488  {
489  const unsigned int n1 = fe1.dofs_per_cell;
490 
491  AssertDimension(fe1.get_fe().n_components(), dim);
492  AssertDimension(fe2.get_fe().n_components(), dim);
497 
498  const double nu1 = int_factor;
499  const double nu2 = (ext_factor < 0) ? int_factor : ext_factor;
500  const double penalty = .5 * pen * (nu1 + nu2);
501 
502 
503  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
504  {
505  const double dx = fe1.JxW(k);
506  const Tensor<1,dim> n = fe1.normal_vector(k);
507 
508  for (unsigned int i=0; i<n1; ++i)
509  for (unsigned int d1=0; d1<dim; ++d1)
510  {
511  const double v1 = fe1.shape_value_component(i,k,d1);
512  const double v2 = fe2.shape_value_component(i,k,d1);
513  const double u1 = input1[d1][k];
514  const double u2 = input2[d1][k];
515 
516  result1(i) += dx * penalty * u1*v1;
517  result1(i) -= dx * penalty * u2*v1;
518  result2(i) -= dx * penalty * u1*v2;
519  result2(i) += dx * penalty * u2*v2;
520 
521  for (unsigned int d2=0; d2<dim; ++d2)
522  {
523  // v . nabla u n
524  result1(i) -= .25*dx* (nu1*Dinput1[d1][k][d2]+nu2*Dinput2[d1][k][d2]) * n[d2] * v1;
525  result2(i) += .25*dx* (nu1*Dinput1[d1][k][d2]+nu2*Dinput2[d1][k][d2]) * n[d2] * v2;
526  // v . (nabla u)^T n
527  result1(i) -= .25*dx* (nu1*Dinput1[d2][k][d1]+nu2*Dinput2[d2][k][d1]) * n[d2] * v1;
528  result2(i) += .25*dx* (nu1*Dinput1[d2][k][d1]+nu2*Dinput2[d2][k][d1]) * n[d2] * v2;
529  // u nabla v n
530  result1(i) -= .25*dx* nu1*fe1.shape_grad_component(i,k,d1)[d2] * n[d2] * (u1-u2);
531  result2(i) -= .25*dx* nu2*fe2.shape_grad_component(i,k,d1)[d2] * n[d2] * (u1-u2);
532  // u (nabla v)^T n
533  result1(i) -= .25*dx* nu1*fe1.shape_grad_component(i,k,d2)[d1] * n[d2] * (u1-u2);
534  result2(i) -= .25*dx* nu2*fe2.shape_grad_component(i,k,d2)[d1] * n[d2] * (u1-u2);
535  }
536  }
537  }
538  }
539  }
540 }
541 
542 DEAL_II_NAMESPACE_CLOSE
543 
544 #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:354
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
const unsigned int dofs_per_cell
Definition: fe_values.h:1459
#define AssertVectorVectorDimension(vec, dim1, dim2)
Definition: exceptions.h:1157
const FiniteElement< dim, spacedim > & get_fe() const
void nitsche_tangential_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
Definition: elasticity.h:160
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:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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:396
const unsigned int n_quadrature_points
Definition: fe_values.h:1452
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:231
void nitsche_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
Definition: elasticity.h:115
void nitsche_tangential_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:279
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:476
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