Reference documentation for deal.II version 9.1.0-pre
divergence.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_divergence_h
17 #define dealii_integrators_divergence_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/integrators/grad_div.h>
29 
30 #include <deal.II/lac/full_matrix.h>
31 
32 #include <deal.II/meshworker/dof_info.h>
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 namespace LocalIntegrators
37 {
46  namespace Divergence
47  {
57  template <int dim>
58  void
60  const FEValuesBase<dim> &fe,
61  const FEValuesBase<dim> &fetest,
62  double factor = 1.)
63  {
64  const unsigned int n_dofs = fe.dofs_per_cell;
65  const unsigned int t_dofs = fetest.dofs_per_cell;
66  AssertDimension(fe.get_fe().n_components(), dim);
67  AssertDimension(fetest.get_fe().n_components(), 1);
68  AssertDimension(M.m(), t_dofs);
69  AssertDimension(M.n(), n_dofs);
70 
71  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
72  {
73  const double dx = fe.JxW(k) * factor;
74  for (unsigned int i = 0; i < t_dofs; ++i)
75  {
76  const double vv = fetest.shape_value(i, k);
77  for (unsigned int d = 0; d < dim; ++d)
78  for (unsigned int j = 0; j < n_dofs; ++j)
79  {
80  const double du = fe.shape_grad_component(j, k, d)[d];
81  M(i, j) += dx * du * vv;
82  }
83  }
84  }
85  }
86 
99  template <int dim, typename number>
100  void
102  Vector<number> & result,
103  const FEValuesBase<dim> & fetest,
104  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &input,
105  const double factor = 1.)
106  {
107  AssertDimension(fetest.get_fe().n_components(), 1);
109  const unsigned int t_dofs = fetest.dofs_per_cell;
110  Assert(result.size() == t_dofs,
111  ExcDimensionMismatch(result.size(), t_dofs));
112 
113  for (unsigned int k = 0; k < fetest.n_quadrature_points; ++k)
114  {
115  const double dx = factor * fetest.JxW(k);
116 
117  for (unsigned int i = 0; i < t_dofs; ++i)
118  for (unsigned int d = 0; d < dim; ++d)
119  result(i) += dx * input[d][k][d] * fetest.shape_value(i, k);
120  }
121  }
122 
123 
136  template <int dim, typename number>
137  void
139  Vector<number> & result,
140  const FEValuesBase<dim> & fetest,
141  const VectorSlice<const std::vector<std::vector<double>>> &input,
142  const double factor = 1.)
143  {
144  AssertDimension(fetest.get_fe().n_components(), 1);
146  const unsigned int t_dofs = fetest.dofs_per_cell;
147  Assert(result.size() == t_dofs,
148  ExcDimensionMismatch(result.size(), t_dofs));
149 
150  for (unsigned int k = 0; k < fetest.n_quadrature_points; ++k)
151  {
152  const double dx = factor * fetest.JxW(k);
153 
154  for (unsigned int i = 0; i < t_dofs; ++i)
155  for (unsigned int d = 0; d < dim; ++d)
156  result(i) -= dx * input[d][k] * fetest.shape_grad(i, k)[d];
157  }
158  }
159 
160 
171  template <int dim>
172  void
174  const FEValuesBase<dim> &fe,
175  const FEValuesBase<dim> &fetest,
176  double factor = 1.)
177  {
178  const unsigned int t_dofs = fetest.dofs_per_cell;
179  const unsigned int n_dofs = fe.dofs_per_cell;
180 
181  AssertDimension(fetest.get_fe().n_components(), dim);
182  AssertDimension(fe.get_fe().n_components(), 1);
183  AssertDimension(M.m(), t_dofs);
184  AssertDimension(M.n(), n_dofs);
185 
186  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
187  {
188  const double dx = fe.JxW(k) * factor;
189  for (unsigned int d = 0; d < dim; ++d)
190  for (unsigned int i = 0; i < t_dofs; ++i)
191  {
192  const double vv = fetest.shape_value_component(i, k, d);
193  for (unsigned int j = 0; j < n_dofs; ++j)
194  {
195  const Tensor<1, dim> &Du = fe.shape_grad(j, k);
196  M(i, j) += dx * vv * Du[d];
197  }
198  }
199  }
200  }
201 
214  template <int dim, typename number>
215  void
217  const FEValuesBase<dim> & fetest,
218  const std::vector<Tensor<1, dim>> &input,
219  const double factor = 1.)
220  {
221  AssertDimension(fetest.get_fe().n_components(), dim);
222  AssertDimension(input.size(), fetest.n_quadrature_points);
223  const unsigned int t_dofs = fetest.dofs_per_cell;
224  Assert(result.size() == t_dofs,
225  ExcDimensionMismatch(result.size(), t_dofs));
226 
227  for (unsigned int k = 0; k < fetest.n_quadrature_points; ++k)
228  {
229  const double dx = factor * fetest.JxW(k);
230 
231  for (unsigned int i = 0; i < t_dofs; ++i)
232  for (unsigned int d = 0; d < dim; ++d)
233  result(i) +=
234  dx * input[k][d] * fetest.shape_value_component(i, k, d);
235  }
236  }
237 
250  template <int dim, typename number>
251  void
253  const FEValuesBase<dim> & fetest,
254  const std::vector<double> &input,
255  const double factor = 1.)
256  {
257  AssertDimension(fetest.get_fe().n_components(), dim);
258  AssertDimension(input.size(), fetest.n_quadrature_points);
259  const unsigned int t_dofs = fetest.dofs_per_cell;
260  Assert(result.size() == t_dofs,
261  ExcDimensionMismatch(result.size(), t_dofs));
262 
263  for (unsigned int k = 0; k < fetest.n_quadrature_points; ++k)
264  {
265  const double dx = factor * fetest.JxW(k);
266 
267  for (unsigned int i = 0; i < t_dofs; ++i)
268  for (unsigned int d = 0; d < dim; ++d)
269  result(i) -=
270  dx * input[k] * fetest.shape_grad_component(i, k, d)[d];
271  }
272  }
273 
282  template <int dim>
283  void
285  const FEValuesBase<dim> &fe,
286  const FEValuesBase<dim> &fetest,
287  double factor = 1.)
288  {
289  const unsigned int n_dofs = fe.dofs_per_cell;
290  const unsigned int t_dofs = fetest.dofs_per_cell;
291 
292  AssertDimension(fe.get_fe().n_components(), dim);
293  AssertDimension(fetest.get_fe().n_components(), 1);
294  AssertDimension(M.m(), t_dofs);
295  AssertDimension(M.n(), n_dofs);
296 
297  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
298  {
299  const Tensor<1, dim> ndx = factor * fe.JxW(k) * fe.normal_vector(k);
300  for (unsigned int i = 0; i < t_dofs; ++i)
301  for (unsigned int j = 0; j < n_dofs; ++j)
302  for (unsigned int d = 0; d < dim; ++d)
303  M(i, j) += ndx[d] * fe.shape_value_component(j, k, d) *
304  fetest.shape_value(i, k);
305  }
306  }
307 
318  template <int dim, typename number>
319  void
321  Vector<number> & result,
322  const FEValuesBase<dim> & fe,
323  const FEValuesBase<dim> & fetest,
324  const VectorSlice<const std::vector<std::vector<double>>> &data,
325  double factor = 1.)
326  {
327  const unsigned int t_dofs = fetest.dofs_per_cell;
328 
329  AssertDimension(fe.get_fe().n_components(), dim);
330  AssertDimension(fetest.get_fe().n_components(), 1);
331  AssertDimension(result.size(), t_dofs);
333 
334  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
335  {
336  const Tensor<1, dim> ndx = factor * fe.normal_vector(k) * fe.JxW(k);
337 
338  for (unsigned int i = 0; i < t_dofs; ++i)
339  for (unsigned int d = 0; d < dim; ++d)
340  result(i) += ndx[d] * fetest.shape_value(i, k) * data[d][k];
341  }
342  }
343 
354  template <int dim, typename number>
355  void
357  const FEValuesBase<dim> & fetest,
358  const std::vector<double> &data,
359  double factor = 1.)
360  {
361  const unsigned int t_dofs = fetest.dofs_per_cell;
362 
363  AssertDimension(fetest.get_fe().n_components(), dim);
364  AssertDimension(result.size(), t_dofs);
365  AssertDimension(data.size(), fetest.n_quadrature_points);
366 
367  for (unsigned int k = 0; k < fetest.n_quadrature_points; ++k)
368  {
369  const Tensor<1, dim> ndx =
370  factor * fetest.normal_vector(k) * fetest.JxW(k);
371 
372  for (unsigned int i = 0; i < t_dofs; ++i)
373  for (unsigned int d = 0; d < dim; ++d)
374  result(i) +=
375  ndx[d] * fetest.shape_value_component(i, k, d) * data[k];
376  }
377  }
378 
391  template <int dim>
392  void
394  FullMatrix<double> & M12,
395  FullMatrix<double> & M21,
396  FullMatrix<double> & M22,
397  const FEValuesBase<dim> &fe1,
398  const FEValuesBase<dim> &fe2,
399  const FEValuesBase<dim> &fetest1,
400  const FEValuesBase<dim> &fetest2,
401  double factor = 1.)
402  {
403  const unsigned int n_dofs = fe1.dofs_per_cell;
404  const unsigned int t_dofs = fetest1.dofs_per_cell;
405 
406  AssertDimension(fe1.get_fe().n_components(), dim);
407  AssertDimension(fe2.get_fe().n_components(), dim);
408  AssertDimension(fetest1.get_fe().n_components(), 1);
409  AssertDimension(fetest2.get_fe().n_components(), 1);
410  AssertDimension(M11.m(), t_dofs);
411  AssertDimension(M11.n(), n_dofs);
412  AssertDimension(M12.m(), t_dofs);
413  AssertDimension(M12.n(), n_dofs);
414  AssertDimension(M21.m(), t_dofs);
415  AssertDimension(M21.n(), n_dofs);
416  AssertDimension(M22.m(), t_dofs);
417  AssertDimension(M22.n(), n_dofs);
418 
419  for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
420  {
421  const double dx = factor * fe1.JxW(k);
422  for (unsigned int i = 0; i < t_dofs; ++i)
423  for (unsigned int j = 0; j < n_dofs; ++j)
424  for (unsigned int d = 0; d < dim; ++d)
425  {
426  const double un1 = fe1.shape_value_component(j, k, d) *
427  fe1.normal_vector(k)[d];
428  const double un2 = -fe2.shape_value_component(j, k, d) *
429  fe1.normal_vector(k)[d];
430  const double v1 = fetest1.shape_value(i, k);
431  const double v2 = fetest2.shape_value(i, k);
432 
433  M11(i, j) += .5 * dx * un1 * v1;
434  M12(i, j) += .5 * dx * un2 * v1;
435  M21(i, j) += .5 * dx * un1 * v2;
436  M22(i, j) += .5 * dx * un2 * v2;
437  }
438  }
439  }
440 
444  template <int dim>
445  DEAL_II_DEPRECATED void
447  const FEValuesBase<dim> &fe,
448  const double factor = 1.);
449 
450  template <int dim>
451  void
453  const FEValuesBase<dim> &fe,
454  const double factor)
455  {
456  GradDiv::cell_matrix(M, fe, factor);
457  }
458 
462  template <int dim, typename number>
463  DEAL_II_DEPRECATED void
465  Vector<number> & result,
466  const FEValuesBase<dim> & fetest,
467  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &input,
468  const double factor = 1.);
469 
470  template <int dim, typename number>
471  void
473  Vector<number> & result,
474  const FEValuesBase<dim> & fetest,
475  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &input,
476  const double factor)
477  {
478  GradDiv::cell_residual(result, fetest, input, factor);
479  }
480 
493  template <int dim>
494  void
496  FullMatrix<double> & M12,
497  FullMatrix<double> & M21,
498  FullMatrix<double> & M22,
499  const FEValuesBase<dim> &fe1,
500  const FEValuesBase<dim> &fe2,
501  double factor = 1.)
502  {
503  const unsigned int n_dofs = fe1.dofs_per_cell;
504 
505  AssertDimension(fe1.get_fe().n_components(), dim);
506  AssertDimension(fe2.get_fe().n_components(), dim);
507  AssertDimension(M11.m(), n_dofs);
508  AssertDimension(M11.n(), n_dofs);
509  AssertDimension(M12.m(), n_dofs);
510  AssertDimension(M12.n(), n_dofs);
511  AssertDimension(M21.m(), n_dofs);
512  AssertDimension(M21.n(), n_dofs);
513  AssertDimension(M22.m(), n_dofs);
514  AssertDimension(M22.n(), n_dofs);
515 
516  for (unsigned int k = 0; k < fe1.n_quadrature_points; ++k)
517  {
518  const double dx = factor * fe1.JxW(k);
519  for (unsigned int i = 0; i < n_dofs; ++i)
520  for (unsigned int j = 0; j < n_dofs; ++j)
521  for (unsigned int d = 0; d < dim; ++d)
522  {
523  const double un1 = fe1.shape_value_component(j, k, d) *
524  fe1.normal_vector(k)[d];
525  const double un2 = -fe2.shape_value_component(j, k, d) *
526  fe1.normal_vector(k)[d];
527  const double vn1 = fe1.shape_value_component(i, k, d) *
528  fe1.normal_vector(k)[d];
529  const double vn2 = -fe2.shape_value_component(i, k, d) *
530  fe1.normal_vector(k)[d];
531 
532  M11(i, j) += dx * un1 * vn1;
533  M12(i, j) += dx * un2 * vn1;
534  M21(i, j) += dx * un1 * vn2;
535  M22(i, j) += dx * un2 * vn2;
536  }
537  }
538  }
539 
551  template <int dim>
552  double
554  const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &Du)
555  {
556  AssertDimension(fe.get_fe().n_components(), dim);
558 
559  double result = 0;
560  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
561  {
562  double div = Du[0][k][0];
563  for (unsigned int d = 1; d < dim; ++d)
564  div += Du[d][k][d];
565  result += div * div * fe.JxW(k);
566  }
567  return result;
568  }
569 
570  } // namespace Divergence
571 } // namespace LocalIntegrators
572 
573 
574 DEAL_II_NAMESPACE_CLOSE
575 
576 #endif
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
const Tensor< 1, spacedim > & shape_grad(const unsigned int function_no, const unsigned int quadrature_point) const
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1366
void u_dot_n_residual(Vector< number > &result, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const VectorSlice< const std::vector< std::vector< double >>> &data, double factor=1.)
Definition: divergence.h:320
const unsigned int dofs_per_cell
Definition: fe_values.h:2033
void grad_div_residual(Vector< number > &result, const FEValuesBase< dim > &fetest, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim >>>> &input, const double factor=1.)
Definition: divergence.h:472
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
size_type size() const
#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 cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.)
Definition: divergence.h:59
void gradient_residual(Vector< number > &result, const FEValuesBase< dim > &fetest, const std::vector< Tensor< 1, dim >> &input, const double factor=1.)
Definition: divergence.h:216
#define Assert(cond, exc)
Definition: exceptions.h:1227
static::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
const FiniteElement< dim, spacedim > & get_fe() const
void gradient_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.)
Definition: divergence.h:173
void u_dot_n_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.)
Definition: divergence.h:284
size_type m() const
void u_dot_n_jump_matrix(FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, double factor=1.)
Definition: divergence.h:495
const unsigned int n_quadrature_points
Definition: fe_values.h:2026
void u_times_n_residual(Vector< number > &result, const FEValuesBase< dim > &fetest, const std::vector< double > &data, double factor=1.)
Definition: divergence.h:356
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
void grad_div_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
Definition: divergence.h:452
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: divergence.h:101
double norm(const FEValuesBase< dim > &fe, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim >>>> &Du)
Definition: divergence.h:553
Tensor< 1, spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
double JxW(const unsigned int quadrature_point) const