Reference documentation for deal.II version 8.4.2
laplace.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2015 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 
202  template <int dim>
204  Vector<double> &result,
205  const FEValuesBase<dim> &fe,
206  const std::vector<double> &input,
207  const std::vector<Tensor<1,dim> > &Dinput,
208  const std::vector<double> &data,
209  double penalty,
210  double factor = 1.)
211  {
212  const unsigned int n_dofs = fe.dofs_per_cell;
213  AssertDimension(input.size(), fe.n_quadrature_points);
214  AssertDimension(Dinput.size(), fe.n_quadrature_points);
215  AssertDimension(data.size(), fe.n_quadrature_points);
216 
217  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
218  {
219  const double dx = factor * fe.JxW(k);
220  const Tensor<1,dim> n = fe.normal_vector(k);
221  for (unsigned int i=0; i<n_dofs; ++i)
222  {
223  const double dnv = fe.shape_grad(i,k) * n;
224  const double dnu = Dinput[k] * n;
225  const double v= fe.shape_value(i,k);
226  const double u= input[k];
227  const double g= data[k];
228 
229  result(i) += dx*(2.*penalty*(u-g)*v - dnv*(u-g) - dnu*v);
230  }
231  }
232  }
233 
251  template <int dim>
253  Vector<double> &result,
254  const FEValuesBase<dim> &fe,
255  const VectorSlice<const std::vector<std::vector<double> > > &input,
256  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput,
257  const VectorSlice<const std::vector<std::vector<double> > > &data,
258  double penalty,
259  double factor = 1.)
260  {
261  const unsigned int n_dofs = fe.dofs_per_cell;
262  const unsigned int n_comp = fe.get_fe().n_components();
266 
267  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
268  {
269  const double dx = factor * fe.JxW(k);
270  const Tensor<1,dim> n = fe.normal_vector(k);
271  for (unsigned int i=0; i<n_dofs; ++i)
272  for (unsigned int d=0; d<n_comp; ++d)
273  {
274  const double dnv = fe.shape_grad_component(i,k,d) * n;
275  const double dnu = Dinput[d][k] * n;
276  const double v= fe.shape_value_component(i,k,d);
277  const double u= input[d][k];
278  const double g= data[d][k];
279 
280  result(i) += dx*(2.*penalty*(u-g)*v - dnv*(u-g) - dnu*v);
281  }
282  }
283  }
284 
304  template <int dim>
305  void ip_matrix (
306  FullMatrix<double> &M11,
307  FullMatrix<double> &M12,
308  FullMatrix<double> &M21,
309  FullMatrix<double> &M22,
310  const FEValuesBase<dim> &fe1,
311  const FEValuesBase<dim> &fe2,
312  double penalty,
313  double factor1 = 1.,
314  double factor2 = -1.)
315  {
316  const unsigned int n_dofs = fe1.dofs_per_cell;
317  AssertDimension(M11.n(), n_dofs);
318  AssertDimension(M11.m(), n_dofs);
319  AssertDimension(M12.n(), n_dofs);
320  AssertDimension(M12.m(), n_dofs);
321  AssertDimension(M21.n(), n_dofs);
322  AssertDimension(M21.m(), n_dofs);
323  AssertDimension(M22.n(), n_dofs);
324  AssertDimension(M22.m(), n_dofs);
325 
326  const double nui = factor1;
327  const double nue = (factor2 < 0) ? factor1 : factor2;
328  const double nu = .5*(nui+nue);
329 
330  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
331  {
332  const double dx = fe1.JxW(k);
333  const Tensor<1,dim> n = fe1.normal_vector(k);
334  for (unsigned int d=0; d<fe1.get_fe().n_components(); ++d)
335  {
336  for (unsigned int i=0; i<n_dofs; ++i)
337  {
338  for (unsigned int j=0; j<n_dofs; ++j)
339  {
340  const double vi = fe1.shape_value_component(i,k,d);
341  const double dnvi = n * fe1.shape_grad_component(i,k,d);
342  const double ve = fe2.shape_value_component(i,k,d);
343  const double dnve = n * fe2.shape_grad_component(i,k,d);
344  const double ui = fe1.shape_value_component(j,k,d);
345  const double dnui = n * fe1.shape_grad_component(j,k,d);
346  const double ue = fe2.shape_value_component(j,k,d);
347  const double dnue = n * fe2.shape_grad_component(j,k,d);
348  M11(i,j) += dx*(-.5*nui*dnvi*ui-.5*nui*dnui*vi+nu*penalty*ui*vi);
349  M12(i,j) += dx*( .5*nui*dnvi*ue-.5*nue*dnue*vi-nu*penalty*vi*ue);
350  M21(i,j) += dx*(-.5*nue*dnve*ui+.5*nui*dnui*ve-nu*penalty*ui*ve);
351  M22(i,j) += dx*( .5*nue*dnve*ue+.5*nue*dnue*ve+nu*penalty*ue*ve);
352  }
353  }
354  }
355  }
356  }
357 
372  template <int dim>
374  FullMatrix<double> &M11,
375  FullMatrix<double> &M12,
376  FullMatrix<double> &M21,
377  FullMatrix<double> &M22,
378  const FEValuesBase<dim> &fe1,
379  const FEValuesBase<dim> &fe2,
380  double penalty,
381  double factor1 = 1.,
382  double factor2 = -1.)
383  {
384  const unsigned int n_dofs = fe1.dofs_per_cell;
385  AssertDimension(fe1.get_fe().n_components(), dim);
386  AssertDimension(fe2.get_fe().n_components(), dim);
387  AssertDimension(M11.n(), n_dofs);
388  AssertDimension(M11.m(), n_dofs);
389  AssertDimension(M12.n(), n_dofs);
390  AssertDimension(M12.m(), n_dofs);
391  AssertDimension(M21.n(), n_dofs);
392  AssertDimension(M21.m(), n_dofs);
393  AssertDimension(M22.n(), n_dofs);
394  AssertDimension(M22.m(), n_dofs);
395 
396  const double nui = factor1;
397  const double nue = (factor2 < 0) ? factor1 : factor2;
398  const double nu = .5*(nui+nue);
399 
400  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
401  {
402  const double dx = fe1.JxW(k);
403  const Tensor<1,dim> n = fe1.normal_vector(k);
404  for (unsigned int i=0; i<n_dofs; ++i)
405  {
406  for (unsigned int j=0; j<n_dofs; ++j)
407  {
408  double u1dotn = 0.;
409  double v1dotn = 0.;
410  double u2dotn = 0.;
411  double v2dotn = 0.;
412 
413  double ngradu1n = 0.;
414  double ngradv1n = 0.;
415  double ngradu2n = 0.;
416  double ngradv2n = 0.;
417 
418  for (unsigned int d=0; d<dim; ++d)
419  {
420  u1dotn += n[d]*fe1.shape_value_component(j,k,d);
421  v1dotn += n[d]*fe1.shape_value_component(i,k,d);
422  u2dotn += n[d]*fe2.shape_value_component(j,k,d);
423  v2dotn += n[d]*fe2.shape_value_component(i,k,d);
424 
425  ngradu1n += n*fe1.shape_grad_component(j,k,d)*n[d];
426  ngradv1n += n*fe1.shape_grad_component(i,k,d)*n[d];
427  ngradu2n += n*fe2.shape_grad_component(j,k,d)*n[d];
428  ngradv2n += n*fe2.shape_grad_component(i,k,d)*n[d];
429  }
430 
431  for (unsigned int d=0; d<fe1.get_fe().n_components(); ++d)
432  {
433  const double vi = fe1.shape_value_component(i,k,d)-v1dotn*n[d];
434  const double dnvi = n * fe1.shape_grad_component(i,k,d)-ngradv1n*n[d];
435 
436  const double ve = fe2.shape_value_component(i,k,d)-v2dotn*n[d];
437  const double dnve = n * fe2.shape_grad_component(i,k,d)-ngradv2n*n[d];
438 
439  const double ui = fe1.shape_value_component(j,k,d)-u1dotn*n[d];
440  const double dnui = n * fe1.shape_grad_component(j,k,d)-ngradu1n*n[d];
441 
442  const double ue = fe2.shape_value_component(j,k,d)-u2dotn*n[d];
443  const double dnue = n * fe2.shape_grad_component(j,k,d)-ngradu2n*n[d];
444 
445  M11(i,j) += dx*(-.5*nui*dnvi*ui-.5*nui*dnui*vi+nu*penalty*ui*vi);
446  M12(i,j) += dx*( .5*nui*dnvi*ue-.5*nue*dnue*vi-nu*penalty*vi*ue);
447  M21(i,j) += dx*(-.5*nue*dnve*ui+.5*nui*dnui*ve-nu*penalty*ui*ve);
448  M22(i,j) += dx*( .5*nue*dnve*ue+.5*nue*dnue*ve+nu*penalty*ue*ve);
449  }
450  }
451  }
452  }
453  }
454 
465  template<int dim>
466  void
468  Vector<double> &result1,
469  Vector<double> &result2,
470  const FEValuesBase<dim> &fe1,
471  const FEValuesBase<dim> &fe2,
472  const std::vector<double> &input1,
473  const std::vector<Tensor<1,dim> > &Dinput1,
474  const std::vector<double> &input2,
475  const std::vector<Tensor<1,dim> > &Dinput2,
476  double pen,
477  double int_factor = 1.,
478  double ext_factor = -1.)
479  {
480  Assert(fe1.get_fe().n_components() == 1,
481  ExcDimensionMismatch(fe1.get_fe().n_components(), 1));
482  Assert(fe2.get_fe().n_components() == 1,
483  ExcDimensionMismatch(fe2.get_fe().n_components(), 1));
484 
485  const double nui = int_factor;
486  const double nue = (ext_factor < 0) ? int_factor : ext_factor;
487  const double penalty = .5 * pen * (nui + nue);
488 
489  const unsigned int n_dofs = fe1.dofs_per_cell;
490 
491  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
492  {
493  const double dx = fe1.JxW(k);
494  const Tensor<1,dim> n = fe1.normal_vector(k);
495 
496  for (unsigned int i=0; i<n_dofs; ++i)
497  {
498  const double vi = fe1.shape_value(i,k);
499  const Tensor<1,dim> &Dvi = fe1.shape_grad(i,k);
500  const double dnvi = Dvi * n;
501  const double ve = fe2.shape_value(i,k);
502  const Tensor<1,dim> &Dve = fe2.shape_grad(i,k);
503  const double dnve = Dve * n;
504 
505  const double ui = input1[k];
506  const Tensor<1,dim> &Dui = Dinput1[k];
507  const double dnui = Dui * n;
508  const double ue = input2[k];
509  const Tensor<1,dim> &Due = Dinput2[k];
510  const double dnue = Due * n;
511 
512  result1(i) += dx*(-.5*nui*dnvi*ui-.5*nui*dnui*vi+penalty*ui*vi);
513  result1(i) += dx*( .5*nui*dnvi*ue-.5*nue*dnue*vi-penalty*vi*ue);
514  result2(i) += dx*(-.5*nue*dnve*ui+.5*nui*dnui*ve-penalty*ui*ve);
515  result2(i) += dx*( .5*nue*dnve*ue+.5*nue*dnue*ve+penalty*ue*ve);
516  }
517  }
518  }
519 
520 
532  template<int dim>
533  void
535  Vector<double> &result1,
536  Vector<double> &result2,
537  const FEValuesBase<dim> &fe1,
538  const FEValuesBase<dim> &fe2,
539  const VectorSlice<const std::vector<std::vector<double> > > &input1,
540  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput1,
541  const VectorSlice<const std::vector<std::vector<double> > > &input2,
542  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput2,
543  double pen,
544  double int_factor = 1.,
545  double ext_factor = -1.)
546  {
547  const unsigned int n_comp = fe1.get_fe().n_components();
548  const unsigned int n1 = fe1.dofs_per_cell;
549 
551  AssertVectorVectorDimension(Dinput1, n_comp, fe1.n_quadrature_points);
553  AssertVectorVectorDimension(Dinput2, n_comp, fe2.n_quadrature_points);
554 
555  const double nui = int_factor;
556  const double nue = (ext_factor < 0) ? int_factor : ext_factor;
557  const double penalty = .5 * pen * (nui + nue);
558 
559 
560  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
561  {
562  const double dx = fe1.JxW(k);
563  const Tensor<1,dim> n = fe1.normal_vector(k);
564 
565  for (unsigned int i=0; i<n1; ++i)
566  for (unsigned int d=0; d<n_comp; ++d)
567  {
568  const double vi = fe1.shape_value_component(i,k,d);
569  const Tensor<1,dim> &Dvi = fe1.shape_grad_component(i,k,d);
570  const double dnvi = Dvi * n;
571  const double ve = fe2.shape_value_component(i,k,d);
572  const Tensor<1,dim> &Dve = fe2.shape_grad_component(i,k,d);
573  const double dnve = Dve * n;
574 
575  const double ui = input1[d][k];
576  const Tensor<1,dim> &Dui = Dinput1[d][k];
577  const double dnui = Dui * n;
578  const double ue = input2[d][k];
579  const Tensor<1,dim> &Due = Dinput2[d][k];
580  const double dnue = Due * n;
581 
582  result1(i) += dx*(-.5*nui*dnvi*ui-.5*nui*dnui*vi+penalty*ui*vi);
583  result1(i) += dx*( .5*nui*dnvi*ue-.5*nue*dnue*vi-penalty*vi*ue);
584  result2(i) += dx*(-.5*nue*dnve*ui+.5*nui*dnui*ve-penalty*ui*ve);
585  result2(i) += dx*( .5*nue*dnve*ue+.5*nue*dnue*ve+penalty*ue*ve);
586  }
587  }
588  }
589 
590 
591 
606  template <int dim, int spacedim, typename number>
610  unsigned int deg1,
611  unsigned int deg2)
612  {
613  const unsigned int normal1 = GeometryInfo<dim>::unit_normal_direction[dinfo1.face_number];
614  const unsigned int normal2 = GeometryInfo<dim>::unit_normal_direction[dinfo2.face_number];
615  const unsigned int deg1sq = (deg1 == 0) ? 1 : deg1 * (deg1+1);
616  const unsigned int deg2sq = (deg2 == 0) ? 1 : deg2 * (deg2+1);
617 
618  double penalty1 = deg1sq / dinfo1.cell->extent_in_direction(normal1);
619  double penalty2 = deg2sq / dinfo2.cell->extent_in_direction(normal2);
620  if (dinfo1.cell->has_children() ^ dinfo2.cell->has_children())
621  {
622  Assert (dinfo1.face == dinfo2.face, ExcInternalError());
623  Assert (dinfo1.face->has_children(), ExcInternalError());
624  penalty1 *= 2;
625  }
626  const double penalty = 0.5*(penalty1 + penalty2);
627  return penalty;
628  }
629  }
630 }
631 
632 
633 DEAL_II_NAMESPACE_CLOSE
634 
635 #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:1052
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:305
const unsigned int dofs_per_cell
Definition: fe_values.h:1464
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:1063
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:203
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:467
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:294
std::size_t size() const
const unsigned int n_quadrature_points
Definition: fe_values.h:1457
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:607
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:373
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