Reference documentation for deal.II version 9.1.0-pre
laplace.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_laplace_h
17 #define dealii_integrators_laplace_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 Laplace
44  {
55  template <int dim>
56  void
58  const FEValuesBase<dim> &fe,
59  const double factor = 1.)
60  {
61  const unsigned int n_dofs = fe.dofs_per_cell;
62  const unsigned int n_components = fe.get_fe().n_components();
63 
64  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
65  {
66  const double dx = fe.JxW(k) * factor;
67  for (unsigned int i = 0; i < n_dofs; ++i)
68  {
69  double Mii = 0.0;
70  for (unsigned int d = 0; d < n_components; ++d)
71  Mii += dx * (fe.shape_grad_component(i, k, d) *
72  fe.shape_grad_component(i, k, d));
73 
74  M(i, i) += Mii;
75 
76  for (unsigned int j = i + 1; j < n_dofs; ++j)
77  {
78  double Mij = 0.0;
79  for (unsigned int d = 0; d < n_components; ++d)
80  Mij += dx * (fe.shape_grad_component(j, k, d) *
81  fe.shape_grad_component(i, k, d));
82 
83  M(i, j) += Mij;
84  M(j, i) += Mij;
85  }
86  }
87  }
88  }
89 
95  template <int dim>
96  inline void
98  const FEValuesBase<dim> & fe,
99  const std::vector<Tensor<1, dim>> &input,
100  double factor = 1.)
101  {
102  const unsigned int nq = fe.n_quadrature_points;
103  const unsigned int n_dofs = fe.dofs_per_cell;
104  Assert(input.size() == nq, ExcDimensionMismatch(input.size(), nq));
105  Assert(result.size() == n_dofs,
106  ExcDimensionMismatch(result.size(), n_dofs));
107 
108  for (unsigned int k = 0; k < nq; ++k)
109  {
110  const double dx = factor * fe.JxW(k);
111  for (unsigned int i = 0; i < n_dofs; ++i)
112  result(i) += dx * (input[k] * fe.shape_grad(i, k));
113  }
114  }
115 
116 
122  template <int dim>
123  inline void
125  Vector<double> & result,
126  const FEValuesBase<dim> & fe,
127  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &input,
128  double factor = 1.)
129  {
130  const unsigned int nq = fe.n_quadrature_points;
131  const unsigned int n_dofs = fe.dofs_per_cell;
132  const unsigned int n_comp = fe.get_fe().n_components();
133 
135  Assert(result.size() == n_dofs,
136  ExcDimensionMismatch(result.size(), n_dofs));
137 
138  for (unsigned int k = 0; k < nq; ++k)
139  {
140  const double dx = factor * fe.JxW(k);
141  for (unsigned int i = 0; i < n_dofs; ++i)
142  for (unsigned int d = 0; d < n_comp; ++d)
143  {
144  result(i) +=
145  dx * (input[d][k] * fe.shape_grad_component(i, k, d));
146  }
147  }
148  }
149 
150 
164  template <int dim>
165  void
167  const FEValuesBase<dim> &fe,
168  double penalty,
169  double factor = 1.)
170  {
171  const unsigned int n_dofs = fe.dofs_per_cell;
172  const unsigned int n_comp = fe.get_fe().n_components();
173 
174  Assert(M.m() == n_dofs, ExcDimensionMismatch(M.m(), n_dofs));
175  Assert(M.n() == n_dofs, ExcDimensionMismatch(M.n(), n_dofs));
176 
177  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
178  {
179  const double dx = fe.JxW(k) * factor;
180  const Tensor<1, dim> n = fe.normal_vector(k);
181  for (unsigned int i = 0; i < n_dofs; ++i)
182  for (unsigned int j = 0; j < n_dofs; ++j)
183  for (unsigned int d = 0; d < n_comp; ++d)
184  M(i, j) += dx * (2. * fe.shape_value_component(i, k, d) *
185  penalty * fe.shape_value_component(j, k, d) -
186  (n * fe.shape_grad_component(i, k, d)) *
187  fe.shape_value_component(j, k, d) -
188  (n * fe.shape_grad_component(j, k, d)) *
189  fe.shape_value_component(i, k, d));
190  }
191  }
192 
208  template <int dim>
209  void
211  const FEValuesBase<dim> &fe,
212  double penalty,
213  double factor = 1.)
214  {
215  const unsigned int n_dofs = fe.dofs_per_cell;
216  AssertDimension(fe.get_fe().n_components(), dim);
217  Assert(M.m() == n_dofs, ExcDimensionMismatch(M.m(), n_dofs));
218  Assert(M.n() == n_dofs, ExcDimensionMismatch(M.n(), n_dofs));
219 
220  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
221  {
222  const double dx = fe.JxW(k) * factor;
223  const Tensor<1, dim> n = fe.normal_vector(k);
224  for (unsigned int i = 0; i < n_dofs; ++i)
225  for (unsigned int j = 0; j < n_dofs; ++j)
226  {
227  double udotn = 0.;
228  double vdotn = 0.;
229  double ngradun = 0.;
230  double ngradvn = 0.;
231 
232  for (unsigned int d = 0; d < dim; ++d)
233  {
234  udotn += n[d] * fe.shape_value_component(j, k, d);
235  vdotn += n[d] * fe.shape_value_component(i, k, d);
236  ngradun += n * fe.shape_grad_component(j, k, d) * n[d];
237  ngradvn += n * fe.shape_grad_component(i, k, d) * n[d];
238  }
239 
240  for (unsigned int d = 0; d < dim; ++d)
241  {
242  const double v_t =
243  fe.shape_value_component(i, k, d) - vdotn * n[d];
244  const double dnv_t =
245  n * fe.shape_grad_component(i, k, d) - ngradvn * n[d];
246  const double u_t =
247  fe.shape_value_component(j, k, d) - udotn * n[d];
248  const double dnu_t =
249  n * fe.shape_grad_component(j, k, d) - ngradun * n[d];
250 
251  M(i, j) += dx * (2. * penalty * u_t * v_t - dnu_t * v_t -
252  dnv_t * u_t);
253  }
254  }
255  }
256  }
257 
274  template <int dim>
275  void
277  const FEValuesBase<dim> & fe,
278  const std::vector<double> & input,
279  const std::vector<Tensor<1, dim>> &Dinput,
280  const std::vector<double> & data,
281  double penalty,
282  double factor = 1.)
283  {
284  const unsigned int n_dofs = fe.dofs_per_cell;
285  AssertDimension(input.size(), fe.n_quadrature_points);
286  AssertDimension(Dinput.size(), fe.n_quadrature_points);
287  AssertDimension(data.size(), fe.n_quadrature_points);
288 
289  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
290  {
291  const double dx = factor * fe.JxW(k);
292  const Tensor<1, dim> n = fe.normal_vector(k);
293  for (unsigned int i = 0; i < n_dofs; ++i)
294  {
295  const double dnv = fe.shape_grad(i, k) * n;
296  const double dnu = Dinput[k] * n;
297  const double v = fe.shape_value(i, k);
298  const double u = input[k];
299  const double g = data[k];
300 
301  result(i) +=
302  dx * (2. * penalty * (u - g) * v - dnv * (u - g) - dnu * v);
303  }
304  }
305  }
306 
324  template <int dim>
325  void
327  Vector<double> & result,
328  const FEValuesBase<dim> & fe,
329  const VectorSlice<const std::vector<std::vector<double>>> & input,
330  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &Dinput,
331  const VectorSlice<const std::vector<std::vector<double>>> & data,
332  double penalty,
333  double factor = 1.)
334  {
335  const unsigned int n_dofs = fe.dofs_per_cell;
336  const unsigned int n_comp = fe.get_fe().n_components();
340 
341  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
342  {
343  const double dx = factor * fe.JxW(k);
344  const Tensor<1, dim> n = fe.normal_vector(k);
345  for (unsigned int i = 0; i < n_dofs; ++i)
346  for (unsigned int d = 0; d < n_comp; ++d)
347  {
348  const double dnv = fe.shape_grad_component(i, k, d) * n;
349  const double dnu = Dinput[d][k] * n;
350  const double v = fe.shape_value_component(i, k, d);
351  const double u = input[d][k];
352  const double g = data[d][k];
353 
354  result(i) +=
355  dx * (2. * penalty * (u - g) * v - dnv * (u - g) - dnu * v);
356  }
357  }
358  }
359 
379  template <int dim>
380  void
382  FullMatrix<double> & M12,
383  FullMatrix<double> & M21,
384  FullMatrix<double> & M22,
385  const FEValuesBase<dim> &fe1,
386  const FEValuesBase<dim> &fe2,
387  double penalty,
388  double factor1 = 1.,
389  double factor2 = -1.)
390  {
391  const unsigned int n_dofs = fe1.dofs_per_cell;
392  AssertDimension(M11.n(), n_dofs);
393  AssertDimension(M11.m(), n_dofs);
394  AssertDimension(M12.n(), n_dofs);
395  AssertDimension(M12.m(), n_dofs);
396  AssertDimension(M21.n(), n_dofs);
397  AssertDimension(M21.m(), n_dofs);
398  AssertDimension(M22.n(), n_dofs);
399  AssertDimension(M22.m(), n_dofs);
400 
401  const double nui = factor1;
402  const double nue = (factor2 < 0) ? factor1 : factor2;
403  const double nu = .5 * (nui + nue);
404 
405  for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
406  {
407  const double dx = fe1.JxW(k);
408  const Tensor<1, dim> n = fe1.normal_vector(k);
409  for (unsigned int d = 0; d < fe1.get_fe().n_components(); ++d)
410  {
411  for (unsigned int i = 0; i < n_dofs; ++i)
412  {
413  for (unsigned int j = 0; j < n_dofs; ++j)
414  {
415  const double vi = fe1.shape_value_component(i, k, d);
416  const double dnvi = n * fe1.shape_grad_component(i, k, d);
417  const double ve = fe2.shape_value_component(i, k, d);
418  const double dnve = n * fe2.shape_grad_component(i, k, d);
419  const double ui = fe1.shape_value_component(j, k, d);
420  const double dnui = n * fe1.shape_grad_component(j, k, d);
421  const double ue = fe2.shape_value_component(j, k, d);
422  const double dnue = n * fe2.shape_grad_component(j, k, d);
423  M11(i, j) +=
424  dx * (-.5 * nui * dnvi * ui - .5 * nui * dnui * vi +
425  nu * penalty * ui * vi);
426  M12(i, j) +=
427  dx * (.5 * nui * dnvi * ue - .5 * nue * dnue * vi -
428  nu * penalty * vi * ue);
429  M21(i, j) +=
430  dx * (-.5 * nue * dnve * ui + .5 * nui * dnui * ve -
431  nu * penalty * ui * ve);
432  M22(i, j) +=
433  dx * (.5 * nue * dnve * ue + .5 * nue * dnue * ve +
434  nu * penalty * ue * ve);
435  }
436  }
437  }
438  }
439  }
440 
455  template <int dim>
456  void
458  FullMatrix<double> & M12,
459  FullMatrix<double> & M21,
460  FullMatrix<double> & M22,
461  const FEValuesBase<dim> &fe1,
462  const FEValuesBase<dim> &fe2,
463  double penalty,
464  double factor1 = 1.,
465  double factor2 = -1.)
466  {
467  const unsigned int n_dofs = fe1.dofs_per_cell;
468  AssertDimension(fe1.get_fe().n_components(), dim);
469  AssertDimension(fe2.get_fe().n_components(), dim);
470  AssertDimension(M11.n(), n_dofs);
471  AssertDimension(M11.m(), n_dofs);
472  AssertDimension(M12.n(), n_dofs);
473  AssertDimension(M12.m(), n_dofs);
474  AssertDimension(M21.n(), n_dofs);
475  AssertDimension(M21.m(), n_dofs);
476  AssertDimension(M22.n(), n_dofs);
477  AssertDimension(M22.m(), n_dofs);
478 
479  const double nui = factor1;
480  const double nue = (factor2 < 0) ? factor1 : factor2;
481  const double nu = .5 * (nui + nue);
482 
483  for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
484  {
485  const double dx = fe1.JxW(k);
486  const Tensor<1, dim> n = fe1.normal_vector(k);
487  for (unsigned int i = 0; i < n_dofs; ++i)
488  {
489  // We compute the tangential component by subtracting
490  // the normal component from the total. Thus, compute
491  // the normal component first
492  for (unsigned int j = 0; j < n_dofs; ++j)
493  {
494  double u1dotn = 0.;
495  double v1dotn = 0.;
496  double u2dotn = 0.;
497  double v2dotn = 0.;
498 
499  double ngradu1n = 0.;
500  double ngradv1n = 0.;
501  double ngradu2n = 0.;
502  double ngradv2n = 0.;
503 
504  for (unsigned int d = 0; d < dim; ++d)
505  {
506  u1dotn += n[d] * fe1.shape_value_component(j, k, d);
507  v1dotn += n[d] * fe1.shape_value_component(i, k, d);
508  u2dotn += n[d] * fe2.shape_value_component(j, k, d);
509  v2dotn += n[d] * fe2.shape_value_component(i, k, d);
510 
511  ngradu1n += n * fe1.shape_grad_component(j, k, d) * n[d];
512  ngradv1n += n * fe1.shape_grad_component(i, k, d) * n[d];
513  ngradu2n += n * fe2.shape_grad_component(j, k, d) * n[d];
514  ngradv2n += n * fe2.shape_grad_component(i, k, d) * n[d];
515  }
516 
517  // The following code is equal to ip_matrix() with
518  // the only exception that all variables introduced
519  // below denote tangential components.
520  for (unsigned int d = 0; d < dim; ++d)
521  {
522  const double vi =
523  fe1.shape_value_component(i, k, d) - v1dotn * n[d];
524  const double dnvi =
525  n * fe1.shape_grad_component(i, k, d) - ngradv1n * n[d];
526 
527  const double ve =
528  fe2.shape_value_component(i, k, d) - v2dotn * n[d];
529  const double dnve =
530  n * fe2.shape_grad_component(i, k, d) - ngradv2n * n[d];
531 
532  const double ui =
533  fe1.shape_value_component(j, k, d) - u1dotn * n[d];
534  const double dnui =
535  n * fe1.shape_grad_component(j, k, d) - ngradu1n * n[d];
536 
537  const double ue =
538  fe2.shape_value_component(j, k, d) - u2dotn * n[d];
539  const double dnue =
540  n * fe2.shape_grad_component(j, k, d) - ngradu2n * n[d];
541 
542  M11(i, j) +=
543  dx * (-.5 * nui * dnvi * ui - .5 * nui * dnui * vi +
544  nu * penalty * ui * vi);
545  M12(i, j) +=
546  dx * (.5 * nui * dnvi * ue - .5 * nue * dnue * vi -
547  nu * penalty * vi * ue);
548  M21(i, j) +=
549  dx * (-.5 * nue * dnve * ui + .5 * nui * dnui * ve -
550  nu * penalty * ui * ve);
551  M22(i, j) +=
552  dx * (.5 * nue * dnve * ue + .5 * nue * dnue * ve +
553  nu * penalty * ue * ve);
554  }
555  }
556  }
557  }
558  }
559 
570  template <int dim>
571  void
573  Vector<double> & result2,
574  const FEValuesBase<dim> & fe1,
575  const FEValuesBase<dim> & fe2,
576  const std::vector<double> & input1,
577  const std::vector<Tensor<1, dim>> &Dinput1,
578  const std::vector<double> & input2,
579  const std::vector<Tensor<1, dim>> &Dinput2,
580  double pen,
581  double int_factor = 1.,
582  double ext_factor = -1.)
583  {
584  Assert(fe1.get_fe().n_components() == 1,
585  ExcDimensionMismatch(fe1.get_fe().n_components(), 1));
586  Assert(fe2.get_fe().n_components() == 1,
587  ExcDimensionMismatch(fe2.get_fe().n_components(), 1));
588 
589  const double nui = int_factor;
590  const double nue = (ext_factor < 0) ? int_factor : ext_factor;
591  const double penalty = .5 * pen * (nui + nue);
592 
593  const unsigned int n_dofs = fe1.dofs_per_cell;
594 
595  for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
596  {
597  const double dx = fe1.JxW(k);
598  const Tensor<1, dim> n = fe1.normal_vector(k);
599 
600  for (unsigned int i = 0; i < n_dofs; ++i)
601  {
602  const double vi = fe1.shape_value(i, k);
603  const Tensor<1, dim> &Dvi = fe1.shape_grad(i, k);
604  const double dnvi = Dvi * n;
605  const double ve = fe2.shape_value(i, k);
606  const Tensor<1, dim> &Dve = fe2.shape_grad(i, k);
607  const double dnve = Dve * n;
608 
609  const double ui = input1[k];
610  const Tensor<1, dim> &Dui = Dinput1[k];
611  const double dnui = Dui * n;
612  const double ue = input2[k];
613  const Tensor<1, dim> &Due = Dinput2[k];
614  const double dnue = Due * n;
615 
616  result1(i) += dx * (-.5 * nui * dnvi * ui - .5 * nui * dnui * vi +
617  penalty * ui * vi);
618  result1(i) += dx * (.5 * nui * dnvi * ue - .5 * nue * dnue * vi -
619  penalty * vi * ue);
620  result2(i) += dx * (-.5 * nue * dnve * ui + .5 * nui * dnui * ve -
621  penalty * ui * ve);
622  result2(i) += dx * (.5 * nue * dnve * ue + .5 * nue * dnue * ve +
623  penalty * ue * ve);
624  }
625  }
626  }
627 
628 
640  template <int dim>
641  void
643  Vector<double> & result1,
644  Vector<double> & result2,
645  const FEValuesBase<dim> & fe1,
646  const FEValuesBase<dim> & fe2,
647  const VectorSlice<const std::vector<std::vector<double>>> &input1,
648  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>>
649  & Dinput1,
650  const VectorSlice<const std::vector<std::vector<double>>> &input2,
651  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>>
652  & Dinput2,
653  double pen,
654  double int_factor = 1.,
655  double ext_factor = -1.)
656  {
657  const unsigned int n_comp = fe1.get_fe().n_components();
658  const unsigned int n1 = fe1.dofs_per_cell;
659 
661  AssertVectorVectorDimension(Dinput1, n_comp, fe1.n_quadrature_points);
663  AssertVectorVectorDimension(Dinput2, n_comp, fe2.n_quadrature_points);
664 
665  const double nui = int_factor;
666  const double nue = (ext_factor < 0) ? int_factor : ext_factor;
667  const double penalty = .5 * pen * (nui + nue);
668 
669 
670  for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
671  {
672  const double dx = fe1.JxW(k);
673  const Tensor<1, dim> n = fe1.normal_vector(k);
674 
675  for (unsigned int i = 0; i < n1; ++i)
676  for (unsigned int d = 0; d < n_comp; ++d)
677  {
678  const double vi = fe1.shape_value_component(i, k, d);
679  const Tensor<1, dim> &Dvi = fe1.shape_grad_component(i, k, d);
680  const double dnvi = Dvi * n;
681  const double ve = fe2.shape_value_component(i, k, d);
682  const Tensor<1, dim> &Dve = fe2.shape_grad_component(i, k, d);
683  const double dnve = Dve * n;
684 
685  const double ui = input1[d][k];
686  const Tensor<1, dim> &Dui = Dinput1[d][k];
687  const double dnui = Dui * n;
688  const double ue = input2[d][k];
689  const Tensor<1, dim> &Due = Dinput2[d][k];
690  const double dnue = Due * n;
691 
692  result1(i) += dx * (-.5 * nui * dnvi * ui -
693  .5 * nui * dnui * vi + penalty * ui * vi);
694  result1(i) += dx * (.5 * nui * dnvi * ue -
695  .5 * nue * dnue * vi - penalty * vi * ue);
696  result2(i) += dx * (-.5 * nue * dnve * ui +
697  .5 * nui * dnui * ve - penalty * ui * ve);
698  result2(i) += dx * (.5 * nue * dnve * ue +
699  .5 * nue * dnue * ve + penalty * ue * ve);
700  }
701  }
702  }
703 
704 
705 
720  template <int dim, int spacedim, typename number>
721  double
724  unsigned int deg1,
725  unsigned int deg2)
726  {
727  const unsigned int normal1 =
729  const unsigned int normal2 =
731  const unsigned int deg1sq = (deg1 == 0) ? 1 : deg1 * (deg1 + 1);
732  const unsigned int deg2sq = (deg2 == 0) ? 1 : deg2 * (deg2 + 1);
733 
734  double penalty1 = deg1sq / dinfo1.cell->extent_in_direction(normal1);
735  double penalty2 = deg2sq / dinfo2.cell->extent_in_direction(normal2);
736  if (dinfo1.cell->has_children() ^ dinfo2.cell->has_children())
737  {
738  Assert(dinfo1.face == dinfo2.face, ExcInternalError());
739  Assert(dinfo1.face->has_children(), ExcInternalError());
740  penalty1 *= 2;
741  }
742  const double penalty = 0.5 * (penalty1 + penalty2);
743  return penalty;
744  }
745  } // namespace Laplace
746 } // namespace LocalIntegrators
747 
748 
749 DEAL_II_NAMESPACE_CLOSE
750 
751 #endif
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
void cell_residual(Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< Tensor< 1, dim >> &input, double factor=1.)
Definition: laplace.h:97
const Tensor< 1, spacedim > & shape_grad(const unsigned int function_no, const unsigned int quadrature_point) const
Triangulation< dim, spacedim >::cell_iterator cell
The current cell.
Definition: dof_info.h:78
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1366
void ip_matrix(FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, double penalty, double factor1=1., double factor2=-1.)
Definition: laplace.h:381
void nitsche_tangential_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
Definition: laplace.h:210
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
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
Definition: laplace.h:57
#define AssertVectorVectorDimension(vec, dim1, dim2)
Definition: exceptions.h:1377
Library of integrals over cells and faces.
Definition: advection.h:34
size_type n() const
#define Assert(cond, exc)
Definition: exceptions.h:1227
static::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
const FiniteElement< dim, spacedim > & get_fe() const
size_type m() const
void ip_residual(Vector< double > &result1, Vector< double > &result2, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const std::vector< double > &input1, const std::vector< Tensor< 1, dim >> &Dinput1, const std::vector< double > &input2, const std::vector< Tensor< 1, dim >> &Dinput2, double pen, double int_factor=1., double ext_factor=-1.)
Definition: laplace.h:572
const unsigned int n_quadrature_points
Definition: fe_values.h:2026
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
void nitsche_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
Definition: laplace.h:166
double compute_penalty(const MeshWorker::DoFInfo< dim, spacedim, number > &dinfo1, const MeshWorker::DoFInfo< dim, spacedim, number > &dinfo2, unsigned int deg1, unsigned int deg2)
Definition: laplace.h:722
void nitsche_residual(Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< double > &input, const std::vector< Tensor< 1, dim >> &Dinput, const std::vector< double > &data, double penalty, double factor=1.)
Definition: laplace.h:276
unsigned int face_number
Definition: dof_info.h:89
Triangulation< dim, spacedim >::face_iterator face
The current face.
Definition: dof_info.h:81
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 ip_tangential_matrix(FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, double penalty, double factor1=1., double factor2=-1.)
Definition: laplace.h:457
static::ExceptionBase & ExcInternalError()