Reference documentation for deal.II version 8.5.1
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 at
12 // the top level of the deal.II distribution.
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 #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 Laplace
40  {
51  template<int dim>
52  void cell_matrix (
54  const FEValuesBase<dim> &fe,
55  const double factor = 1.)
56  {
57  const unsigned int n_dofs = fe.dofs_per_cell;
58  const unsigned int n_components = fe.get_fe().n_components();
59 
60  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
61  {
62  const double dx = fe.JxW(k) * factor;
63  for (unsigned int i=0; i<n_dofs; ++i)
64  {
65  double Mii = 0.0;
66  for (unsigned int d=0; d<n_components; ++d)
67  Mii += dx *
68  (fe.shape_grad_component(i,k,d) * fe.shape_grad_component(i,k,d));
69 
70  M(i,i) += Mii;
71 
72  for (unsigned int j=i+1; j<n_dofs; ++j)
73  {
74  double Mij = 0.0;
75  for (unsigned int d=0; d<n_components; ++d)
76  Mij += dx *
77  (fe.shape_grad_component(j,k,d) * fe.shape_grad_component(i,k,d));
78 
79  M(i,j) += Mij;
80  M(j,i) += Mij;
81  }
82  }
83  }
84 
85  }
86 
92  template <int dim>
93  inline void
95  Vector<double> &result,
96  const FEValuesBase<dim> &fe,
97  const std::vector<Tensor<1,dim> > &input,
98  double factor = 1.)
99  {
100  const unsigned int nq = fe.n_quadrature_points;
101  const unsigned int n_dofs = fe.dofs_per_cell;
102  Assert(input.size() == nq, ExcDimensionMismatch(input.size(), nq));
103  Assert(result.size() == n_dofs, ExcDimensionMismatch(result.size(), n_dofs));
104 
105  for (unsigned int k=0; k<nq; ++k)
106  {
107  const double dx = factor * fe.JxW(k);
108  for (unsigned int i=0; i<n_dofs; ++i)
109  result(i) += dx * (input[k] * fe.shape_grad(i,k));
110  }
111  }
112 
113 
119  template <int dim>
120  inline void
122  Vector<double> &result,
123  const FEValuesBase<dim> &fe,
124  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &input,
125  double factor = 1.)
126  {
127  const unsigned int nq = fe.n_quadrature_points;
128  const unsigned int n_dofs = fe.dofs_per_cell;
129  const unsigned int n_comp = fe.get_fe().n_components();
130 
132  Assert(result.size() == n_dofs, ExcDimensionMismatch(result.size(), n_dofs));
133 
134  for (unsigned int k=0; k<nq; ++k)
135  {
136  const double dx = factor * fe.JxW(k);
137  for (unsigned int i=0; i<n_dofs; ++i)
138  for (unsigned int d=0; d<n_comp; ++d)
139  {
140 
141  result(i) += dx * (input[d][k] * fe.shape_grad_component(i,k,d));
142  }
143  }
144  }
145 
146 
160  template <int dim>
163  const FEValuesBase<dim> &fe,
164  double penalty,
165  double factor = 1.)
166  {
167  const unsigned int n_dofs = fe.dofs_per_cell;
168  const unsigned int n_comp = fe.get_fe().n_components();
169 
170  Assert (M.m() == n_dofs, ExcDimensionMismatch(M.m(), n_dofs));
171  Assert (M.n() == n_dofs, ExcDimensionMismatch(M.n(), n_dofs));
172 
173  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
174  {
175  const double dx = fe.JxW(k) * factor;
176  const Tensor<1,dim> n = fe.normal_vector(k);
177  for (unsigned int i=0; i<n_dofs; ++i)
178  for (unsigned int j=0; j<n_dofs; ++j)
179  for (unsigned int d=0; d<n_comp; ++d)
180  M(i,j) += dx *
181  (2. * fe.shape_value_component(i,k,d) * penalty * fe.shape_value_component(j,k,d)
182  - (n * fe.shape_grad_component(i,k,d)) * fe.shape_value_component(j,k,d)
183  - (n * fe.shape_grad_component(j,k,d)) * fe.shape_value_component(i,k,d));
184  }
185  }
186 
201  template <int dim>
204  const FEValuesBase<dim> &fe,
205  double penalty,
206  double factor = 1.)
207  {
208  const unsigned int n_dofs = fe.dofs_per_cell;
209  AssertDimension(fe.get_fe().n_components(), dim);
210  Assert (M.m() == n_dofs, ExcDimensionMismatch(M.m(), n_dofs));
211  Assert (M.n() == n_dofs, ExcDimensionMismatch(M.n(), n_dofs));
212 
213  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
214  {
215  const double dx = fe.JxW(k) * factor;
216  const Tensor<1,dim> n = fe.normal_vector(k);
217  for (unsigned int i=0; i<n_dofs; ++i)
218  for (unsigned int j=0; j<n_dofs; ++j)
219  {
220  double udotn = 0.;
221  double vdotn = 0.;
222  double ngradun = 0.;
223  double ngradvn = 0.;
224 
225  for (unsigned int d=0; d<dim; ++d)
226  {
227  udotn += n[d]*fe.shape_value_component(j,k,d);
228  vdotn += n[d]*fe.shape_value_component(i,k,d);
229  ngradun += n*fe.shape_grad_component(j,k,d)*n[d];
230  ngradvn += n*fe.shape_grad_component(i,k,d)*n[d];
231  }
232 
233  for (unsigned int d=0; d<dim; ++d)
234  {
235  const double v_t = fe.shape_value_component(i,k,d)-vdotn*n[d];
236  const double dnv_t = n * fe.shape_grad_component(i,k,d)-ngradvn*n[d];
237  const double u_t = fe.shape_value_component(j,k,d)-udotn*n[d];
238  const double dnu_t = n * fe.shape_grad_component(j,k,d)-ngradun*n[d];
239 
240  M(i,j) += dx *(2.*penalty*u_t*v_t-dnu_t*v_t-dnv_t*u_t);
241  }
242  }
243  }
244  }
245 
261  template <int dim>
263  Vector<double> &result,
264  const FEValuesBase<dim> &fe,
265  const std::vector<double> &input,
266  const std::vector<Tensor<1,dim> > &Dinput,
267  const std::vector<double> &data,
268  double penalty,
269  double factor = 1.)
270  {
271  const unsigned int n_dofs = fe.dofs_per_cell;
272  AssertDimension(input.size(), fe.n_quadrature_points);
273  AssertDimension(Dinput.size(), fe.n_quadrature_points);
274  AssertDimension(data.size(), fe.n_quadrature_points);
275 
276  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
277  {
278  const double dx = factor * fe.JxW(k);
279  const Tensor<1,dim> n = fe.normal_vector(k);
280  for (unsigned int i=0; i<n_dofs; ++i)
281  {
282  const double dnv = fe.shape_grad(i,k) * n;
283  const double dnu = Dinput[k] * n;
284  const double v= fe.shape_value(i,k);
285  const double u= input[k];
286  const double g= data[k];
287 
288  result(i) += dx*(2.*penalty*(u-g)*v - dnv*(u-g) - dnu*v);
289  }
290  }
291  }
292 
310  template <int dim>
312  Vector<double> &result,
313  const FEValuesBase<dim> &fe,
314  const VectorSlice<const std::vector<std::vector<double> > > &input,
315  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput,
316  const VectorSlice<const std::vector<std::vector<double> > > &data,
317  double penalty,
318  double factor = 1.)
319  {
320  const unsigned int n_dofs = fe.dofs_per_cell;
321  const unsigned int n_comp = fe.get_fe().n_components();
325 
326  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
327  {
328  const double dx = factor * fe.JxW(k);
329  const Tensor<1,dim> n = fe.normal_vector(k);
330  for (unsigned int i=0; i<n_dofs; ++i)
331  for (unsigned int d=0; d<n_comp; ++d)
332  {
333  const double dnv = fe.shape_grad_component(i,k,d) * n;
334  const double dnu = Dinput[d][k] * n;
335  const double v= fe.shape_value_component(i,k,d);
336  const double u= input[d][k];
337  const double g= data[d][k];
338 
339  result(i) += dx*(2.*penalty*(u-g)*v - dnv*(u-g) - dnu*v);
340  }
341  }
342  }
343 
363  template <int dim>
364  void ip_matrix (
365  FullMatrix<double> &M11,
366  FullMatrix<double> &M12,
367  FullMatrix<double> &M21,
368  FullMatrix<double> &M22,
369  const FEValuesBase<dim> &fe1,
370  const FEValuesBase<dim> &fe2,
371  double penalty,
372  double factor1 = 1.,
373  double factor2 = -1.)
374  {
375  const unsigned int n_dofs = fe1.dofs_per_cell;
376  AssertDimension(M11.n(), n_dofs);
377  AssertDimension(M11.m(), n_dofs);
378  AssertDimension(M12.n(), n_dofs);
379  AssertDimension(M12.m(), n_dofs);
380  AssertDimension(M21.n(), n_dofs);
381  AssertDimension(M21.m(), n_dofs);
382  AssertDimension(M22.n(), n_dofs);
383  AssertDimension(M22.m(), n_dofs);
384 
385  const double nui = factor1;
386  const double nue = (factor2 < 0) ? factor1 : factor2;
387  const double nu = .5*(nui+nue);
388 
389  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
390  {
391  const double dx = fe1.JxW(k);
392  const Tensor<1,dim> n = fe1.normal_vector(k);
393  for (unsigned int d=0; d<fe1.get_fe().n_components(); ++d)
394  {
395  for (unsigned int i=0; i<n_dofs; ++i)
396  {
397  for (unsigned int j=0; j<n_dofs; ++j)
398  {
399  const double vi = fe1.shape_value_component(i,k,d);
400  const double dnvi = n * fe1.shape_grad_component(i,k,d);
401  const double ve = fe2.shape_value_component(i,k,d);
402  const double dnve = n * fe2.shape_grad_component(i,k,d);
403  const double ui = fe1.shape_value_component(j,k,d);
404  const double dnui = n * fe1.shape_grad_component(j,k,d);
405  const double ue = fe2.shape_value_component(j,k,d);
406  const double dnue = n * fe2.shape_grad_component(j,k,d);
407  M11(i,j) += dx*(-.5*nui*dnvi*ui-.5*nui*dnui*vi+nu*penalty*ui*vi);
408  M12(i,j) += dx*( .5*nui*dnvi*ue-.5*nue*dnue*vi-nu*penalty*vi*ue);
409  M21(i,j) += dx*(-.5*nue*dnve*ui+.5*nui*dnui*ve-nu*penalty*ui*ve);
410  M22(i,j) += dx*( .5*nue*dnve*ue+.5*nue*dnue*ve+nu*penalty*ue*ve);
411  }
412  }
413  }
414  }
415  }
416 
431  template <int dim>
433  FullMatrix<double> &M11,
434  FullMatrix<double> &M12,
435  FullMatrix<double> &M21,
436  FullMatrix<double> &M22,
437  const FEValuesBase<dim> &fe1,
438  const FEValuesBase<dim> &fe2,
439  double penalty,
440  double factor1 = 1.,
441  double factor2 = -1.)
442  {
443  const unsigned int n_dofs = fe1.dofs_per_cell;
444  AssertDimension(fe1.get_fe().n_components(), dim);
445  AssertDimension(fe2.get_fe().n_components(), dim);
446  AssertDimension(M11.n(), n_dofs);
447  AssertDimension(M11.m(), n_dofs);
448  AssertDimension(M12.n(), n_dofs);
449  AssertDimension(M12.m(), n_dofs);
450  AssertDimension(M21.n(), n_dofs);
451  AssertDimension(M21.m(), n_dofs);
452  AssertDimension(M22.n(), n_dofs);
453  AssertDimension(M22.m(), n_dofs);
454 
455  const double nui = factor1;
456  const double nue = (factor2 < 0) ? factor1 : factor2;
457  const double nu = .5*(nui+nue);
458 
459  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
460  {
461  const double dx = fe1.JxW(k);
462  const Tensor<1,dim> n = fe1.normal_vector(k);
463  for (unsigned int i=0; i<n_dofs; ++i)
464  {
465  // We compute the tangential component by subtracting
466  // the normal component from the total. Thus, compute
467  // the normal component first
468  for (unsigned int j=0; j<n_dofs; ++j)
469  {
470  double u1dotn = 0.;
471  double v1dotn = 0.;
472  double u2dotn = 0.;
473  double v2dotn = 0.;
474 
475  double ngradu1n = 0.;
476  double ngradv1n = 0.;
477  double ngradu2n = 0.;
478  double ngradv2n = 0.;
479 
480  for (unsigned int d=0; d<dim; ++d)
481  {
482  u1dotn += n[d]*fe1.shape_value_component(j,k,d);
483  v1dotn += n[d]*fe1.shape_value_component(i,k,d);
484  u2dotn += n[d]*fe2.shape_value_component(j,k,d);
485  v2dotn += n[d]*fe2.shape_value_component(i,k,d);
486 
487  ngradu1n += n*fe1.shape_grad_component(j,k,d)*n[d];
488  ngradv1n += n*fe1.shape_grad_component(i,k,d)*n[d];
489  ngradu2n += n*fe2.shape_grad_component(j,k,d)*n[d];
490  ngradv2n += n*fe2.shape_grad_component(i,k,d)*n[d];
491  }
492 
493  // The following code is equal to ip_matrix() with
494  // the only exception that all variables introduced
495  // below denote tangential components.
496  for (unsigned int d=0; d<dim; ++d)
497  {
498  const double vi = fe1.shape_value_component(i,k,d)-v1dotn*n[d];
499  const double dnvi = n * fe1.shape_grad_component(i,k,d)-ngradv1n*n[d];
500 
501  const double ve = fe2.shape_value_component(i,k,d)-v2dotn*n[d];
502  const double dnve = n * fe2.shape_grad_component(i,k,d)-ngradv2n*n[d];
503 
504  const double ui = fe1.shape_value_component(j,k,d)-u1dotn*n[d];
505  const double dnui = n * fe1.shape_grad_component(j,k,d)-ngradu1n*n[d];
506 
507  const double ue = fe2.shape_value_component(j,k,d)-u2dotn*n[d];
508  const double dnue = n * fe2.shape_grad_component(j,k,d)-ngradu2n*n[d];
509 
510  M11(i,j) += dx*(-.5*nui*dnvi*ui-.5*nui*dnui*vi+nu*penalty*ui*vi);
511  M12(i,j) += dx*( .5*nui*dnvi*ue-.5*nue*dnue*vi-nu*penalty*vi*ue);
512  M21(i,j) += dx*(-.5*nue*dnve*ui+.5*nui*dnui*ve-nu*penalty*ui*ve);
513  M22(i,j) += dx*( .5*nue*dnve*ue+.5*nue*dnue*ve+nu*penalty*ue*ve);
514  }
515  }
516  }
517  }
518  }
519 
530  template<int dim>
531  void
533  Vector<double> &result1,
534  Vector<double> &result2,
535  const FEValuesBase<dim> &fe1,
536  const FEValuesBase<dim> &fe2,
537  const std::vector<double> &input1,
538  const std::vector<Tensor<1,dim> > &Dinput1,
539  const std::vector<double> &input2,
540  const std::vector<Tensor<1,dim> > &Dinput2,
541  double pen,
542  double int_factor = 1.,
543  double ext_factor = -1.)
544  {
545  Assert(fe1.get_fe().n_components() == 1,
546  ExcDimensionMismatch(fe1.get_fe().n_components(), 1));
547  Assert(fe2.get_fe().n_components() == 1,
548  ExcDimensionMismatch(fe2.get_fe().n_components(), 1));
549 
550  const double nui = int_factor;
551  const double nue = (ext_factor < 0) ? int_factor : ext_factor;
552  const double penalty = .5 * pen * (nui + nue);
553 
554  const unsigned int n_dofs = fe1.dofs_per_cell;
555 
556  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
557  {
558  const double dx = fe1.JxW(k);
559  const Tensor<1,dim> n = fe1.normal_vector(k);
560 
561  for (unsigned int i=0; i<n_dofs; ++i)
562  {
563  const double vi = fe1.shape_value(i,k);
564  const Tensor<1,dim> &Dvi = fe1.shape_grad(i,k);
565  const double dnvi = Dvi * n;
566  const double ve = fe2.shape_value(i,k);
567  const Tensor<1,dim> &Dve = fe2.shape_grad(i,k);
568  const double dnve = Dve * n;
569 
570  const double ui = input1[k];
571  const Tensor<1,dim> &Dui = Dinput1[k];
572  const double dnui = Dui * n;
573  const double ue = input2[k];
574  const Tensor<1,dim> &Due = Dinput2[k];
575  const double dnue = Due * n;
576 
577  result1(i) += dx*(-.5*nui*dnvi*ui-.5*nui*dnui*vi+penalty*ui*vi);
578  result1(i) += dx*( .5*nui*dnvi*ue-.5*nue*dnue*vi-penalty*vi*ue);
579  result2(i) += dx*(-.5*nue*dnve*ui+.5*nui*dnui*ve-penalty*ui*ve);
580  result2(i) += dx*( .5*nue*dnve*ue+.5*nue*dnue*ve+penalty*ue*ve);
581  }
582  }
583  }
584 
585 
597  template<int dim>
598  void
600  Vector<double> &result1,
601  Vector<double> &result2,
602  const FEValuesBase<dim> &fe1,
603  const FEValuesBase<dim> &fe2,
604  const VectorSlice<const std::vector<std::vector<double> > > &input1,
605  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput1,
606  const VectorSlice<const std::vector<std::vector<double> > > &input2,
607  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput2,
608  double pen,
609  double int_factor = 1.,
610  double ext_factor = -1.)
611  {
612  const unsigned int n_comp = fe1.get_fe().n_components();
613  const unsigned int n1 = fe1.dofs_per_cell;
614 
616  AssertVectorVectorDimension(Dinput1, n_comp, fe1.n_quadrature_points);
618  AssertVectorVectorDimension(Dinput2, n_comp, fe2.n_quadrature_points);
619 
620  const double nui = int_factor;
621  const double nue = (ext_factor < 0) ? int_factor : ext_factor;
622  const double penalty = .5 * pen * (nui + nue);
623 
624 
625  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
626  {
627  const double dx = fe1.JxW(k);
628  const Tensor<1,dim> n = fe1.normal_vector(k);
629 
630  for (unsigned int i=0; i<n1; ++i)
631  for (unsigned int d=0; d<n_comp; ++d)
632  {
633  const double vi = fe1.shape_value_component(i,k,d);
634  const Tensor<1,dim> &Dvi = fe1.shape_grad_component(i,k,d);
635  const double dnvi = Dvi * n;
636  const double ve = fe2.shape_value_component(i,k,d);
637  const Tensor<1,dim> &Dve = fe2.shape_grad_component(i,k,d);
638  const double dnve = Dve * n;
639 
640  const double ui = input1[d][k];
641  const Tensor<1,dim> &Dui = Dinput1[d][k];
642  const double dnui = Dui * n;
643  const double ue = input2[d][k];
644  const Tensor<1,dim> &Due = Dinput2[d][k];
645  const double dnue = Due * n;
646 
647  result1(i) += dx*(-.5*nui*dnvi*ui-.5*nui*dnui*vi+penalty*ui*vi);
648  result1(i) += dx*( .5*nui*dnvi*ue-.5*nue*dnue*vi-penalty*vi*ue);
649  result2(i) += dx*(-.5*nue*dnve*ui+.5*nui*dnui*ve-penalty*ui*ve);
650  result2(i) += dx*( .5*nue*dnve*ue+.5*nue*dnue*ve+penalty*ue*ve);
651  }
652  }
653  }
654 
655 
656 
671  template <int dim, int spacedim, typename number>
675  unsigned int deg1,
676  unsigned int deg2)
677  {
678  const unsigned int normal1 = GeometryInfo<dim>::unit_normal_direction[dinfo1.face_number];
679  const unsigned int normal2 = GeometryInfo<dim>::unit_normal_direction[dinfo2.face_number];
680  const unsigned int deg1sq = (deg1 == 0) ? 1 : deg1 * (deg1+1);
681  const unsigned int deg2sq = (deg2 == 0) ? 1 : deg2 * (deg2+1);
682 
683  double penalty1 = deg1sq / dinfo1.cell->extent_in_direction(normal1);
684  double penalty2 = deg2sq / dinfo2.cell->extent_in_direction(normal2);
685  if (dinfo1.cell->has_children() ^ dinfo2.cell->has_children())
686  {
687  Assert (dinfo1.face == dinfo2.face, ExcInternalError());
688  Assert (dinfo1.face->has_children(), ExcInternalError());
689  penalty1 *= 2;
690  }
691  const double penalty = 0.5*(penalty1 + penalty2);
692  return penalty;
693  }
694  }
695 }
696 
697 
698 DEAL_II_NAMESPACE_CLOSE
699 
700 #endif
size_type m() const
Triangulation< dim, spacedim >::cell_iterator cell
The current cell.
Definition: dof_info.h:72
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
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:364
void nitsche_tangential_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
Definition: laplace.h:202
const unsigned int dofs_per_cell
Definition: fe_values.h:1459
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
Definition: laplace.h:52
#define AssertVectorVectorDimension(vec, dim1, dim2)
Definition: exceptions.h:1157
const FiniteElement< dim, spacedim > & get_fe() const
void cell_residual(Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< Tensor< 1, dim > > &input, double factor=1.)
Definition: laplace.h:94
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 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:262
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:532
size_type n() const
const double & shape_value(const unsigned int function_no, const unsigned int point_no) 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
const unsigned int n_quadrature_points
Definition: fe_values.h:1452
void nitsche_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
Definition: laplace.h:161
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:672
double JxW(const unsigned int quadrature_point) const
unsigned int face_number
Definition: dof_info.h:83
Triangulation< dim, spacedim >::face_iterator face
The current face.
Definition: dof_info.h:75
const Tensor< 1, spacedim > & shape_grad(const unsigned int function_no, 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:432
static ::ExceptionBase & ExcInternalError()
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