Reference documentation for deal.II version 9.1.0-pre
elasticity.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2017 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.md at
12 // the top level directory of deal.II.
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 
22 #include <deal.II/base/exceptions.h>
23 #include <deal.II/base/quadrature.h>
24 
25 #include <deal.II/fe/fe_values.h>
26 #include <deal.II/fe/mapping.h>
27 
28 #include <deal.II/lac/full_matrix.h>
29 
30 #include <deal.II/meshworker/dof_info.h>
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 namespace LocalIntegrators
35 {
43  namespace Elasticity
44  {
51  template <int dim>
52  inline void
54  const FEValuesBase<dim> &fe,
55  const double factor = 1.)
56  {
57  const unsigned int n_dofs = fe.dofs_per_cell;
58 
59  AssertDimension(fe.get_fe().n_components(), dim);
60  AssertDimension(M.m(), n_dofs);
61  AssertDimension(M.n(), n_dofs);
62 
63  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
64  {
65  const double dx = factor * fe.JxW(k);
66  for (unsigned int i = 0; i < n_dofs; ++i)
67  for (unsigned int j = 0; j < n_dofs; ++j)
68  for (unsigned int d1 = 0; d1 < dim; ++d1)
69  for (unsigned int d2 = 0; d2 < dim; ++d2)
70  M(i, j) += dx * .25 *
71  (fe.shape_grad_component(j, k, d1)[d2] +
72  fe.shape_grad_component(j, k, d2)[d1]) *
73  (fe.shape_grad_component(i, k, d1)[d2] +
74  fe.shape_grad_component(i, k, d2)[d1]);
75  }
76  }
77 
78 
84  template <int dim, typename number>
85  inline void
87  Vector<number> & result,
88  const FEValuesBase<dim> & fe,
89  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &input,
90  double factor = 1.)
91  {
92  const unsigned int nq = fe.n_quadrature_points;
93  const unsigned int n_dofs = fe.dofs_per_cell;
94  AssertDimension(fe.get_fe().n_components(), dim);
95 
97  Assert(result.size() == n_dofs,
98  ExcDimensionMismatch(result.size(), n_dofs));
99 
100  for (unsigned int k = 0; k < nq; ++k)
101  {
102  const double dx = factor * fe.JxW(k);
103  for (unsigned int i = 0; i < n_dofs; ++i)
104  for (unsigned int d1 = 0; d1 < dim; ++d1)
105  for (unsigned int d2 = 0; d2 < dim; ++d2)
106  {
107  result(i) += dx * .25 *
108  (input[d1][k][d2] + input[d2][k][d1]) *
109  (fe.shape_grad_component(i, k, d1)[d2] +
110  fe.shape_grad_component(i, k, d2)[d1]);
111  }
112  }
113  }
114 
115 
124  template <int dim>
125  inline void
127  const FEValuesBase<dim> &fe,
128  double penalty,
129  double factor = 1.)
130  {
131  const unsigned int n_dofs = fe.dofs_per_cell;
132 
133  AssertDimension(fe.get_fe().n_components(), dim);
134  AssertDimension(M.m(), n_dofs);
135  AssertDimension(M.n(), n_dofs);
136 
137  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
138  {
139  const double dx = factor * fe.JxW(k);
140  const Tensor<1, dim> n = fe.normal_vector(k);
141  for (unsigned int i = 0; i < n_dofs; ++i)
142  for (unsigned int j = 0; j < n_dofs; ++j)
143  for (unsigned int d1 = 0; d1 < dim; ++d1)
144  {
145  const double u = fe.shape_value_component(j, k, d1);
146  const double v = fe.shape_value_component(i, k, d1);
147  M(i, j) += dx * 2. * penalty * u * v;
148  for (unsigned int d2 = 0; d2 < dim; ++d2)
149  {
150  // v . nabla u n
151  M(i, j) -= .5 * dx *
152  fe.shape_grad_component(j, k, d1)[d2] * n[d2] *
153  v;
154  // v (nabla u)^T n
155  M(i, j) -= .5 * dx *
156  fe.shape_grad_component(j, k, d2)[d1] * n[d2] *
157  v;
158  // u nabla v n
159  M(i, j) -= .5 * dx *
160  fe.shape_grad_component(i, k, d1)[d2] * n[d2] *
161  u;
162  // u (nabla v)^T n
163  M(i, j) -= .5 * dx *
164  fe.shape_grad_component(i, k, d2)[d1] * n[d2] *
165  u;
166  }
167  }
168  }
169  }
170 
179  template <int dim>
180  inline void
182  const FEValuesBase<dim> &fe,
183  double penalty,
184  double factor = 1.)
185  {
186  const unsigned int n_dofs = fe.dofs_per_cell;
187 
188  AssertDimension(fe.get_fe().n_components(), dim);
189  AssertDimension(M.m(), n_dofs);
190  AssertDimension(M.n(), n_dofs);
191 
192  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
193  {
194  const double dx = factor * fe.JxW(k);
195  const Tensor<1, dim> n = fe.normal_vector(k);
196  for (unsigned int i = 0; i < n_dofs; ++i)
197  for (unsigned int j = 0; j < n_dofs; ++j)
198  {
199  double udotn = 0.;
200  double vdotn = 0.;
201  double ngradun = 0.;
202  double ngradvn = 0.;
203 
204  for (unsigned int d = 0; d < dim; ++d)
205  {
206  udotn += n[d] * fe.shape_value_component(j, k, d);
207  vdotn += n[d] * fe.shape_value_component(i, k, d);
208  ngradun += n * fe.shape_grad_component(j, k, d) * n[d];
209  ngradvn += n * fe.shape_grad_component(i, k, d) * n[d];
210  }
211  for (unsigned int d1 = 0; d1 < dim; ++d1)
212  {
213  const double u =
214  fe.shape_value_component(j, k, d1) - udotn * n[d1];
215  const double v =
216  fe.shape_value_component(i, k, d1) - vdotn * n[d1];
217  M(i, j) += dx * 2. * penalty * u * v;
218  // Correct the gradients below and subtract normal component
219  M(i, j) += dx * (ngradun * v + ngradvn * u);
220  for (unsigned int d2 = 0; d2 < dim; ++d2)
221  {
222  // v . nabla u n
223  M(i, j) -= .5 * dx *
224  fe.shape_grad_component(j, k, d1)[d2] *
225  n[d2] * v;
226  // v (nabla u)^T n
227  M(i, j) -= .5 * dx *
228  fe.shape_grad_component(j, k, d2)[d1] *
229  n[d2] * v;
230  // u nabla v n
231  M(i, j) -= .5 * dx *
232  fe.shape_grad_component(i, k, d1)[d2] *
233  n[d2] * u;
234  // u (nabla v)^T n
235  M(i, j) -= .5 * dx *
236  fe.shape_grad_component(i, k, d2)[d1] *
237  n[d2] * u;
238  }
239  }
240  }
241  }
242  }
243 
261  template <int dim, typename number>
262  void
264  Vector<number> & result,
265  const FEValuesBase<dim> & fe,
266  const VectorSlice<const std::vector<std::vector<double>>> & input,
267  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &Dinput,
268  const VectorSlice<const std::vector<std::vector<double>>> & data,
269  double penalty,
270  double factor = 1.)
271  {
272  const unsigned int n_dofs = fe.dofs_per_cell;
276 
277  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
278  {
279  const double dx = factor * fe.JxW(k);
280  const Tensor<1, dim> n = fe.normal_vector(k);
281  for (unsigned int i = 0; i < n_dofs; ++i)
282  for (unsigned int d1 = 0; d1 < dim; ++d1)
283  {
284  const double u = input[d1][k];
285  const double v = fe.shape_value_component(i, k, d1);
286  const double g = data[d1][k];
287  result(i) += dx * 2. * penalty * (u - g) * v;
288 
289  for (unsigned int d2 = 0; d2 < dim; ++d2)
290  {
291  // v . nabla u n
292  result(i) -= .5 * dx * v * Dinput[d1][k][d2] * n[d2];
293  // v . (nabla u)^T n
294  result(i) -= .5 * dx * v * Dinput[d2][k][d1] * n[d2];
295  // u nabla v n
296  result(i) -= .5 * dx * (u - g) *
297  fe.shape_grad_component(i, k, d1)[d2] * n[d2];
298  // u (nabla v)^T n
299  result(i) -= .5 * dx * (u - g) *
300  fe.shape_grad_component(i, k, d2)[d1] * n[d2];
301  }
302  }
303  }
304  }
305 
314  template <int dim, typename number>
315  inline void
317  Vector<number> & result,
318  const FEValuesBase<dim> & fe,
319  const VectorSlice<const std::vector<std::vector<double>>> & input,
320  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &Dinput,
321  const VectorSlice<const std::vector<std::vector<double>>> & data,
322  double penalty,
323  double factor = 1.)
324  {
325  const unsigned int n_dofs = fe.dofs_per_cell;
329 
330  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
331  {
332  const double dx = factor * fe.JxW(k);
333  const Tensor<1, dim> n = fe.normal_vector(k);
334  for (unsigned int i = 0; i < n_dofs; ++i)
335  {
336  double udotn = 0.;
337  double gdotn = 0.;
338  double vdotn = 0.;
339  double ngradun = 0.;
340  double ngradvn = 0.;
341 
342  for (unsigned int d = 0; d < dim; ++d)
343  {
344  udotn += n[d] * input[d][k];
345  gdotn += n[d] * data[d][k];
346  vdotn += n[d] * fe.shape_value_component(i, k, d);
347  ngradun += n * Dinput[d][k] * n[d];
348  ngradvn += n * fe.shape_grad_component(i, k, d) * n[d];
349  }
350  for (unsigned int d1 = 0; d1 < dim; ++d1)
351  {
352  const double u = input[d1][k] - udotn * n[d1];
353  const double v =
354  fe.shape_value_component(i, k, d1) - vdotn * n[d1];
355  const double g = data[d1][k] - gdotn * n[d1];
356  result(i) += dx * 2. * penalty * (u - g) * v;
357  // Correct the gradients below and subtract normal component
358  result(i) += dx * (ngradun * v + ngradvn * (u - g));
359  for (unsigned int d2 = 0; d2 < dim; ++d2)
360  {
361  // v . nabla u n
362  result(i) -= .5 * dx * Dinput[d1][k][d2] * n[d2] * v;
363  // v (nabla u)^T n
364  result(i) -= .5 * dx * Dinput[d2][k][d1] * n[d2] * v;
365  // u nabla v n
366  result(i) -= .5 * dx * (u - g) *
367  fe.shape_grad_component(i, k, d1)[d2] *
368  n[d2];
369  // u (nabla v)^T n
370  result(i) -= .5 * dx * (u - g) *
371  fe.shape_grad_component(i, k, d2)[d1] *
372  n[d2];
373  }
374  }
375  }
376  }
377  }
378 
395  template <int dim, typename number>
396  void
398  Vector<number> & result,
399  const FEValuesBase<dim> & fe,
400  const VectorSlice<const std::vector<std::vector<double>>> & input,
401  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &Dinput,
402  double penalty,
403  double factor = 1.)
404  {
405  const unsigned int n_dofs = fe.dofs_per_cell;
408 
409  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
410  {
411  const double dx = factor * fe.JxW(k);
412  const Tensor<1, dim> n = fe.normal_vector(k);
413  for (unsigned int i = 0; i < n_dofs; ++i)
414  for (unsigned int d1 = 0; d1 < dim; ++d1)
415  {
416  const double u = input[d1][k];
417  const double v = fe.shape_value_component(i, k, d1);
418  result(i) += dx * 2. * penalty * u * v;
419 
420  for (unsigned int d2 = 0; d2 < dim; ++d2)
421  {
422  // v . nabla u n
423  result(i) -= .5 * dx * v * Dinput[d1][k][d2] * n[d2];
424  // v . (nabla u)^T n
425  result(i) -= .5 * dx * v * Dinput[d2][k][d1] * n[d2];
426  // u nabla v n
427  result(i) -= .5 * dx * u *
428  fe.shape_grad_component(i, k, d1)[d2] * n[d2];
429  // u (nabla v)^T n
430  result(i) -= .5 * dx * u *
431  fe.shape_grad_component(i, k, d2)[d1] * n[d2];
432  }
433  }
434  }
435  }
436 
440  template <int dim>
441  inline void
443  FullMatrix<double> & M12,
444  FullMatrix<double> & M21,
445  FullMatrix<double> & M22,
446  const FEValuesBase<dim> &fe1,
447  const FEValuesBase<dim> &fe2,
448  const double pen,
449  const double int_factor = 1.,
450  const double ext_factor = -1.)
451  {
452  const unsigned int n_dofs = fe1.dofs_per_cell;
453 
454  AssertDimension(fe1.get_fe().n_components(), dim);
455  AssertDimension(fe2.get_fe().n_components(), dim);
456  AssertDimension(M11.m(), n_dofs);
457  AssertDimension(M11.n(), n_dofs);
458  AssertDimension(M12.m(), n_dofs);
459  AssertDimension(M12.n(), n_dofs);
460  AssertDimension(M21.m(), n_dofs);
461  AssertDimension(M21.n(), n_dofs);
462  AssertDimension(M22.m(), n_dofs);
463  AssertDimension(M22.n(), n_dofs);
464 
465  const double nu1 = int_factor;
466  const double nu2 = (ext_factor < 0) ? int_factor : ext_factor;
467  const double penalty = .5 * pen * (nu1 + nu2);
468 
469  for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
470  {
471  const double dx = fe1.JxW(k);
472  const Tensor<1, dim> n = fe1.normal_vector(k);
473  for (unsigned int i = 0; i < n_dofs; ++i)
474  for (unsigned int j = 0; j < n_dofs; ++j)
475  for (unsigned int d1 = 0; d1 < dim; ++d1)
476  {
477  const double u1 = fe1.shape_value_component(j, k, d1);
478  const double u2 = fe2.shape_value_component(j, k, d1);
479  const double v1 = fe1.shape_value_component(i, k, d1);
480  const double v2 = fe2.shape_value_component(i, k, d1);
481 
482  M11(i, j) += dx * penalty * u1 * v1;
483  M12(i, j) -= dx * penalty * u2 * v1;
484  M21(i, j) -= dx * penalty * u1 * v2;
485  M22(i, j) += dx * penalty * u2 * v2;
486 
487  for (unsigned int d2 = 0; d2 < dim; ++d2)
488  {
489  // v . nabla u n
490  M11(i, j) -= .25 * dx * nu1 *
491  fe1.shape_grad_component(j, k, d1)[d2] *
492  n[d2] * v1;
493  M12(i, j) -= .25 * dx * nu2 *
494  fe2.shape_grad_component(j, k, d1)[d2] *
495  n[d2] * v1;
496  M21(i, j) += .25 * dx * nu1 *
497  fe1.shape_grad_component(j, k, d1)[d2] *
498  n[d2] * v2;
499  M22(i, j) += .25 * dx * nu2 *
500  fe2.shape_grad_component(j, k, d1)[d2] *
501  n[d2] * v2;
502  // v (nabla u)^T n
503  M11(i, j) -= .25 * dx * nu1 *
504  fe1.shape_grad_component(j, k, d2)[d1] *
505  n[d2] * v1;
506  M12(i, j) -= .25 * dx * nu2 *
507  fe2.shape_grad_component(j, k, d2)[d1] *
508  n[d2] * v1;
509  M21(i, j) += .25 * dx * nu1 *
510  fe1.shape_grad_component(j, k, d2)[d1] *
511  n[d2] * v2;
512  M22(i, j) += .25 * dx * nu2 *
513  fe2.shape_grad_component(j, k, d2)[d1] *
514  n[d2] * v2;
515  // u nabla v n
516  M11(i, j) -= .25 * dx * nu1 *
517  fe1.shape_grad_component(i, k, d1)[d2] *
518  n[d2] * u1;
519  M12(i, j) += .25 * dx * nu1 *
520  fe1.shape_grad_component(i, k, d1)[d2] *
521  n[d2] * u2;
522  M21(i, j) -= .25 * dx * nu2 *
523  fe2.shape_grad_component(i, k, d1)[d2] *
524  n[d2] * u1;
525  M22(i, j) += .25 * dx * nu2 *
526  fe2.shape_grad_component(i, k, d1)[d2] *
527  n[d2] * u2;
528  // u (nabla v)^T n
529  M11(i, j) -= .25 * dx * nu1 *
530  fe1.shape_grad_component(i, k, d2)[d1] *
531  n[d2] * u1;
532  M12(i, j) += .25 * dx * nu1 *
533  fe1.shape_grad_component(i, k, d2)[d1] *
534  n[d2] * u2;
535  M21(i, j) -= .25 * dx * nu2 *
536  fe2.shape_grad_component(i, k, d2)[d1] *
537  n[d2] * u1;
538  M22(i, j) += .25 * dx * nu2 *
539  fe2.shape_grad_component(i, k, d2)[d1] *
540  n[d2] * u2;
541  }
542  }
543  }
544  }
551  template <int dim, typename number>
552  void
554  Vector<number> & result1,
555  Vector<number> & result2,
556  const FEValuesBase<dim> & fe1,
557  const FEValuesBase<dim> & fe2,
558  const VectorSlice<const std::vector<std::vector<double>>> &input1,
559  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>>
560  & Dinput1,
561  const VectorSlice<const std::vector<std::vector<double>>> &input2,
562  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>>
563  & Dinput2,
564  double pen,
565  double int_factor = 1.,
566  double ext_factor = -1.)
567  {
568  const unsigned int n1 = fe1.dofs_per_cell;
569 
570  AssertDimension(fe1.get_fe().n_components(), dim);
571  AssertDimension(fe2.get_fe().n_components(), dim);
576 
577  const double nu1 = int_factor;
578  const double nu2 = (ext_factor < 0) ? int_factor : ext_factor;
579  const double penalty = .5 * pen * (nu1 + nu2);
580 
581 
582  for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
583  {
584  const double dx = fe1.JxW(k);
585  const Tensor<1, dim> n = fe1.normal_vector(k);
586 
587  for (unsigned int i = 0; i < n1; ++i)
588  for (unsigned int d1 = 0; d1 < dim; ++d1)
589  {
590  const double v1 = fe1.shape_value_component(i, k, d1);
591  const double v2 = fe2.shape_value_component(i, k, d1);
592  const double u1 = input1[d1][k];
593  const double u2 = input2[d1][k];
594 
595  result1(i) += dx * penalty * u1 * v1;
596  result1(i) -= dx * penalty * u2 * v1;
597  result2(i) -= dx * penalty * u1 * v2;
598  result2(i) += dx * penalty * u2 * v2;
599 
600  for (unsigned int d2 = 0; d2 < dim; ++d2)
601  {
602  // v . nabla u n
603  result1(i) -=
604  .25 * dx *
605  (nu1 * Dinput1[d1][k][d2] + nu2 * Dinput2[d1][k][d2]) *
606  n[d2] * v1;
607  result2(i) +=
608  .25 * dx *
609  (nu1 * Dinput1[d1][k][d2] + nu2 * Dinput2[d1][k][d2]) *
610  n[d2] * v2;
611  // v . (nabla u)^T n
612  result1(i) -=
613  .25 * dx *
614  (nu1 * Dinput1[d2][k][d1] + nu2 * Dinput2[d2][k][d1]) *
615  n[d2] * v1;
616  result2(i) +=
617  .25 * dx *
618  (nu1 * Dinput1[d2][k][d1] + nu2 * Dinput2[d2][k][d1]) *
619  n[d2] * v2;
620  // u nabla v n
621  result1(i) -= .25 * dx * nu1 *
622  fe1.shape_grad_component(i, k, d1)[d2] *
623  n[d2] * (u1 - u2);
624  result2(i) -= .25 * dx * nu2 *
625  fe2.shape_grad_component(i, k, d1)[d2] *
626  n[d2] * (u1 - u2);
627  // u (nabla v)^T n
628  result1(i) -= .25 * dx * nu1 *
629  fe1.shape_grad_component(i, k, d2)[d1] *
630  n[d2] * (u1 - u2);
631  result2(i) -= .25 * dx * nu2 *
632  fe2.shape_grad_component(i, k, d2)[d1] *
633  n[d2] * (u1 - u2);
634  }
635  }
636  }
637  }
638  } // namespace Elasticity
639 } // namespace LocalIntegrators
640 
641 DEAL_II_NAMESPACE_CLOSE
642 
643 #endif
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1366
const unsigned int dofs_per_cell
Definition: fe_values.h:2033
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const
size_type size() const
#define AssertVectorVectorDimension(vec, dim1, dim2)
Definition: exceptions.h:1377
void nitsche_tangential_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
Definition: elasticity.h:181
Library of integrals over cells and faces.
Definition: advection.h:34
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:397
size_type n() const
#define Assert(cond, exc)
Definition: exceptions.h:1227
static::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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:86
const FiniteElement< dim, spacedim > & get_fe() 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:442
size_type m() const
const unsigned int n_quadrature_points
Definition: fe_values.h:2026
void nitsche_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
Definition: elasticity.h:126
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:263
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:553
Tensor< 1, spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
double JxW(const unsigned int quadrature_point) const
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:316
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
Definition: elasticity.h:53