Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
laplace.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2010 - 2021 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
24
26#include <deal.II/fe/mapping.h>
27
29
31
33
34namespace LocalIntegrators
35{
41 namespace Laplace
42 {
50 template <int dim>
51 void
53 const FEValuesBase<dim> &fe,
54 const double factor = 1.)
55 {
56 const unsigned int n_dofs = fe.dofs_per_cell;
57 const unsigned int n_components = fe.get_fe().n_components();
58
59 for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
60 {
61 const double dx = fe.JxW(k) * factor;
62 for (unsigned int i = 0; i < n_dofs; ++i)
63 {
64 double Mii = 0.0;
65 for (unsigned int d = 0; d < n_components; ++d)
66 Mii += dx * (fe.shape_grad_component(i, k, d) *
67 fe.shape_grad_component(i, k, d));
68
69 M(i, i) += Mii;
70
71 for (unsigned int j = i + 1; j < n_dofs; ++j)
72 {
73 double Mij = 0.0;
74 for (unsigned int d = 0; d < n_components; ++d)
75 Mij += dx * (fe.shape_grad_component(j, k, d) *
76 fe.shape_grad_component(i, k, d));
77
78 M(i, j) += Mij;
79 M(j, i) += Mij;
80 }
81 }
82 }
83 }
84
90 template <int dim>
91 inline void
93 const FEValuesBase<dim> & fe,
94 const std::vector<Tensor<1, dim>> &input,
95 double factor = 1.)
96 {
97 const unsigned int nq = fe.n_quadrature_points;
98 const unsigned int n_dofs = fe.dofs_per_cell;
99 Assert(input.size() == nq, ExcDimensionMismatch(input.size(), nq));
100 Assert(result.size() == n_dofs,
101 ExcDimensionMismatch(result.size(), n_dofs));
102
103 for (unsigned int k = 0; k < nq; ++k)
104 {
105 const double dx = factor * fe.JxW(k);
106 for (unsigned int i = 0; i < n_dofs; ++i)
107 result(i) += dx * (input[k] * fe.shape_grad(i, k));
108 }
109 }
110
111
117 template <int dim>
118 inline void
120 const FEValuesBase<dim> & fe,
121 const ArrayView<const std::vector<Tensor<1, dim>>> &input,
122 double factor = 1.)
123 {
124 const unsigned int nq = fe.n_quadrature_points;
125 const unsigned int n_dofs = fe.dofs_per_cell;
126 const unsigned int n_comp = fe.get_fe().n_components();
127
129 Assert(result.size() == n_dofs,
130 ExcDimensionMismatch(result.size(), n_dofs));
131
132 for (unsigned int k = 0; k < nq; ++k)
133 {
134 const double dx = factor * fe.JxW(k);
135 for (unsigned int i = 0; i < n_dofs; ++i)
136 for (unsigned int d = 0; d < n_comp; ++d)
137 {
138 result(i) +=
139 dx * (input[d][k] * fe.shape_grad_component(i, k, d));
140 }
141 }
142 }
143
144
155 template <int dim>
156 void
158 const FEValuesBase<dim> &fe,
159 double penalty,
160 double factor = 1.)
161 {
162 const unsigned int n_dofs = fe.dofs_per_cell;
163 const unsigned int n_comp = fe.get_fe().n_components();
164
165 Assert(M.m() == n_dofs, ExcDimensionMismatch(M.m(), n_dofs));
166 Assert(M.n() == n_dofs, ExcDimensionMismatch(M.n(), n_dofs));
167
168 for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
169 {
170 const double dx = fe.JxW(k) * factor;
171 const Tensor<1, dim> &n = fe.normal_vector(k);
172 for (unsigned int i = 0; i < n_dofs; ++i)
173 for (unsigned int j = 0; j < n_dofs; ++j)
174 for (unsigned int d = 0; d < n_comp; ++d)
175 M(i, j) += dx * (2. * fe.shape_value_component(i, k, d) *
176 penalty * fe.shape_value_component(j, k, d) -
177 (n * fe.shape_grad_component(i, k, d)) *
178 fe.shape_value_component(j, k, d) -
179 (n * fe.shape_grad_component(j, k, d)) *
180 fe.shape_value_component(i, k, d));
181 }
182 }
183
196 template <int dim>
197 void
199 const FEValuesBase<dim> &fe,
200 double penalty,
201 double factor = 1.)
202 {
203 const unsigned int n_dofs = fe.dofs_per_cell;
205 Assert(M.m() == n_dofs, ExcDimensionMismatch(M.m(), n_dofs));
206 Assert(M.n() == n_dofs, ExcDimensionMismatch(M.n(), n_dofs));
207
208 for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
209 {
210 const double dx = fe.JxW(k) * factor;
211 const Tensor<1, dim> n = fe.normal_vector(k);
212 for (unsigned int i = 0; i < n_dofs; ++i)
213 for (unsigned int j = 0; j < n_dofs; ++j)
214 {
215 double udotn = 0.;
216 double vdotn = 0.;
217 double ngradun = 0.;
218 double ngradvn = 0.;
219
220 for (unsigned int d = 0; d < dim; ++d)
221 {
222 udotn += n[d] * fe.shape_value_component(j, k, d);
223 vdotn += n[d] * fe.shape_value_component(i, k, d);
224 ngradun += n * fe.shape_grad_component(j, k, d) * n[d];
225 ngradvn += n * fe.shape_grad_component(i, k, d) * n[d];
226 }
227
228 for (unsigned int d = 0; d < dim; ++d)
229 {
230 const double v_t =
231 fe.shape_value_component(i, k, d) - vdotn * n[d];
232 const double dnv_t =
233 n * fe.shape_grad_component(i, k, d) - ngradvn * n[d];
234 const double u_t =
235 fe.shape_value_component(j, k, d) - udotn * n[d];
236 const double dnu_t =
237 n * fe.shape_grad_component(j, k, d) - ngradun * n[d];
238
239 M(i, j) += dx * (2. * penalty * u_t * v_t - dnu_t * v_t -
240 dnv_t * u_t);
241 }
242 }
243 }
244 }
245
259 template <int dim>
260 void
262 const FEValuesBase<dim> & fe,
263 const std::vector<double> & input,
264 const std::vector<Tensor<1, dim>> &Dinput,
265 const std::vector<double> & data,
266 double penalty,
267 double factor = 1.)
268 {
269 const unsigned int n_dofs = fe.dofs_per_cell;
270 AssertDimension(input.size(), fe.n_quadrature_points);
271 AssertDimension(Dinput.size(), fe.n_quadrature_points);
272 AssertDimension(data.size(), fe.n_quadrature_points);
273
274 for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
275 {
276 const double dx = factor * fe.JxW(k);
277 const Tensor<1, dim> n = fe.normal_vector(k);
278 for (unsigned int i = 0; i < n_dofs; ++i)
279 {
280 const double dnv = fe.shape_grad(i, k) * n;
281 const double dnu = Dinput[k] * n;
282 const double v = fe.shape_value(i, k);
283 const double u = input[k];
284 const double g = data[k];
285
286 result(i) +=
287 dx * (2. * penalty * (u - g) * v - dnv * (u - g) - dnu * v);
288 }
289 }
290 }
291
306 template <int dim>
307 void
309 const FEValuesBase<dim> & fe,
310 const ArrayView<const std::vector<double>> & input,
311 const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput,
312 const ArrayView<const std::vector<double>> & data,
313 double penalty,
314 double factor = 1.)
315 {
316 const unsigned int n_dofs = fe.dofs_per_cell;
317 const unsigned int n_comp = fe.get_fe().n_components();
321
322 for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
323 {
324 const double dx = factor * fe.JxW(k);
325 const Tensor<1, dim> n = fe.normal_vector(k);
326 for (unsigned int i = 0; i < n_dofs; ++i)
327 for (unsigned int d = 0; d < n_comp; ++d)
328 {
329 const double dnv = fe.shape_grad_component(i, k, d) * n;
330 const double dnu = Dinput[d][k] * n;
331 const double v = fe.shape_value_component(i, k, d);
332 const double u = input[d][k];
333 const double g = data[d][k];
334
335 result(i) +=
336 dx * (2. * penalty * (u - g) * v - dnv * (u - g) - dnu * v);
337 }
338 }
339 }
340
357 template <int dim>
358 void
360 FullMatrix<double> & M12,
361 FullMatrix<double> & M21,
362 FullMatrix<double> & M22,
363 const FEValuesBase<dim> &fe1,
364 const FEValuesBase<dim> &fe2,
365 double penalty,
366 double factor1 = 1.,
367 double factor2 = -1.)
368 {
369 const unsigned int n_dofs = fe1.dofs_per_cell;
370 AssertDimension(M11.n(), n_dofs);
371 AssertDimension(M11.m(), n_dofs);
372 AssertDimension(M12.n(), n_dofs);
373 AssertDimension(M12.m(), n_dofs);
374 AssertDimension(M21.n(), n_dofs);
375 AssertDimension(M21.m(), n_dofs);
376 AssertDimension(M22.n(), n_dofs);
377 AssertDimension(M22.m(), n_dofs);
378
379 const double nui = factor1;
380 const double nue = (factor2 < 0) ? factor1 : factor2;
381 const double nu = .5 * (nui + nue);
382
383 for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
384 {
385 const double dx = fe1.JxW(k);
386 const Tensor<1, dim> &n = fe1.normal_vector(k);
387 for (unsigned int d = 0; d < fe1.get_fe().n_components(); ++d)
388 {
389 for (unsigned int i = 0; i < n_dofs; ++i)
390 {
391 for (unsigned int j = 0; j < n_dofs; ++j)
392 {
393 const double vi = fe1.shape_value_component(i, k, d);
394 const double dnvi = n * fe1.shape_grad_component(i, k, d);
395 const double ve = fe2.shape_value_component(i, k, d);
396 const double dnve = n * fe2.shape_grad_component(i, k, d);
397 const double ui = fe1.shape_value_component(j, k, d);
398 const double dnui = n * fe1.shape_grad_component(j, k, d);
399 const double ue = fe2.shape_value_component(j, k, d);
400 const double dnue = n * fe2.shape_grad_component(j, k, d);
401 M11(i, j) +=
402 dx * (-.5 * nui * dnvi * ui - .5 * nui * dnui * vi +
403 nu * penalty * ui * vi);
404 M12(i, j) +=
405 dx * (.5 * nui * dnvi * ue - .5 * nue * dnue * vi -
406 nu * penalty * vi * ue);
407 M21(i, j) +=
408 dx * (-.5 * nue * dnve * ui + .5 * nui * dnui * ve -
409 nu * penalty * ui * ve);
410 M22(i, j) +=
411 dx * (.5 * nue * dnve * ue + .5 * nue * dnue * ve +
412 nu * penalty * ue * ve);
413 }
414 }
415 }
416 }
417 }
418
430 template <int dim>
431 void
433 FullMatrix<double> & M12,
434 FullMatrix<double> & M21,
435 FullMatrix<double> & M22,
436 const FEValuesBase<dim> &fe1,
437 const FEValuesBase<dim> &fe2,
438 double penalty,
439 double factor1 = 1.,
440 double factor2 = -1.)
441 {
442 const unsigned int n_dofs = fe1.dofs_per_cell;
443 AssertDimension(fe1.get_fe().n_components(), dim);
444 AssertDimension(fe2.get_fe().n_components(), dim);
445 AssertDimension(M11.n(), n_dofs);
446 AssertDimension(M11.m(), n_dofs);
447 AssertDimension(M12.n(), n_dofs);
448 AssertDimension(M12.m(), n_dofs);
449 AssertDimension(M21.n(), n_dofs);
450 AssertDimension(M21.m(), n_dofs);
451 AssertDimension(M22.n(), n_dofs);
452 AssertDimension(M22.m(), n_dofs);
453
454 const double nui = factor1;
455 const double nue = (factor2 < 0) ? factor1 : factor2;
456 const double nu = .5 * (nui + nue);
457
458 for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
459 {
460 const double dx = fe1.JxW(k);
461 const Tensor<1, dim> n = fe1.normal_vector(k);
462 for (unsigned int i = 0; i < n_dofs; ++i)
463 {
464 // We compute the tangential component by subtracting
465 // the normal component from the total. Thus, compute
466 // the normal component first
467 for (unsigned int j = 0; j < n_dofs; ++j)
468 {
469 double u1dotn = 0.;
470 double v1dotn = 0.;
471 double u2dotn = 0.;
472 double v2dotn = 0.;
473
474 double ngradu1n = 0.;
475 double ngradv1n = 0.;
476 double ngradu2n = 0.;
477 double ngradv2n = 0.;
478
479 for (unsigned int d = 0; d < dim; ++d)
480 {
481 u1dotn += n[d] * fe1.shape_value_component(j, k, d);
482 v1dotn += n[d] * fe1.shape_value_component(i, k, d);
483 u2dotn += n[d] * fe2.shape_value_component(j, k, d);
484 v2dotn += n[d] * fe2.shape_value_component(i, k, d);
485
486 ngradu1n += n * fe1.shape_grad_component(j, k, d) * n[d];
487 ngradv1n += n * fe1.shape_grad_component(i, k, d) * n[d];
488 ngradu2n += n * fe2.shape_grad_component(j, k, d) * n[d];
489 ngradv2n += n * fe2.shape_grad_component(i, k, d) * n[d];
490 }
491
492 // The following code is equal to ip_matrix() with
493 // the only exception that all variables introduced
494 // below denote tangential components.
495 for (unsigned int d = 0; d < dim; ++d)
496 {
497 const double vi =
498 fe1.shape_value_component(i, k, d) - v1dotn * n[d];
499 const double dnvi =
500 n * fe1.shape_grad_component(i, k, d) - ngradv1n * n[d];
501
502 const double ve =
503 fe2.shape_value_component(i, k, d) - v2dotn * n[d];
504 const double dnve =
505 n * fe2.shape_grad_component(i, k, d) - ngradv2n * n[d];
506
507 const double ui =
508 fe1.shape_value_component(j, k, d) - u1dotn * n[d];
509 const double dnui =
510 n * fe1.shape_grad_component(j, k, d) - ngradu1n * n[d];
511
512 const double ue =
513 fe2.shape_value_component(j, k, d) - u2dotn * n[d];
514 const double dnue =
515 n * fe2.shape_grad_component(j, k, d) - ngradu2n * n[d];
516
517 M11(i, j) +=
518 dx * (-.5 * nui * dnvi * ui - .5 * nui * dnui * vi +
519 nu * penalty * ui * vi);
520 M12(i, j) +=
521 dx * (.5 * nui * dnvi * ue - .5 * nue * dnue * vi -
522 nu * penalty * vi * ue);
523 M21(i, j) +=
524 dx * (-.5 * nue * dnve * ui + .5 * nui * dnui * ve -
525 nu * penalty * ui * ve);
526 M22(i, j) +=
527 dx * (.5 * nue * dnve * ue + .5 * nue * dnue * ve +
528 nu * penalty * ue * ve);
529 }
530 }
531 }
532 }
533 }
534
542 template <int dim>
543 void
545 Vector<double> & result2,
546 const FEValuesBase<dim> & fe1,
547 const FEValuesBase<dim> & fe2,
548 const std::vector<double> & input1,
549 const std::vector<Tensor<1, dim>> &Dinput1,
550 const std::vector<double> & input2,
551 const std::vector<Tensor<1, dim>> &Dinput2,
552 double pen,
553 double int_factor = 1.,
554 double ext_factor = -1.)
555 {
556 Assert(fe1.get_fe().n_components() == 1,
558 Assert(fe2.get_fe().n_components() == 1,
560
561 const double nui = int_factor;
562 const double nue = (ext_factor < 0) ? int_factor : ext_factor;
563 const double penalty = .5 * pen * (nui + nue);
564
565 const unsigned int n_dofs = fe1.dofs_per_cell;
566
567 for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
568 {
569 const double dx = fe1.JxW(k);
570 const Tensor<1, dim> n = fe1.normal_vector(k);
571
572 for (unsigned int i = 0; i < n_dofs; ++i)
573 {
574 const double vi = fe1.shape_value(i, k);
575 const Tensor<1, dim> &Dvi = fe1.shape_grad(i, k);
576 const double dnvi = Dvi * n;
577 const double ve = fe2.shape_value(i, k);
578 const Tensor<1, dim> &Dve = fe2.shape_grad(i, k);
579 const double dnve = Dve * n;
580
581 const double ui = input1[k];
582 const Tensor<1, dim> &Dui = Dinput1[k];
583 const double dnui = Dui * n;
584 const double ue = input2[k];
585 const Tensor<1, dim> &Due = Dinput2[k];
586 const double dnue = Due * n;
587
588 result1(i) += dx * (-.5 * nui * dnvi * ui - .5 * nui * dnui * vi +
589 penalty * ui * vi);
590 result1(i) += dx * (.5 * nui * dnvi * ue - .5 * nue * dnue * vi -
591 penalty * vi * ue);
592 result2(i) += dx * (-.5 * nue * dnve * ui + .5 * nui * dnui * ve -
593 penalty * ui * ve);
594 result2(i) += dx * (.5 * nue * dnve * ue + .5 * nue * dnue * ve +
595 penalty * ue * ve);
596 }
597 }
598 }
599
600
609 template <int dim>
610 void
612 Vector<double> & result2,
613 const FEValuesBase<dim> & fe1,
614 const FEValuesBase<dim> & fe2,
615 const ArrayView<const std::vector<double>> & input1,
616 const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput1,
617 const ArrayView<const std::vector<double>> & input2,
618 const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput2,
619 double pen,
620 double int_factor = 1.,
621 double ext_factor = -1.)
622 {
623 const unsigned int n_comp = fe1.get_fe().n_components();
624 const unsigned int n1 = fe1.dofs_per_cell;
625
630
631 const double nui = int_factor;
632 const double nue = (ext_factor < 0) ? int_factor : ext_factor;
633 const double penalty = .5 * pen * (nui + nue);
634
635
636 for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
637 {
638 const double dx = fe1.JxW(k);
639 const Tensor<1, dim> n = fe1.normal_vector(k);
640
641 for (unsigned int i = 0; i < n1; ++i)
642 for (unsigned int d = 0; d < n_comp; ++d)
643 {
644 const double vi = fe1.shape_value_component(i, k, d);
645 const Tensor<1, dim> &Dvi = fe1.shape_grad_component(i, k, d);
646 const double dnvi = Dvi * n;
647 const double ve = fe2.shape_value_component(i, k, d);
648 const Tensor<1, dim> &Dve = fe2.shape_grad_component(i, k, d);
649 const double dnve = Dve * n;
650
651 const double ui = input1[d][k];
652 const Tensor<1, dim> &Dui = Dinput1[d][k];
653 const double dnui = Dui * n;
654 const double ue = input2[d][k];
655 const Tensor<1, dim> &Due = Dinput2[d][k];
656 const double dnue = Due * n;
657
658 result1(i) += dx * (-.5 * nui * dnvi * ui -
659 .5 * nui * dnui * vi + penalty * ui * vi);
660 result1(i) += dx * (.5 * nui * dnvi * ue -
661 .5 * nue * dnue * vi - penalty * vi * ue);
662 result2(i) += dx * (-.5 * nue * dnve * ui +
663 .5 * nui * dnui * ve - penalty * ui * ve);
664 result2(i) += dx * (.5 * nue * dnve * ue +
665 .5 * nue * dnue * ve + penalty * ue * ve);
666 }
667 }
668 }
669
670
671
683 template <int dim, int spacedim, typename number>
684 double
687 unsigned int deg1,
688 unsigned int deg2)
689 {
690 const unsigned int normal1 =
692 const unsigned int normal2 =
694 const unsigned int deg1sq = (deg1 == 0) ? 1 : deg1 * (deg1 + 1);
695 const unsigned int deg2sq = (deg2 == 0) ? 1 : deg2 * (deg2 + 1);
696
697 double penalty1 = deg1sq / dinfo1.cell->extent_in_direction(normal1);
698 double penalty2 = deg2sq / dinfo2.cell->extent_in_direction(normal2);
699 if (dinfo1.cell->has_children() ^ dinfo2.cell->has_children())
700 {
701 Assert(dinfo1.face == dinfo2.face, ExcInternalError());
702 Assert(dinfo1.face->has_children(), ExcInternalError());
703 penalty1 *= 2;
704 }
705 const double penalty = 0.5 * (penalty1 + penalty2);
706 return penalty;
707 }
708 } // namespace Laplace
709} // namespace LocalIntegrators
710
711
713
714#endif
const unsigned int dofs_per_cell
Definition fe_values.h:2451
const Tensor< 1, spacedim > & normal_vector(const unsigned int q_point) const
double shape_value_component(const unsigned int i, const unsigned int q_point, const unsigned int component) const
const unsigned int n_quadrature_points
Definition fe_values.h:2433
Tensor< 1, spacedim > shape_grad_component(const unsigned int i, const unsigned int q_point, const unsigned int component) const
const Tensor< 1, spacedim > & shape_grad(const unsigned int i, const unsigned int q_point) const
const FiniteElement< dim, spacedim > & get_fe() const
double JxW(const unsigned int q_point) const
const double & shape_value(const unsigned int i, const unsigned int q_point) const
unsigned int n_components() const
size_type n() const
size_type m() const
Triangulation< dim, spacedim >::face_iterator face
The current face.
Definition dof_info.h:82
unsigned int face_number
Definition dof_info.h:90
Triangulation< dim, spacedim >::cell_iterator cell
The current cell.
Definition dof_info.h:79
size_type size() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
#define Assert(cond, exc)
#define AssertVectorVectorDimension(VEC, DIM1, DIM2)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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:432
void cell_residual(Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< Tensor< 1, dim > > &input, double factor=1.)
Definition laplace.h:92
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
Definition laplace.h:52
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:685
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:261
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:359
void nitsche_tangential_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
Definition laplace.h:198
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:544
void nitsche_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
Definition laplace.h:157
Library of integrals over cells and faces.
Definition advection.h:35