Reference documentation for deal.II version 8.5.1
advection.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2016 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_advection_h
17 #define dealii__integrators_advection_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 
31 {
48  namespace Advection
49  {
74  template<int dim>
75  void cell_matrix (
77  const FEValuesBase<dim> &fe,
78  const FEValuesBase<dim> &fetest,
79  const VectorSlice<const std::vector<std::vector<double> > > &velocity,
80  const double factor = 1.)
81  {
82  const unsigned int n_dofs = fe.dofs_per_cell;
83  const unsigned int t_dofs = fetest.dofs_per_cell;
84  const unsigned int n_components = fe.get_fe().n_components();
85 
86  AssertDimension(velocity.size(), dim);
87  // If the size of the
88  // velocity vectors is one,
89  // then do not increment
90  // between quadrature points.
91  const unsigned int v_increment = (velocity[0].size() == 1) ? 0 : 1;
92 
93  if (v_increment == 1)
94  {
96  }
97 
98  AssertDimension(M.n(), n_dofs);
99  AssertDimension(M.m(), t_dofs);
100 
101  for (unsigned k=0; k<fe.n_quadrature_points; ++k)
102  {
103  const double dx = factor * fe.JxW(k);
104  const unsigned int vindex = k * v_increment;
105 
106  for (unsigned j=0; j<n_dofs; ++j)
107  for (unsigned i=0; i<t_dofs; ++i)
108  for (unsigned int c=0; c<n_components; ++c)
109  {
110  double wgradv = velocity[0][vindex]
111  * fe.shape_grad_component(i,k,c)[0];
112  for (unsigned int d=1; d<dim; ++d)
113  wgradv += velocity[d][vindex]
114  * fe.shape_grad_component(i,k,c)[d];
115  M(i,j) -= dx * wgradv * fe.shape_value_component(j,k,c);
116  }
117  }
118  }
119 
120 
121 
130  template <int dim>
131  inline void
133  Vector<double> &result,
134  const FEValuesBase<dim> &fe,
135  const std::vector<Tensor<1,dim> > &input,
136  const VectorSlice<const std::vector<std::vector<double> > > &velocity,
137  double factor = 1.)
138  {
139  const unsigned int nq = fe.n_quadrature_points;
140  const unsigned int n_dofs = fe.dofs_per_cell;
141  Assert(input.size() == nq, ExcDimensionMismatch(input.size(), nq));
142  Assert(result.size() == n_dofs, ExcDimensionMismatch(result.size(), n_dofs));
143 
144  AssertDimension(velocity.size(), dim);
145  const unsigned int v_increment = (velocity[0].size() == 1) ? 0 : 1;
146  if (v_increment == 1)
147  {
149  }
150 
151  for (unsigned k=0; k<nq; ++k)
152  {
153  const double dx = factor * fe.JxW(k);
154  for (unsigned i=0; i<n_dofs; ++i)
155  for (unsigned int d=0; d<dim; ++d)
156  result(i) += dx * input[k][d]
157  * fe.shape_value(i,k) * velocity[d][k * v_increment];
158  }
159  }
160 
161 
162 
173  template <int dim>
174  inline void
176  Vector<double> &result,
177  const FEValuesBase<dim> &fe,
178  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &input,
179  const VectorSlice<const std::vector<std::vector<double> > > &velocity,
180  double factor = 1.)
181  {
182  const unsigned int nq = fe.n_quadrature_points;
183  const unsigned int n_dofs = fe.dofs_per_cell;
184  const unsigned int n_comp = fe.get_fe().n_components();
185 
187  Assert(result.size() == n_dofs, ExcDimensionMismatch(result.size(), n_dofs));
188 
189  AssertDimension(velocity.size(), dim);
190  const unsigned int v_increment = (velocity[0].size() == 1) ? 0 : 1;
191  if (v_increment == 1)
192  {
194  }
195 
196  for (unsigned k=0; k<nq; ++k)
197  {
198  const double dx = factor * fe.JxW(k);
199  for (unsigned i=0; i<n_dofs; ++i)
200  for (unsigned int c=0; c<n_comp; ++c)
201  for (unsigned int d=0; d<dim; ++d)
202  result(i) += dx * input[c][k][d]
203  * fe.shape_value_component(i,k,c) * velocity[d][k * v_increment];
204  }
205  }
206 
207 
208 
214  template <int dim>
215  inline void
217  Vector<double> &result,
218  const FEValuesBase<dim> &fe,
219  const std::vector<double> &input,
220  const VectorSlice<const std::vector<std::vector<double> > > &velocity,
221  double factor = 1.)
222  {
223  const unsigned int nq = fe.n_quadrature_points;
224  const unsigned int n_dofs = fe.dofs_per_cell;
225  Assert(input.size() == nq, ExcDimensionMismatch(input.size(), nq));
226  Assert(result.size() == n_dofs, ExcDimensionMismatch(result.size(), n_dofs));
227 
228  AssertDimension(velocity.size(), dim);
229  const unsigned int v_increment = (velocity[0].size() == 1) ? 0 : 1;
230  if (v_increment == 1)
231  {
233  }
234 
235  for (unsigned k=0; k<nq; ++k)
236  {
237  const double dx = factor * fe.JxW(k);
238  for (unsigned i=0; i<n_dofs; ++i)
239  for (unsigned int d=0; d<dim; ++d)
240  result(i) -= dx * input[k]
241  * fe.shape_grad(i,k)[d] * velocity[d][k * v_increment];
242  }
243  }
244 
245 
246 
254  template <int dim>
255  inline void
257  Vector<double> &result,
258  const FEValuesBase<dim> &fe,
259  const VectorSlice<const std::vector<std::vector<double> > > &input,
260  const VectorSlice<const std::vector<std::vector<double> > > &velocity,
261  double factor = 1.)
262  {
263  const unsigned int nq = fe.n_quadrature_points;
264  const unsigned int n_dofs = fe.dofs_per_cell;
265  const unsigned int n_comp = fe.get_fe().n_components();
266 
268  Assert(result.size() == n_dofs, ExcDimensionMismatch(result.size(), n_dofs));
269 
270  AssertDimension(velocity.size(), dim);
271  const unsigned int v_increment = (velocity[0].size() == 1) ? 0 : 1;
272  if (v_increment == 1)
273  {
275  }
276 
277  for (unsigned k=0; k<nq; ++k)
278  {
279  const double dx = factor * fe.JxW(k);
280  for (unsigned i=0; i<n_dofs; ++i)
281  for (unsigned int c=0; c<n_comp; ++c)
282  for (unsigned int d=0; d<dim; ++d)
283  result(i) -= dx * input[c][k]
284  * fe.shape_grad_component(i,k,c)[d] * velocity[d][k * v_increment];
285  }
286  }
287 
288 
289 
307  template <int dim>
310  const FEValuesBase<dim> &fe,
311  const FEValuesBase<dim> &fetest,
312  const VectorSlice<const std::vector<std::vector<double> > > &velocity,
313  double factor = 1.)
314  {
315  const unsigned int n_dofs = fe.dofs_per_cell;
316  const unsigned int t_dofs = fetest.dofs_per_cell;
317  unsigned int n_components = fe.get_fe().n_components();
318  AssertDimension (M.m(), n_dofs);
319  AssertDimension (M.n(), n_dofs);
320 
321  AssertDimension(velocity.size(), dim);
322  const unsigned int v_increment = (velocity[0].size() == 1) ? 0 : 1;
323  if (v_increment == 1)
324  {
326  }
327 
328  for (unsigned k=0; k<fe.n_quadrature_points; ++k)
329  {
330  const double dx = factor * fe.JxW(k);
331 
332  double nv = 0.;
333  for (unsigned int d=0; d<dim; ++d)
334  nv += fe.normal_vector(k)[d] * velocity[d][k * v_increment];
335 
336  if (nv > 0)
337  {
338  for (unsigned i=0; i<t_dofs; ++i)
339  for (unsigned j=0; j<n_dofs; ++j)
340  {
341  if (fe.get_fe().is_primitive())
342  M(i,j) += dx * nv * fe.shape_value(i,k) * fe.shape_value(j,k);
343  else
344  for (unsigned int c=0; c<n_components; ++c)
345  M(i,j) += dx * nv * fetest.shape_value_component(i,k,c)
346  * fe.shape_value_component(j,k,c);
347  }
348  }
349  }
350  }
351 
352 
353 
378  template <int dim>
379  inline void
381  Vector<double> &result,
382  const FEValuesBase<dim> &fe,
383  const std::vector<double> &input,
384  const std::vector<double> &data,
385  const VectorSlice<const std::vector<std::vector<double> > > &velocity,
386  double factor = 1.)
387  {
388  const unsigned int n_dofs = fe.dofs_per_cell;
389 
390  AssertDimension(input.size(), fe.n_quadrature_points);
391  AssertDimension(data.size(), fe.n_quadrature_points);
392 
393  AssertDimension(velocity.size(), dim);
394  const unsigned int v_increment = (velocity[0].size() == 1) ? 0 : 1;
395  if (v_increment == 1)
396  {
398  }
399 
400 
401  for (unsigned k=0; k<fe.n_quadrature_points; ++k)
402  {
403  const double dx = factor * fe.JxW(k);
404 
405  double nv = 0.;
406  for (unsigned int d=0; d<dim; ++d)
407  nv += fe.normal_vector(k)[d] * velocity[d][k * v_increment];
408 
409  // Always use the upwind value
410  const double val = (nv > 0.) ? input[k] : -data[k];
411 
412  for (unsigned i=0; i<n_dofs; ++i)
413  {
414  const double v= fe.shape_value(i,k);
415  result(i) += dx * nv * val *v;
416  }
417  }
418  }
419 
420 
421 
446  template <int dim>
447  inline void
449  Vector<double> &result,
450  const FEValuesBase<dim> &fe,
451  const VectorSlice<const std::vector<std::vector<double> > > &input,
452  const VectorSlice<const std::vector<std::vector<double> > > &data,
453  const VectorSlice<const std::vector<std::vector<double> > > &velocity,
454  double factor = 1.)
455  {
456  const unsigned int n_dofs = fe.dofs_per_cell;
457  const unsigned int n_comp = fe.get_fe().n_components();
458 
461 
462  AssertDimension(velocity.size(), dim);
463  const unsigned int v_increment = (velocity[0].size() == 1) ? 0 : 1;
464  if (v_increment == 1)
465  {
467  }
468 
469 
470  for (unsigned k=0; k<fe.n_quadrature_points; ++k)
471  {
472  const double dx = factor * fe.JxW(k);
473 
474  double nv = 0.;
475  for (unsigned int d=0; d<dim; ++d)
476  nv += fe.normal_vector(k)[d] * velocity[d][k * v_increment];
477 
478  std::vector<double> val(n_comp);
479 
480  for (unsigned int d=0; d<n_comp; ++d)
481  {
482  val[d] = (nv > 0.) ? input[d][k] : -data[d][k];
483  for (unsigned i=0; i<n_dofs; ++i)
484  {
485  const double v= fe.shape_value_component(i,k,d);
486  result(i) += dx * nv * val[d] *v;
487  }
488  }
489  }
490  }
491 
492 
493 
514  template <int dim>
516  FullMatrix<double> &M11,
517  FullMatrix<double> &M12,
518  FullMatrix<double> &M21,
519  FullMatrix<double> &M22,
520  const FEValuesBase<dim> &fe1,
521  const FEValuesBase<dim> &fe2,
522  const FEValuesBase<dim> &fetest1,
523  const FEValuesBase<dim> &fetest2,
524  const VectorSlice<const std::vector<std::vector<double> > > &velocity,
525  const double factor = 1.)
526  {
527  const unsigned int n1 = fe1.dofs_per_cell;
528  // Multiply the quadrature point
529  // index below with this factor to
530  // have simpler data for constant
531  // velocities.
532  AssertDimension(velocity.size(), dim);
533  const unsigned int v_increment = (velocity[0].size() == 1) ? 0 : 1;
534  if (v_increment == 1)
535  {
537  }
538 
539  for (unsigned k=0; k<fe1.n_quadrature_points; ++k)
540  {
541  double nbeta = fe1.normal_vector(k)[0] * velocity[0][k * v_increment];
542  for (unsigned int d=1; d<dim; ++d)
543  nbeta += fe1.normal_vector(k)[d] * velocity[d][k * v_increment];
544  const double dx_nbeta = factor * std::abs(nbeta) * fe1.JxW(k);
545  FullMatrix<double> &M1 = nbeta > 0. ? M11 : M22;
546  FullMatrix<double> &M2 = nbeta > 0. ? M21 : M12;
547  const FEValuesBase<dim> &fe = nbeta > 0. ? fe1 : fe2;
548  const FEValuesBase<dim> &fetest = nbeta > 0. ? fetest1 : fetest2;
549  const FEValuesBase<dim> &fetestn = nbeta > 0. ? fetest2 : fetest1;
550  for (unsigned i=0; i<n1; ++i)
551  for (unsigned j=0; j<n1; ++j)
552  {
553  if (fe1.get_fe().is_primitive())
554  {
555  M1(i,j) += dx_nbeta*fe.shape_value(j,k)*fetest.shape_value(i,k);
556  M2(i,j) -= dx_nbeta*fe.shape_value(j,k)*fetestn.shape_value(i,k);
557  }
558  else
559  {
560  for (unsigned int d=0; d<fe1.get_fe().n_components(); ++d)
561  {
562  M1(i,j) += dx_nbeta*fe.shape_value_component(j,k,d)*fetest.shape_value_component(i,k,d);
563  M2(i,j) -= dx_nbeta*fe.shape_value_component(j,k,d)*fetestn.shape_value_component(i,k,d);
564  }
565  }
566  }
567  }
568  }
569 
570 
571 
592  template <int dim>
594  Vector<double> &result1,
595  Vector<double> &result2,
596  const FEValuesBase<dim> &fe1,
597  const FEValuesBase<dim> &fe2,
598  const std::vector<double> &input1,
599  const std::vector<double> &input2,
600  const VectorSlice<const std::vector<std::vector<double> > > &velocity,
601  const double factor = 1.)
602  {
603  Assert(fe1.get_fe().n_components() == 1,
604  ExcDimensionMismatch(fe1.get_fe().n_components(), 1));
605  Assert(fe2.get_fe().n_components() == 1,
606  ExcDimensionMismatch(fe2.get_fe().n_components(), 1));
607 
608  const unsigned int n1 = fe1.dofs_per_cell;
609  // Multiply the quadrature point
610  // index below with this factor to
611  // have simpler data for constant
612  // velocities.
613  AssertDimension(velocity.size(), dim);
614  const unsigned int v_increment = (velocity[0].size() == 1) ? 0 : 1;
615  if (v_increment == 1)
616  {
618  }
619 
620  for (unsigned k=0; k<fe1.n_quadrature_points; ++k)
621  {
622  double nbeta = fe1.normal_vector(k)[0] * velocity[0][k * v_increment];
623  for (unsigned int d=1; d<dim; ++d)
624  nbeta += fe1.normal_vector(k)[d] * velocity[d][k * v_increment];
625  const double dx_nbeta = factor * nbeta * fe1.JxW(k);
626 
627  for (unsigned i=0; i<n1; ++i)
628  {
629  const double v1 = fe1.shape_value(i,k);
630  const double v2 = fe2.shape_value(i,k);
631  const double u1 = input1[k];
632  const double u2 = input2[k];
633  if (nbeta > 0)
634  {
635  result1(i) += dx_nbeta*u1*v1;
636  result2(i) -= dx_nbeta*u1*v2;
637  }
638  else
639  {
640 
641  result1(i) += dx_nbeta*u2*v1;
642  result2(i) -= dx_nbeta*u2*v2;
643  }
644  }
645  }
646  }
647 
648 
649 
670  template <int dim>
672  Vector<double> &result1,
673  Vector<double> &result2,
674  const FEValuesBase<dim> &fe1,
675  const FEValuesBase<dim> &fe2,
676  const VectorSlice<const std::vector<std::vector<double> > > &input1,
677  const VectorSlice<const std::vector<std::vector<double> > > &input2,
678  const VectorSlice<const std::vector<std::vector<double> > > &velocity,
679  const double factor = 1.)
680  {
681  const unsigned int n_comp = fe1.get_fe().n_components();
682  const unsigned int n1 = fe1.dofs_per_cell;
685 
686  // Multiply the quadrature point
687  // index below with this factor to
688  // have simpler data for constant
689  // velocities.
690  AssertDimension(velocity.size(), dim);
691  const unsigned int v_increment = (velocity[0].size() == 1) ? 0 : 1;
692  if (v_increment == 1)
693  {
695  }
696 
697  for (unsigned k=0; k<fe1.n_quadrature_points; ++k)
698  {
699  double nbeta = fe1.normal_vector(k)[0] * velocity[0][k * v_increment];
700  for (unsigned int d=1; d<dim; ++d)
701  nbeta += fe1.normal_vector(k)[d] * velocity[d][k * v_increment];
702  const double dx_nbeta = factor * nbeta * fe1.JxW(k);
703 
704  for (unsigned i=0; i<n1; ++i)
705  for (unsigned int d=0; d<n_comp; ++d)
706  {
707  const double v1 = fe1.shape_value_component(i,k,d);
708  const double v2 = fe2.shape_value_component(i,k,d);
709  const double u1 = input1[d][k];
710  const double u2 = input2[d][k];
711  if (nbeta > 0)
712  {
713  result1(i) += dx_nbeta*u1*v1;
714  result2(i) -= dx_nbeta*u1*v2;
715  }
716  else
717  {
718 
719  result1(i) += dx_nbeta*u2*v1;
720  result2(i) -= dx_nbeta*u2*v2;
721  }
722  }
723  }
724  }
725 
726  }
727 }
728 
729 
730 DEAL_II_NAMESPACE_CLOSE
731 
732 #endif
size_type m() const
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
const unsigned int dofs_per_cell
Definition: fe_values.h:1459
#define AssertVectorVectorDimension(vec, dim1, dim2)
Definition: exceptions.h:1157
const FiniteElement< dim, spacedim > & get_fe() const
void upwind_value_residual(Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< double > &input, const std::vector< double > &data, const VectorSlice< const std::vector< std::vector< double > > > &velocity, double factor=1.)
Definition: advection.h:380
void cell_residual(Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< Tensor< 1, dim > > &input, const VectorSlice< const std::vector< std::vector< double > > > &velocity, double factor=1.)
Definition: advection.h:132
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 upwind_face_residual(Vector< double > &result1, Vector< double > &result2, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const std::vector< double > &input1, const std::vector< double > &input2, const VectorSlice< const std::vector< std::vector< double > > > &velocity, const double factor=1.)
Definition: advection.h:593
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 upwind_value_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const VectorSlice< const std::vector< std::vector< double > > > &velocity, double factor=1.)
Definition: advection.h:308
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const VectorSlice< const std::vector< std::vector< double > > > &velocity, const double factor=1.)
Definition: advection.h:75
double JxW(const unsigned int quadrature_point) const
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