Reference documentation for deal.II version 8.4.2
divergence.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_divergence_h
17 #define dealii__integrators_divergence_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 {
40  namespace Divergence
41  {
52  template <int dim>
55  const Tensor<2,dim> &h0,
56  const Tensor<2,dim> &h1,
57  const Tensor<2,dim> &h2)
58  {
59  Tensor<1,dim> result;
60  for (unsigned int d=0; d<dim; ++d)
61  {
62  result[d] += h0[d][0];
63  if (dim >=2) result[d] += h1[d][1];
64  if (dim >=3) result[d] += h2[d][2];
65  }
66  return result;
67  }
68 
69 
79  template <int dim>
80  void cell_matrix (
82  const FEValuesBase<dim> &fe,
83  const FEValuesBase<dim> &fetest,
84  double factor = 1.)
85  {
86  unsigned int fecomp = fe.get_fe().n_components();
87  const unsigned int n_dofs = fe.dofs_per_cell;
88  const unsigned int t_dofs = fetest.dofs_per_cell;
89  AssertDimension(fecomp, dim);
90  AssertDimension(fetest.get_fe().n_components(), 1);
91  AssertDimension(M.m(), t_dofs);
92  AssertDimension(M.n(), n_dofs);
93 
94  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
95  {
96  const double dx = fe.JxW(k) * factor;
97  for (unsigned int i=0; i<t_dofs; ++i)
98  {
99  const double vv = fetest.shape_value(i,k);
100  for (unsigned int d=0; d<dim; ++d)
101  for (unsigned int j=0; j<n_dofs; ++j)
102  {
103  const double du = fe.shape_grad_component(j,k,d)[d];
104  M(i,j) += dx * du * vv;
105  }
106  }
107  }
108  }
109 
122  template <int dim, typename number>
124  Vector<number> &result,
125  const FEValuesBase<dim> &fetest,
126  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &input,
127  const double factor = 1.)
128  {
129  AssertDimension(fetest.get_fe().n_components(), 1);
131  const unsigned int t_dofs = fetest.dofs_per_cell;
132  Assert (result.size() == t_dofs, ExcDimensionMismatch(result.size(), t_dofs));
133 
134  for (unsigned int k=0; k<fetest.n_quadrature_points; ++k)
135  {
136  const double dx = factor * fetest.JxW(k);
137 
138  for (unsigned int i=0; i<t_dofs; ++i)
139  for (unsigned int d=0; d<dim; ++d)
140  result(i) += dx * input[d][k][d] * fetest.shape_value(i,k);
141  }
142  }
143 
144 
157  template <int dim, typename number>
159  Vector<number> &result,
160  const FEValuesBase<dim> &fetest,
161  const VectorSlice<const std::vector<std::vector<double> > > &input,
162  const double factor = 1.)
163  {
164  AssertDimension(fetest.get_fe().n_components(), 1);
166  const unsigned int t_dofs = fetest.dofs_per_cell;
167  Assert (result.size() == t_dofs, ExcDimensionMismatch(result.size(), t_dofs));
168 
169  for (unsigned int k=0; k<fetest.n_quadrature_points; ++k)
170  {
171  const double dx = factor * fetest.JxW(k);
172 
173  for (unsigned int i=0; i<t_dofs; ++i)
174  for (unsigned int d=0; d<dim; ++d)
175  result(i) -= dx * input[d][k] * fetest.shape_grad(i,k)[d];
176  }
177  }
178 
179 
190  template <int dim>
193  const FEValuesBase<dim> &fe,
194  const FEValuesBase<dim> &fetest,
195  double factor = 1.)
196  {
197  unsigned int fecomp = fetest.get_fe().n_components();
198  const unsigned int t_dofs = fetest.dofs_per_cell;
199  const unsigned int n_dofs = fe.dofs_per_cell;
200 
201  AssertDimension(fecomp, dim);
202  AssertDimension(fe.get_fe().n_components(), 1);
203  AssertDimension(M.m(), t_dofs);
204  AssertDimension(M.n(), n_dofs);
205 
206  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
207  {
208  const double dx = fe.JxW(k) * factor;
209  for (unsigned int d=0; d<dim; ++d)
210  for (unsigned int i=0; i<t_dofs; ++i)
211  {
212  const double vv = fetest.shape_value_component(i,k,d);
213  for (unsigned int j=0; j<n_dofs; ++j)
214  {
215  const Tensor<1,dim> &Du = fe.shape_grad(j,k);
216  M(i,j) += dx * vv * Du[d];
217  }
218  }
219  }
220  }
221 
234  template <int dim, typename number>
236  Vector<number> &result,
237  const FEValuesBase<dim> &fetest,
238  const std::vector<Tensor<1,dim> > &input,
239  const double factor = 1.)
240  {
241  AssertDimension(fetest.get_fe().n_components(), dim);
242  AssertDimension(input.size(), fetest.n_quadrature_points);
243  const unsigned int t_dofs = fetest.dofs_per_cell;
244  Assert (result.size() == t_dofs, ExcDimensionMismatch(result.size(), t_dofs));
245 
246  for (unsigned int k=0; k<fetest.n_quadrature_points; ++k)
247  {
248  const double dx = factor * fetest.JxW(k);
249 
250  for (unsigned int i=0; i<t_dofs; ++i)
251  for (unsigned int d=0; d<dim; ++d)
252  result(i) += dx * input[k][d] * fetest.shape_value_component(i,k,d);
253  }
254  }
255 
268  template <int dim, typename number>
270  Vector<number> &result,
271  const FEValuesBase<dim> &fetest,
272  const std::vector<double> &input,
273  const double factor = 1.)
274  {
275  AssertDimension(fetest.get_fe().n_components(), dim);
276  AssertDimension(input.size(), fetest.n_quadrature_points);
277  const unsigned int t_dofs = fetest.dofs_per_cell;
278  Assert (result.size() == t_dofs, ExcDimensionMismatch(result.size(), t_dofs));
279 
280  for (unsigned int k=0; k<fetest.n_quadrature_points; ++k)
281  {
282  const double dx = factor * fetest.JxW(k);
283 
284  for (unsigned int i=0; i<t_dofs; ++i)
285  for (unsigned int d=0; d<dim; ++d)
286  result(i) -= dx * input[k] * fetest.shape_grad_component(i,k,d)[d];
287  }
288  }
289 
298  template<int dim>
299  void
302  const FEValuesBase<dim> &fe,
303  const FEValuesBase<dim> &fetest,
304  double factor = 1.)
305  {
306  unsigned int fecomp = fe.get_fe().n_components();
307  const unsigned int n_dofs = fe.dofs_per_cell;
308  const unsigned int t_dofs = fetest.dofs_per_cell;
309 
310  AssertDimension(fecomp, dim);
311  AssertDimension(fetest.get_fe().n_components(), 1);
312  AssertDimension(M.m(), t_dofs);
313  AssertDimension(M.n(), n_dofs);
314 
315  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
316  {
317  const Tensor<1,dim> ndx = factor * fe.JxW(k) * fe.normal_vector(k);
318  for (unsigned int i=0; i<t_dofs; ++i)
319  for (unsigned int j=0; j<n_dofs; ++j)
320  for (unsigned int d=0; d<dim; ++d)
321  M(i,j) += ndx[d] * fe.shape_value_component(j,k,d)
322  * fetest.shape_value(i,k);
323  }
324  }
325 
336  template<int dim, typename number>
337  void
339  Vector<number> &result,
340  const FEValuesBase<dim> &fe,
341  const FEValuesBase<dim> &fetest,
342  const VectorSlice<const std::vector<std::vector<double> > > &data,
343  double factor = 1.)
344  {
345  unsigned int fecomp = fe.get_fe().n_components();
346  const unsigned int t_dofs = fetest.dofs_per_cell;
347 
348  AssertDimension(fecomp, dim);
349  AssertDimension(fetest.get_fe().n_components(), 1);
350  AssertDimension(result.size(), t_dofs);
352 
353  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
354  {
355  const Tensor<1,dim> ndx = factor * fe.normal_vector(k) * fe.JxW(k);
356 
357  for (unsigned int i=0; i<t_dofs; ++i)
358  for (unsigned int d=0; d<dim; ++d)
359  result(i) += ndx[d] * fetest.shape_value(i,k) * data[d][k];
360  }
361  }
362 
373  template<int dim, typename number>
374  void
376  Vector<number> &result,
377  const FEValuesBase<dim> &fetest,
378  const std::vector<double> &data,
379  double factor = 1.)
380  {
381  const unsigned int t_dofs = fetest.dofs_per_cell;
382 
383  AssertDimension(fetest.get_fe().n_components(), dim);
384  AssertDimension(result.size(), t_dofs);
385  AssertDimension(data.size(), fetest.n_quadrature_points);
386 
387  for (unsigned int k=0; k<fetest.n_quadrature_points; ++k)
388  {
389  const Tensor<1,dim> ndx = factor * fetest.normal_vector(k) * fetest.JxW(k);
390 
391  for (unsigned int i=0; i<t_dofs; ++i)
392  for (unsigned int d=0; d<dim; ++d)
393  result(i) += ndx[d] * fetest.shape_value_component(i,k,d) * data[k];
394  }
395  }
396 
408  template<int dim>
409  void
411  FullMatrix<double> &M11,
412  FullMatrix<double> &M12,
413  FullMatrix<double> &M21,
414  FullMatrix<double> &M22,
415  const FEValuesBase<dim> &fe1,
416  const FEValuesBase<dim> &fe2,
417  const FEValuesBase<dim> &fetest1,
418  const FEValuesBase<dim> &fetest2,
419  double factor = 1.)
420  {
421  const unsigned int n_dofs = fe1.dofs_per_cell;
422  const unsigned int t_dofs = fetest1.dofs_per_cell;
423 
424  AssertDimension(fe1.get_fe().n_components(), dim);
425  AssertDimension(fe2.get_fe().n_components(), dim);
426  AssertDimension(fetest1.get_fe().n_components(), 1);
427  AssertDimension(fetest2.get_fe().n_components(), 1);
428  AssertDimension(M11.m(), t_dofs);
429  AssertDimension(M11.n(), n_dofs);
430  AssertDimension(M12.m(), t_dofs);
431  AssertDimension(M12.n(), n_dofs);
432  AssertDimension(M21.m(), t_dofs);
433  AssertDimension(M21.n(), n_dofs);
434  AssertDimension(M22.m(), t_dofs);
435  AssertDimension(M22.n(), n_dofs);
436 
437  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
438  {
439  const double dx = factor * fe1.JxW(k);
440  for (unsigned int i=0; i<t_dofs; ++i)
441  for (unsigned int j=0; j<n_dofs; ++j)
442  for (unsigned int d=0; d<dim; ++d)
443  {
444  const double un1 = fe1.shape_value_component(j,k,d) * fe1.normal_vector(k)[d];
445  const double un2 =-fe2.shape_value_component(j,k,d) * fe1.normal_vector(k)[d];
446  const double v1 = fetest1.shape_value(i,k);
447  const double v2 = fetest2.shape_value(i,k);
448 
449  M11(i,j) += .5 * dx * un1 * v1;
450  M12(i,j) += .5 * dx * un2 * v1;
451  M21(i,j) += .5 * dx * un1 * v2;
452  M22(i,j) += .5 * dx * un2 * v2;
453  }
454  }
455  }
456 
466  template <int dim>
469  const FEValuesBase<dim> &fe,
470  double factor = 1.)
471  {
472  const unsigned int n_dofs = fe.dofs_per_cell;
473 
474  AssertDimension(fe.get_fe().n_components(), dim);
475  AssertDimension(M.m(), n_dofs);
476  AssertDimension(M.n(), n_dofs);
477 
478  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
479  {
480  const double dx = factor * fe.JxW(k);
481  for (unsigned int i=0; i<n_dofs; ++i)
482  for (unsigned int j=0; j<n_dofs; ++j)
483  {
484  double dv = 0.;
485  double du = 0.;
486  for (unsigned int d=0; d<dim; ++d)
487  {
488  dv += fe.shape_grad_component(i,k,d)[d];
489  du += fe.shape_grad_component(j,k,d)[d];
490  }
491 
492  M(i,j) += dx * du * dv;
493  }
494  }
495  }
496 
506  template <int dim, typename number>
508  Vector<number> &result,
509  const FEValuesBase<dim> &fetest,
510  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &input,
511  const double factor = 1.)
512  {
513  const unsigned int n_dofs = fetest.dofs_per_cell;
514 
515  AssertDimension(fetest.get_fe().n_components(), dim);
517 
518  for (unsigned int k=0; k<fetest.n_quadrature_points; ++k)
519  {
520  const double dx = factor * fetest.JxW(k);
521  for (unsigned int i=0; i<n_dofs; ++i)
522  {
523  double dv = 0.;
524  double du = 0.;
525  for (unsigned int d=0; d<dim; ++d)
526  {
527  dv += fetest.shape_grad_component(i,k,d)[d];
528  du += input[d][k][d];
529  }
530 
531  result(i) += dx * du * dv;
532  }
533  }
534  }
535 
548  template<int dim>
549  void
551  FullMatrix<double> &M11,
552  FullMatrix<double> &M12,
553  FullMatrix<double> &M21,
554  FullMatrix<double> &M22,
555  const FEValuesBase<dim> &fe1,
556  const FEValuesBase<dim> &fe2,
557  double factor = 1.)
558  {
559  const unsigned int n_dofs = fe1.dofs_per_cell;
560 
561  AssertDimension(fe1.get_fe().n_components(), dim);
562  AssertDimension(fe2.get_fe().n_components(), dim);
563  AssertDimension(M11.m(), n_dofs);
564  AssertDimension(M11.n(), n_dofs);
565  AssertDimension(M12.m(), n_dofs);
566  AssertDimension(M12.n(), n_dofs);
567  AssertDimension(M21.m(), n_dofs);
568  AssertDimension(M21.n(), n_dofs);
569  AssertDimension(M22.m(), n_dofs);
570  AssertDimension(M22.n(), n_dofs);
571 
572  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
573  {
574  const double dx = factor * fe1.JxW(k);
575  for (unsigned int i=0; i<n_dofs; ++i)
576  for (unsigned int j=0; j<n_dofs; ++j)
577  for (unsigned int d=0; d<dim; ++d)
578  {
579  const double un1 = fe1.shape_value_component(j,k,d) * fe1.normal_vector(k)[d];
580  const double un2 =-fe2.shape_value_component(j,k,d) * fe1.normal_vector(k)[d];
581  const double vn1 = fe1.shape_value_component(i,k,d) * fe1.normal_vector(k)[d];
582  const double vn2 =-fe2.shape_value_component(i,k,d) * fe1.normal_vector(k)[d];
583 
584  M11(i,j) += dx * un1 * vn1;
585  M12(i,j) += dx * un2 * vn1;
586  M21(i,j) += dx * un1 * vn2;
587  M22(i,j) += dx * un2 * vn2;
588  }
589  }
590  }
591 
603  template <int dim>
604  double norm(const FEValuesBase<dim> &fe,
605  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Du)
606  {
607  unsigned int fecomp = fe.get_fe().n_components();
608 
609  AssertDimension(fecomp, dim);
611 
612  double result = 0;
613  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
614  {
615  double div = Du[0][k][0];
616  for (unsigned int d=1; d<dim; ++d)
617  div += Du[d][k][d];
618  result += div*div*fe.JxW(k);
619  }
620  return result;
621  }
622 
623  }
624 }
625 
626 
627 DEAL_II_NAMESPACE_CLOSE
628 
629 #endif
size_type m() const
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1052
const unsigned int dofs_per_cell
Definition: fe_values.h:1464
#define AssertVectorVectorDimension(vec, dim1, dim2)
Definition: exceptions.h:1063
const FiniteElement< dim, spacedim > & get_fe() const
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:123
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:507
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 cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.)
Definition: divergence.h:80
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
void gradient_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.)
Definition: divergence.h:191
void u_dot_n_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.)
Definition: divergence.h:300
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:550
const unsigned int n_quadrature_points
Definition: fe_values.h:1457
void u_times_n_residual(Vector< number > &result, const FEValuesBase< dim > &fetest, const std::vector< double > &data, double factor=1.)
Definition: divergence.h:375
double norm(const FEValuesBase< dim > &fe, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim > > > > &Du)
Definition: divergence.h:604
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:338
void gradient_residual(Vector< number > &result, const FEValuesBase< dim > &fetest, const std::vector< Tensor< 1, dim > > &input, const double factor=1.)
Definition: divergence.h:235
double JxW(const unsigned int quadrature_point) const
Tensor< 1, dim > grad_div(const Tensor< 2, dim > &h0, const Tensor< 2, dim > &h1, const Tensor< 2, dim > &h2)
Definition: divergence.h:54
void grad_div_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double factor=1.)
Definition: divergence.h:467
const Tensor< 1, spacedim > & shape_grad(const unsigned int function_no, const unsigned int quadrature_point) const
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