deal.II version GIT relicensing-2248-g6a8f632045 2024-12-12 07:30:00+00:00
\(\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// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2010 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_integrators_laplace_h
16#define dealii_integrators_laplace_h
17
18
19#include <deal.II/base/config.h>
20
23
25#include <deal.II/fe/mapping.h>
26
28
30
32
33namespace LocalIntegrators
34{
40 namespace Laplace
41 {
49 template <int dim>
52 const FEValuesBase<dim> &fe,
53 const double factor = 1.)
54 {
55 const unsigned int n_dofs = fe.dofs_per_cell;
56 const unsigned int n_components = fe.get_fe().n_components();
57
58 for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
59 {
60 const double dx = fe.JxW(k) * factor;
61 for (unsigned int i = 0; i < n_dofs; ++i)
62 {
63 double Mii = 0.0;
64 for (unsigned int d = 0; d < n_components; ++d)
65 Mii += dx * (fe.shape_grad_component(i, k, d) *
66 fe.shape_grad_component(i, k, d));
67
68 M(i, i) += Mii;
69
70 for (unsigned int j = i + 1; j < n_dofs; ++j)
71 {
72 double Mij = 0.0;
73 for (unsigned int d = 0; d < n_components; ++d)
74 Mij += dx * (fe.shape_grad_component(j, k, d) *
75 fe.shape_grad_component(i, k, d));
76
77 M(i, j) += Mij;
78 M(j, i) += Mij;
79 }
80 }
81 }
82 }
83
89 template <int dim>
90 DEAL_II_DEPRECATED inline void
92 const FEValuesBase<dim> &fe,
93 const std::vector<Tensor<1, dim>> &input,
94 double factor = 1.)
95 {
96 const unsigned int nq = fe.n_quadrature_points;
97 const unsigned int n_dofs = fe.dofs_per_cell;
98 Assert(input.size() == nq, ExcDimensionMismatch(input.size(), nq));
99 Assert(result.size() == n_dofs,
100 ExcDimensionMismatch(result.size(), n_dofs));
101
102 for (unsigned int k = 0; k < nq; ++k)
103 {
104 const double dx = factor * fe.JxW(k);
105 for (unsigned int i = 0; i < n_dofs; ++i)
106 result(i) += dx * (input[k] * fe.shape_grad(i, k));
107 }
108 }
109
110
116 template <int dim>
117 DEAL_II_DEPRECATED inline void
119 const FEValuesBase<dim> &fe,
120 const ArrayView<const std::vector<Tensor<1, dim>>> &input,
121 double factor = 1.)
122 {
123 const unsigned int nq = fe.n_quadrature_points;
124 const unsigned int n_dofs = fe.dofs_per_cell;
125 const unsigned int n_comp = fe.get_fe().n_components();
126
128 Assert(result.size() == n_dofs,
129 ExcDimensionMismatch(result.size(), n_dofs));
130
131 for (unsigned int k = 0; k < nq; ++k)
132 {
133 const double dx = factor * fe.JxW(k);
134 for (unsigned int i = 0; i < n_dofs; ++i)
135 for (unsigned int d = 0; d < n_comp; ++d)
136 {
137 result(i) +=
138 dx * (input[d][k] * fe.shape_grad_component(i, k, d));
139 }
140 }
141 }
142
143
154 template <int dim>
157 const FEValuesBase<dim> &fe,
158 double penalty,
159 double factor = 1.)
160 {
161 const unsigned int n_dofs = fe.dofs_per_cell;
162 const unsigned int n_comp = fe.get_fe().n_components();
163
164 Assert(M.m() == n_dofs, ExcDimensionMismatch(M.m(), n_dofs));
165 Assert(M.n() == n_dofs, ExcDimensionMismatch(M.n(), n_dofs));
166
167 for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
168 {
169 const double dx = fe.JxW(k) * factor;
170 const Tensor<1, dim> &n = fe.normal_vector(k);
171 for (unsigned int i = 0; i < n_dofs; ++i)
172 for (unsigned int j = 0; j < n_dofs; ++j)
173 for (unsigned int d = 0; d < n_comp; ++d)
174 M(i, j) += dx * (2. * fe.shape_value_component(i, k, d) *
175 penalty * fe.shape_value_component(j, k, d) -
176 (n * fe.shape_grad_component(i, k, d)) *
177 fe.shape_value_component(j, k, d) -
178 (n * fe.shape_grad_component(j, k, d)) *
179 fe.shape_value_component(i, k, d));
180 }
181 }
182
195 template <int dim>
198 const FEValuesBase<dim> &fe,
199 double penalty,
200 double factor = 1.)
201 {
202 const unsigned int n_dofs = fe.dofs_per_cell;
204 Assert(M.m() == n_dofs, ExcDimensionMismatch(M.m(), n_dofs));
205 Assert(M.n() == n_dofs, ExcDimensionMismatch(M.n(), n_dofs));
206
207 for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
208 {
209 const double dx = fe.JxW(k) * factor;
210 const Tensor<1, dim> n = fe.normal_vector(k);
211 for (unsigned int i = 0; i < n_dofs; ++i)
212 for (unsigned int j = 0; j < n_dofs; ++j)
213 {
214 double udotn = 0.;
215 double vdotn = 0.;
216 double ngradun = 0.;
217 double ngradvn = 0.;
218
219 for (unsigned int d = 0; d < dim; ++d)
220 {
221 udotn += n[d] * fe.shape_value_component(j, k, d);
222 vdotn += n[d] * fe.shape_value_component(i, k, d);
223 ngradun += n * fe.shape_grad_component(j, k, d) * n[d];
224 ngradvn += n * fe.shape_grad_component(i, k, d) * n[d];
225 }
226
227 for (unsigned int d = 0; d < dim; ++d)
228 {
229 const double v_t =
230 fe.shape_value_component(i, k, d) - vdotn * n[d];
231 const double dnv_t =
232 n * fe.shape_grad_component(i, k, d) - ngradvn * n[d];
233 const double u_t =
234 fe.shape_value_component(j, k, d) - udotn * n[d];
235 const double dnu_t =
236 n * fe.shape_grad_component(j, k, d) - ngradun * n[d];
237
238 M(i, j) += dx * (2. * penalty * u_t * v_t - dnu_t * v_t -
239 dnv_t * u_t);
240 }
241 }
242 }
243 }
244
258 template <int dim>
261 const FEValuesBase<dim> &fe,
262 const std::vector<double> &input,
263 const std::vector<Tensor<1, dim>> &Dinput,
264 const std::vector<double> &data,
265 double penalty,
266 double factor = 1.)
267 {
268 const unsigned int n_dofs = fe.dofs_per_cell;
269 AssertDimension(input.size(), fe.n_quadrature_points);
270 AssertDimension(Dinput.size(), fe.n_quadrature_points);
272
273 for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
274 {
275 const double dx = factor * fe.JxW(k);
276 const Tensor<1, dim> n = fe.normal_vector(k);
277 for (unsigned int i = 0; i < n_dofs; ++i)
278 {
279 const double dnv = fe.shape_grad(i, k) * n;
280 const double dnu = Dinput[k] * n;
281 const double v = fe.shape_value(i, k);
282 const double u = input[k];
283 const double g = data[k];
284
285 result(i) +=
286 dx * (2. * penalty * (u - g) * v - dnv * (u - g) - dnu * v);
287 }
288 }
289 }
290
305 template <int dim>
308 const FEValuesBase<dim> &fe,
309 const ArrayView<const std::vector<double>> &input,
310 const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput,
311 const ArrayView<const std::vector<double>> &data,
312 double penalty,
313 double factor = 1.)
314 {
315 const unsigned int n_dofs = fe.dofs_per_cell;
316 const unsigned int n_comp = fe.get_fe().n_components();
320
321 for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
322 {
323 const double dx = factor * fe.JxW(k);
324 const Tensor<1, dim> n = fe.normal_vector(k);
325 for (unsigned int i = 0; i < n_dofs; ++i)
326 for (unsigned int d = 0; d < n_comp; ++d)
327 {
328 const double dnv = fe.shape_grad_component(i, k, d) * n;
329 const double dnu = Dinput[d][k] * n;
330 const double v = fe.shape_value_component(i, k, d);
331 const double u = input[d][k];
332 const double g = data[d][k];
333
334 result(i) +=
335 dx * (2. * penalty * (u - g) * v - dnv * (u - g) - dnu * v);
336 }
337 }
338 }
339
356 template <int dim>
362 const FEValuesBase<dim> &fe1,
363 const FEValuesBase<dim> &fe2,
364 double penalty,
365 double factor1 = 1.,
366 double factor2 = -1.)
367 {
368 const unsigned int n_dofs = fe1.dofs_per_cell;
369 AssertDimension(M11.n(), n_dofs);
370 AssertDimension(M11.m(), n_dofs);
371 AssertDimension(M12.n(), n_dofs);
372 AssertDimension(M12.m(), n_dofs);
373 AssertDimension(M21.n(), n_dofs);
374 AssertDimension(M21.m(), n_dofs);
375 AssertDimension(M22.n(), n_dofs);
376 AssertDimension(M22.m(), n_dofs);
377
378 const double nui = factor1;
379 const double nue = (factor2 < 0) ? factor1 : factor2;
380 const double nu = .5 * (nui + nue);
381
382 for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
383 {
384 const double dx = fe1.JxW(k);
385 const Tensor<1, dim> &n = fe1.normal_vector(k);
386 for (unsigned int d = 0; d < fe1.get_fe().n_components(); ++d)
387 {
388 for (unsigned int i = 0; i < n_dofs; ++i)
389 {
390 for (unsigned int j = 0; j < n_dofs; ++j)
391 {
392 const double vi = fe1.shape_value_component(i, k, d);
393 const double dnvi = n * fe1.shape_grad_component(i, k, d);
394 const double ve = fe2.shape_value_component(i, k, d);
395 const double dnve = n * fe2.shape_grad_component(i, k, d);
396 const double ui = fe1.shape_value_component(j, k, d);
397 const double dnui = n * fe1.shape_grad_component(j, k, d);
398 const double ue = fe2.shape_value_component(j, k, d);
399 const double dnue = n * fe2.shape_grad_component(j, k, d);
400 M11(i, j) +=
401 dx * (-.5 * nui * dnvi * ui - .5 * nui * dnui * vi +
402 nu * penalty * ui * vi);
403 M12(i, j) +=
404 dx * (.5 * nui * dnvi * ue - .5 * nue * dnue * vi -
405 nu * penalty * vi * ue);
406 M21(i, j) +=
407 dx * (-.5 * nue * dnve * ui + .5 * nui * dnui * ve -
408 nu * penalty * ui * ve);
409 M22(i, j) +=
410 dx * (.5 * nue * dnve * ue + .5 * nue * dnue * ve +
411 nu * penalty * ue * ve);
412 }
413 }
414 }
415 }
416 }
417
429 template <int dim>
435 const FEValuesBase<dim> &fe1,
436 const FEValuesBase<dim> &fe2,
437 double penalty,
438 double factor1 = 1.,
439 double factor2 = -1.)
440 {
441 const unsigned int n_dofs = fe1.dofs_per_cell;
442 AssertDimension(fe1.get_fe().n_components(), dim);
443 AssertDimension(fe2.get_fe().n_components(), dim);
444 AssertDimension(M11.n(), n_dofs);
445 AssertDimension(M11.m(), n_dofs);
446 AssertDimension(M12.n(), n_dofs);
447 AssertDimension(M12.m(), n_dofs);
448 AssertDimension(M21.n(), n_dofs);
449 AssertDimension(M21.m(), n_dofs);
450 AssertDimension(M22.n(), n_dofs);
451 AssertDimension(M22.m(), n_dofs);
452
453 const double nui = factor1;
454 const double nue = (factor2 < 0) ? factor1 : factor2;
455 const double nu = .5 * (nui + nue);
456
457 for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
458 {
459 const double dx = fe1.JxW(k);
460 const Tensor<1, dim> n = fe1.normal_vector(k);
461 for (unsigned int i = 0; i < n_dofs; ++i)
462 {
463 // We compute the tangential component by subtracting
464 // the normal component from the total. Thus, compute
465 // the normal component first
466 for (unsigned int j = 0; j < n_dofs; ++j)
467 {
468 double u1dotn = 0.;
469 double v1dotn = 0.;
470 double u2dotn = 0.;
471 double v2dotn = 0.;
472
473 double ngradu1n = 0.;
474 double ngradv1n = 0.;
475 double ngradu2n = 0.;
476 double ngradv2n = 0.;
477
478 for (unsigned int d = 0; d < dim; ++d)
479 {
480 u1dotn += n[d] * fe1.shape_value_component(j, k, d);
481 v1dotn += n[d] * fe1.shape_value_component(i, k, d);
482 u2dotn += n[d] * fe2.shape_value_component(j, k, d);
483 v2dotn += n[d] * fe2.shape_value_component(i, k, d);
484
485 ngradu1n += n * fe1.shape_grad_component(j, k, d) * n[d];
486 ngradv1n += n * fe1.shape_grad_component(i, k, d) * n[d];
487 ngradu2n += n * fe2.shape_grad_component(j, k, d) * n[d];
488 ngradv2n += n * fe2.shape_grad_component(i, k, d) * n[d];
489 }
490
491 // The following code is equal to ip_matrix() with
492 // the only exception that all variables introduced
493 // below denote tangential components.
494 for (unsigned int d = 0; d < dim; ++d)
495 {
496 const double vi =
497 fe1.shape_value_component(i, k, d) - v1dotn * n[d];
498 const double dnvi =
499 n * fe1.shape_grad_component(i, k, d) - ngradv1n * n[d];
500
501 const double ve =
502 fe2.shape_value_component(i, k, d) - v2dotn * n[d];
503 const double dnve =
504 n * fe2.shape_grad_component(i, k, d) - ngradv2n * n[d];
505
506 const double ui =
507 fe1.shape_value_component(j, k, d) - u1dotn * n[d];
508 const double dnui =
509 n * fe1.shape_grad_component(j, k, d) - ngradu1n * n[d];
510
511 const double ue =
512 fe2.shape_value_component(j, k, d) - u2dotn * n[d];
513 const double dnue =
514 n * fe2.shape_grad_component(j, k, d) - ngradu2n * n[d];
515
516 M11(i, j) +=
517 dx * (-.5 * nui * dnvi * ui - .5 * nui * dnui * vi +
518 nu * penalty * ui * vi);
519 M12(i, j) +=
520 dx * (.5 * nui * dnvi * ue - .5 * nue * dnue * vi -
521 nu * penalty * vi * ue);
522 M21(i, j) +=
523 dx * (-.5 * nue * dnve * ui + .5 * nui * dnui * ve -
524 nu * penalty * ui * ve);
525 M22(i, j) +=
526 dx * (.5 * nue * dnve * ue + .5 * nue * dnue * ve +
527 nu * penalty * ue * ve);
528 }
529 }
530 }
531 }
532 }
533
541 template <int dim>
544 Vector<double> &result2,
545 const FEValuesBase<dim> &fe1,
546 const FEValuesBase<dim> &fe2,
547 const std::vector<double> &input1,
548 const std::vector<Tensor<1, dim>> &Dinput1,
549 const std::vector<double> &input2,
550 const std::vector<Tensor<1, dim>> &Dinput2,
551 double pen,
552 double int_factor = 1.,
553 double ext_factor = -1.)
554 {
555 Assert(fe1.get_fe().n_components() == 1,
557 Assert(fe2.get_fe().n_components() == 1,
559
560 const double nui = int_factor;
561 const double nue = (ext_factor < 0) ? int_factor : ext_factor;
562 const double penalty = .5 * pen * (nui + nue);
563
564 const unsigned int n_dofs = fe1.dofs_per_cell;
565
566 for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
567 {
568 const double dx = fe1.JxW(k);
569 const Tensor<1, dim> n = fe1.normal_vector(k);
570
571 for (unsigned int i = 0; i < n_dofs; ++i)
572 {
573 const double vi = fe1.shape_value(i, k);
574 const Tensor<1, dim> &Dvi = fe1.shape_grad(i, k);
575 const double dnvi = Dvi * n;
576 const double ve = fe2.shape_value(i, k);
577 const Tensor<1, dim> &Dve = fe2.shape_grad(i, k);
578 const double dnve = Dve * n;
579
580 const double ui = input1[k];
581 const Tensor<1, dim> &Dui = Dinput1[k];
582 const double dnui = Dui * n;
583 const double ue = input2[k];
584 const Tensor<1, dim> &Due = Dinput2[k];
585 const double dnue = Due * n;
586
587 result1(i) += dx * (-.5 * nui * dnvi * ui - .5 * nui * dnui * vi +
588 penalty * ui * vi);
589 result1(i) += dx * (.5 * nui * dnvi * ue - .5 * nue * dnue * vi -
590 penalty * vi * ue);
591 result2(i) += dx * (-.5 * nue * dnve * ui + .5 * nui * dnui * ve -
592 penalty * ui * ve);
593 result2(i) += dx * (.5 * nue * dnve * ue + .5 * nue * dnue * ve +
594 penalty * ue * ve);
595 }
596 }
597 }
598
599
608 template <int dim>
611 Vector<double> &result2,
612 const FEValuesBase<dim> &fe1,
613 const FEValuesBase<dim> &fe2,
614 const ArrayView<const std::vector<double>> &input1,
615 const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput1,
616 const ArrayView<const std::vector<double>> &input2,
617 const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput2,
618 double pen,
619 double int_factor = 1.,
620 double ext_factor = -1.)
621 {
622 const unsigned int n_comp = fe1.get_fe().n_components();
623 const unsigned int n1 = fe1.dofs_per_cell;
624
629
630 const double nui = int_factor;
631 const double nue = (ext_factor < 0) ? int_factor : ext_factor;
632 const double penalty = .5 * pen * (nui + nue);
633
634
635 for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
636 {
637 const double dx = fe1.JxW(k);
638 const Tensor<1, dim> n = fe1.normal_vector(k);
639
640 for (unsigned int i = 0; i < n1; ++i)
641 for (unsigned int d = 0; d < n_comp; ++d)
642 {
643 const double vi = fe1.shape_value_component(i, k, d);
644 const Tensor<1, dim> &Dvi = fe1.shape_grad_component(i, k, d);
645 const double dnvi = Dvi * n;
646 const double ve = fe2.shape_value_component(i, k, d);
647 const Tensor<1, dim> &Dve = fe2.shape_grad_component(i, k, d);
648 const double dnve = Dve * n;
649
650 const double ui = input1[d][k];
651 const Tensor<1, dim> &Dui = Dinput1[d][k];
652 const double dnui = Dui * n;
653 const double ue = input2[d][k];
654 const Tensor<1, dim> &Due = Dinput2[d][k];
655 const double dnue = Due * n;
656
657 result1(i) += dx * (-.5 * nui * dnvi * ui -
658 .5 * nui * dnui * vi + penalty * ui * vi);
659 result1(i) += dx * (.5 * nui * dnvi * ue -
660 .5 * nue * dnue * vi - penalty * vi * ue);
661 result2(i) += dx * (-.5 * nue * dnve * ui +
662 .5 * nui * dnui * ve - penalty * ui * ve);
663 result2(i) += dx * (.5 * nue * dnve * ue +
664 .5 * nue * dnue * ve + penalty * ue * ve);
665 }
666 }
667 }
668
669
670
682 template <int dim, int spacedim, typename number>
683 DEAL_II_DEPRECATED double
686 unsigned int deg1,
687 unsigned int deg2)
688 {
689 const unsigned int normal1 =
691 const unsigned int normal2 =
693 const unsigned int deg1sq = (deg1 == 0) ? 1 : deg1 * (deg1 + 1);
694 const unsigned int deg2sq = (deg2 == 0) ? 1 : deg2 * (deg2 + 1);
695
696 double penalty1 = deg1sq / dinfo1.cell->extent_in_direction(normal1);
697 double penalty2 = deg2sq / dinfo2.cell->extent_in_direction(normal2);
698 if (dinfo1.cell->has_children() ^ dinfo2.cell->has_children())
699 {
700 Assert(dinfo1.face == dinfo2.face, ExcInternalError());
701 Assert(dinfo1.face->has_children(), ExcInternalError());
702 penalty1 *= 2;
703 }
704 const double penalty = 0.5 * (penalty1 + penalty2);
705 return penalty;
706 }
707 } // namespace Laplace
708} // namespace LocalIntegrators
709
710
712
713#endif
const unsigned int dofs_per_cell
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
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:81
unsigned int face_number
Definition dof_info.h:89
Triangulation< dim, spacedim >::cell_iterator cell
The current cell.
Definition dof_info.h:78
virtual size_type size() const override
#define DEAL_II_DEPRECATED
Definition config.h:205
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
#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)
std::vector< index_type > data
Definition mpi.cc:735
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:431
void cell_residual(Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< Tensor< 1, dim > > &input, double factor=1.)
Definition laplace.h:91
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
Definition laplace.h:51
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:684
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:260
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:358
void nitsche_tangential_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
Definition laplace.h:197
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:543
void nitsche_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
Definition laplace.h:156
Library of integrals over cells and faces.
Definition advection.h:34