Reference documentation for deal.II version 9.1.0-pre
grad_div.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_grad_div_h
17 #define dealii_integrators_grad_div_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 {
44  namespace GradDiv
45  {
55  template <int dim>
56  void
58  const FEValuesBase<dim> &fe,
59  double factor = 1.)
60  {
61  const unsigned int n_dofs = fe.dofs_per_cell;
62 
63  AssertDimension(fe.get_fe().n_components(), dim);
64  AssertDimension(M.m(), n_dofs);
65  AssertDimension(M.n(), n_dofs);
66 
67  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
68  {
69  const double dx = factor * fe.JxW(k);
70  for (unsigned int i = 0; i < n_dofs; ++i)
71  for (unsigned int j = 0; j < n_dofs; ++j)
72  {
73  const double divu =
74  fe[FEValuesExtractors::Vector(0)].divergence(j, k);
75  const double divv =
76  fe[FEValuesExtractors::Vector(0)].divergence(i, k);
77 
78  M(i, j) += dx * divu * divv;
79  }
80  }
81  }
82 
92  template <int dim, typename number>
93  void
95  Vector<number> & result,
96  const FEValuesBase<dim> & fetest,
97  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &input,
98  const double factor = 1.)
99  {
100  const unsigned int n_dofs = fetest.dofs_per_cell;
101 
102  AssertDimension(fetest.get_fe().n_components(), dim);
104 
105  for (unsigned int k = 0; k < fetest.n_quadrature_points; ++k)
106  {
107  const double dx = factor * fetest.JxW(k);
108  for (unsigned int i = 0; i < n_dofs; ++i)
109  {
110  const double divv =
111  fetest[FEValuesExtractors::Vector(0)].divergence(i, k);
112  double du = 0.;
113  for (unsigned int d = 0; d < dim; ++d)
114  du += input[d][k][d];
115 
116  result(i) += dx * du * divv;
117  }
118  }
119  }
120 
129  template <int dim>
130  inline void
132  const FEValuesBase<dim> &fe,
133  double penalty,
134  double factor = 1.)
135  {
136  const unsigned int n_dofs = fe.dofs_per_cell;
137 
138  AssertDimension(fe.get_fe().n_components(), dim);
139  AssertDimension(M.m(), n_dofs);
140  AssertDimension(M.n(), n_dofs);
141 
142  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
143  {
144  const double dx = factor * fe.JxW(k);
145  const Tensor<1, dim> n = fe.normal_vector(k);
146  for (unsigned int i = 0; i < n_dofs; ++i)
147  for (unsigned int j = 0; j < n_dofs; ++j)
148  {
149  const double divu =
150  fe[FEValuesExtractors::Vector(0)].divergence(j, k);
151  const double divv =
152  fe[FEValuesExtractors::Vector(0)].divergence(i, k);
153  double un = 0., vn = 0.;
154  for (unsigned int d = 0; d < dim; ++d)
155  {
156  un += fe.shape_value_component(j, k, d) * n[d];
157  vn += fe.shape_value_component(i, k, d) * n[d];
158  }
159 
160  M(i, j) += dx * 2. * penalty * un * vn;
161  M(i, j) -= dx * (divu * vn + divv * un);
162  }
163  }
164  }
165 
184  template <int dim>
185  void
187  Vector<double> & result,
188  const FEValuesBase<dim> & fe,
189  const VectorSlice<const std::vector<std::vector<double>>> & input,
190  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &Dinput,
191  const VectorSlice<const std::vector<std::vector<double>>> & data,
192  double penalty,
193  double factor = 1.)
194  {
195  const unsigned int n_dofs = fe.dofs_per_cell;
196  AssertDimension(fe.get_fe().n_components(), dim)
200 
201  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
202  {
203  const double dx = factor * fe.JxW(k);
204  const Tensor<1, dim> n = fe.normal_vector(k);
205 
206  double umgn = 0.;
207  double divu = 0.;
208  for (unsigned int d = 0; d < dim; ++d)
209  {
210  umgn += (input[d][k] - data[d][k]) * n[d];
211  divu += Dinput[d][k][d];
212  }
213 
214  for (unsigned int i = 0; i < n_dofs; ++i)
215  {
216  double vn = 0.;
217  const double divv =
218  fe[FEValuesExtractors::Vector(0)].divergence(i, k);
219  for (unsigned int d = 0; d < dim; ++d)
220  vn += fe.shape_value_component(i, k, d) * n[d];
221 
222  result(i) +=
223  dx * (2. * penalty * umgn * vn - divv * umgn - divu * vn);
224  }
225  }
226  }
227 
236  template <int dim>
237  void
239  FullMatrix<double> & M12,
240  FullMatrix<double> & M21,
241  FullMatrix<double> & M22,
242  const FEValuesBase<dim> &fe1,
243  const FEValuesBase<dim> &fe2,
244  double penalty,
245  double factor1 = 1.,
246  double factor2 = -1.)
247  {
248  const unsigned int n_dofs = fe1.dofs_per_cell;
249  AssertDimension(M11.n(), n_dofs);
250  AssertDimension(M11.m(), n_dofs);
251  AssertDimension(M12.n(), n_dofs);
252  AssertDimension(M12.m(), n_dofs);
253  AssertDimension(M21.n(), n_dofs);
254  AssertDimension(M21.m(), n_dofs);
255  AssertDimension(M22.n(), n_dofs);
256  AssertDimension(M22.m(), n_dofs);
257 
258  const double fi = factor1;
259  const double fe = (factor2 < 0) ? factor1 : factor2;
260  const double f = .5 * (fi + fe);
261 
262  for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
263  {
264  const double dx = fe1.JxW(k);
265  const Tensor<1, dim> n = fe1.normal_vector(k);
266  for (unsigned int i = 0; i < n_dofs; ++i)
267  for (unsigned int j = 0; j < n_dofs; ++j)
268  {
269  double uni = 0.;
270  double une = 0.;
271  double vni = 0.;
272  double vne = 0.;
273  const double divui =
274  fe1[FEValuesExtractors::Vector(0)].divergence(j, k);
275  const double divue =
276  fe2[FEValuesExtractors::Vector(0)].divergence(j, k);
277  const double divvi =
278  fe1[FEValuesExtractors::Vector(0)].divergence(i, k);
279  const double divve =
280  fe2[FEValuesExtractors::Vector(0)].divergence(i, k);
281 
282  for (unsigned int d = 0; d < dim; ++d)
283  {
284  uni += fe1.shape_value_component(j, k, d) * n[d];
285  une += fe2.shape_value_component(j, k, d) * n[d];
286  vni += fe1.shape_value_component(i, k, d) * n[d];
287  vne += fe2.shape_value_component(i, k, d) * n[d];
288  }
289  M11(i, j) +=
290  dx * (-.5 * fi * divvi * uni - .5 * fi * divui * vni +
291  f * penalty * uni * vni);
292  M12(i, j) +=
293  dx * (.5 * fi * divvi * une - .5 * fe * divue * vni -
294  f * penalty * vni * une);
295  M21(i, j) +=
296  dx * (-.5 * fe * divve * uni + .5 * fi * divui * vne -
297  f * penalty * uni * vne);
298  M22(i, j) +=
299  dx * (.5 * fe * divve * une + .5 * fe * divue * vne +
300  f * penalty * une * vne);
301  }
302  }
303  }
304 
319  template <int dim>
320  void
322  Vector<double> & result1,
323  Vector<double> & result2,
324  const FEValuesBase<dim> & fe1,
325  const FEValuesBase<dim> & fe2,
326  const VectorSlice<const std::vector<std::vector<double>>> &input1,
327  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>>
328  & Dinput1,
329  const VectorSlice<const std::vector<std::vector<double>>> &input2,
330  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>>
331  & Dinput2,
332  double pen,
333  double int_factor = 1.,
334  double ext_factor = -1.)
335  {
336  const unsigned int n1 = fe1.dofs_per_cell;
337 
338  AssertDimension(fe1.get_fe().n_components(), dim);
343 
344  const double fi = int_factor;
345  const double fe = (ext_factor < 0) ? int_factor : ext_factor;
346  const double penalty = .5 * pen * (fi + fe);
347 
348 
349  for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
350  {
351  const double dx = fe1.JxW(k);
352  const Tensor<1, dim> n = fe1.normal_vector(k);
353  double uni = 0.;
354  double une = 0.;
355  double divui = 0.;
356  double divue = 0.;
357  for (unsigned int d = 0; d < dim; ++d)
358  {
359  uni += input1[d][k] * n[d];
360  une += input2[d][k] * n[d];
361  divui += Dinput1[d][k][d];
362  divue += Dinput2[d][k][d];
363  }
364 
365  for (unsigned int i = 0; i < n1; ++i)
366  {
367  double vni = 0.;
368  double vne = 0.;
369  const double divvi =
370  fe1[FEValuesExtractors::Vector(0)].divergence(i, k);
371  const double divve =
372  fe2[FEValuesExtractors::Vector(0)].divergence(i, k);
373  for (unsigned int d = 0; d < dim; ++d)
374  {
375  vni += fe1.shape_value_component(i, k, d) * n[d];
376  vne += fe2.shape_value_component(i, k, d) * n[d];
377  }
378 
379  result1(i) += dx * (-.5 * fi * divvi * uni -
380  .5 * fi * divui * vni + penalty * uni * vni);
381  result1(i) += dx * (.5 * fi * divvi * une -
382  .5 * fe * divue * vni - penalty * vni * une);
383  result2(i) += dx * (-.5 * fe * divve * uni +
384  .5 * fi * divui * vne - penalty * uni * vne);
385  result2(i) += dx * (.5 * fe * divve * une +
386  .5 * fe * divue * vne + penalty * une * vne);
387  }
388  }
389  }
390  } // namespace GradDiv
391 } // namespace LocalIntegrators
392 
393 DEAL_II_NAMESPACE_CLOSE
394 
395 
396 #endif
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
void nitsche_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
Definition: grad_div.h:131
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1366
const unsigned int dofs_per_cell
Definition: fe_values.h:2033
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double factor=1.)
Definition: grad_div.h:57
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: grad_div.h:238
#define AssertVectorVectorDimension(vec, dim1, dim2)
Definition: exceptions.h:1377
Library of integrals over cells and faces.
Definition: advection.h:34
void cell_residual(Vector< number > &result, const FEValuesBase< dim > &fetest, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim >>>> &input, const double factor=1.)
Definition: grad_div.h:94
size_type n() const
void nitsche_residual(Vector< double > &result, const FEValuesBase< dim > &fe, const VectorSlice< const std::vector< std::vector< double >>> &input, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim >>>> &Dinput, const VectorSlice< const std::vector< std::vector< double >>> &data, double penalty, double factor=1.)
Definition: grad_div.h:186
const FiniteElement< dim, spacedim > & get_fe() const
size_type m() const
const unsigned int n_quadrature_points
Definition: fe_values.h:2026
void ip_residual(Vector< double > &result1, Vector< double > &result2, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const VectorSlice< const std::vector< std::vector< double >>> &input1, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim >>>> &Dinput1, const VectorSlice< const std::vector< std::vector< double >>> &input2, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim >>>> &Dinput2, double pen, double int_factor=1., double ext_factor=-1.)
Definition: grad_div.h:321
double JxW(const unsigned int quadrature_point) const