Reference documentation for deal.II version 8.5.1
grad_div.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_grad_div_h
17 #define dealii__integrators_grad_div_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 GradDiv
41  {
51  template <int dim>
52  void cell_matrix (
54  const FEValuesBase<dim> &fe,
55  double factor = 1.)
56  {
57  const unsigned int n_dofs = fe.dofs_per_cell;
58 
59  AssertDimension(fe.get_fe().n_components(), dim);
60  AssertDimension(M.m(), n_dofs);
61  AssertDimension(M.n(), n_dofs);
62 
63  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
64  {
65  const double dx = factor * fe.JxW(k);
66  for (unsigned int i=0; i<n_dofs; ++i)
67  for (unsigned int j=0; j<n_dofs; ++j)
68  {
69  const double divu = fe[FEValuesExtractors::Vector(0)].divergence(j,k);
70  const double divv = fe[FEValuesExtractors::Vector(0)].divergence(i,k);
71 
72  M(i,j) += dx * divu * divv;
73  }
74  }
75  }
76 
86  template <int dim, typename number>
88  Vector<number> &result,
89  const FEValuesBase<dim> &fetest,
90  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &input,
91  const double factor = 1.)
92  {
93  const unsigned int n_dofs = fetest.dofs_per_cell;
94 
95  AssertDimension(fetest.get_fe().n_components(), dim);
97 
98  for (unsigned int k=0; k<fetest.n_quadrature_points; ++k)
99  {
100  const double dx = factor * fetest.JxW(k);
101  for (unsigned int i=0; i<n_dofs; ++i)
102  {
103  const double divv = fetest[FEValuesExtractors::Vector(0)].divergence(i,k);
104  double du = 0.;
105  for (unsigned int d=0; d<dim; ++d)
106  du += input[d][k][d];
107 
108  result(i) += dx * du * divv;
109  }
110  }
111  }
112 
120  template <int dim>
121  inline void nitsche_matrix (
123  const FEValuesBase<dim> &fe,
124  double penalty,
125  double factor = 1.)
126  {
127  const unsigned int n_dofs = fe.dofs_per_cell;
128 
129  AssertDimension(fe.get_fe().n_components(), dim);
130  AssertDimension(M.m(), n_dofs);
131  AssertDimension(M.n(), n_dofs);
132 
133  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
134  {
135  const double dx = factor * fe.JxW(k);
136  const Tensor<1,dim> n = fe.normal_vector(k);
137  for (unsigned int i=0; i<n_dofs; ++i)
138  for (unsigned int j=0; j<n_dofs; ++j)
139  {
140  const double divu = fe[FEValuesExtractors::Vector(0)].divergence(j,k);
141  const double divv = fe[FEValuesExtractors::Vector(0)].divergence(i,k);
142  double un = 0., vn = 0.;
143  for (unsigned int d=0; d<dim; ++d)
144  {
145  un += fe.shape_value_component(j,k,d) * n[d];
146  vn += fe.shape_value_component(i,k,d) * n[d];
147  }
148 
149  M(i,j) += dx * 2. * penalty * un * vn;
150  M(i,j) -= dx*(divu*vn+divv*un);
151  }
152  }
153  }
154 
173  template <int dim>
175  Vector<double> &result,
176  const FEValuesBase<dim> &fe,
177  const VectorSlice<const std::vector<std::vector<double> > > &input,
178  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput,
179  const VectorSlice<const std::vector<std::vector<double> > > &data,
180  double penalty,
181  double factor = 1.)
182  {
183  const unsigned int n_dofs = fe.dofs_per_cell;
184  AssertDimension(fe.get_fe().n_components(), dim)
188 
189  for (unsigned int k=0; k<fe.n_quadrature_points; ++k)
190  {
191  const double dx = factor * fe.JxW(k);
192  const Tensor<1,dim> n = fe.normal_vector(k);
193 
194  double umgn = 0.;
195  double divu = 0.;
196  for (unsigned int d=0; d<dim; ++d)
197  {
198  umgn += (input[d][k] - data[d][k]) * n[d];
199  divu += Dinput[d][k][d];
200  }
201 
202  for (unsigned int i=0; i<n_dofs; ++i)
203  {
204  double vn = 0.;
205  const double divv = fe[FEValuesExtractors::Vector(0)].divergence(i,k);
206  for (unsigned int d=0; d<dim; ++d)
207  vn += fe.shape_value_component(i,k,d) * n[d];
208 
209  result(i) += dx*(2.*penalty*umgn*vn - divv*umgn - divu*vn);
210  }
211  }
212  }
213 
222  template <int dim>
223  void ip_matrix (
224  FullMatrix<double> &M11,
225  FullMatrix<double> &M12,
226  FullMatrix<double> &M21,
227  FullMatrix<double> &M22,
228  const FEValuesBase<dim> &fe1,
229  const FEValuesBase<dim> &fe2,
230  double penalty,
231  double factor1 = 1.,
232  double factor2 = -1.)
233  {
234  const unsigned int n_dofs = fe1.dofs_per_cell;
235  AssertDimension(M11.n(), n_dofs);
236  AssertDimension(M11.m(), n_dofs);
237  AssertDimension(M12.n(), n_dofs);
238  AssertDimension(M12.m(), n_dofs);
239  AssertDimension(M21.n(), n_dofs);
240  AssertDimension(M21.m(), n_dofs);
241  AssertDimension(M22.n(), n_dofs);
242  AssertDimension(M22.m(), n_dofs);
243 
244  const double fi = factor1;
245  const double fe = (factor2 < 0) ? factor1 : factor2;
246  const double f = .5*(fi+fe);
247 
248  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
249  {
250  const double dx = fe1.JxW(k);
251  const Tensor<1,dim> n = fe1.normal_vector(k);
252  for (unsigned int i=0; i<n_dofs; ++i)
253  for (unsigned int j=0; j<n_dofs; ++j)
254  {
255  double uni = 0.;
256  double une = 0.;
257  double vni = 0.;
258  double vne = 0.;
259  const double divui = fe1[FEValuesExtractors::Vector(0)].divergence(j,k);
260  const double divue = fe2[FEValuesExtractors::Vector(0)].divergence(j,k);
261  const double divvi = fe1[FEValuesExtractors::Vector(0)].divergence(i,k);
262  const double divve = fe2[FEValuesExtractors::Vector(0)].divergence(i,k);
263 
264  for (unsigned int d=0; d<dim; ++d)
265  {
266  uni += fe1.shape_value_component(j,k,d) * n[d];
267  une += fe2.shape_value_component(j,k,d) * n[d];
268  vni += fe1.shape_value_component(i,k,d) * n[d];
269  vne += fe2.shape_value_component(i,k,d) * n[d];
270  }
271  M11(i,j) += dx*(-.5*fi*divvi*uni-.5*fi*divui*vni+f*penalty*uni*vni);
272  M12(i,j) += dx*( .5*fi*divvi*une-.5*fe*divue*vni-f*penalty*vni*une);
273  M21(i,j) += dx*(-.5*fe*divve*uni+.5*fi*divui*vne-f*penalty*uni*vne);
274  M22(i,j) += dx*( .5*fe*divve*une+.5*fe*divue*vne+f*penalty*une*vne);
275  }
276  }
277  }
278 
293  template<int dim>
294  void
296  Vector<double> &result1,
297  Vector<double> &result2,
298  const FEValuesBase<dim> &fe1,
299  const FEValuesBase<dim> &fe2,
300  const VectorSlice<const std::vector<std::vector<double> > > &input1,
301  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput1,
302  const VectorSlice<const std::vector<std::vector<double> > > &input2,
303  const VectorSlice<const std::vector<std::vector<Tensor<1,dim> > > > &Dinput2,
304  double pen,
305  double int_factor = 1.,
306  double ext_factor = -1.)
307  {
308  const unsigned int n1 = fe1.dofs_per_cell;
309 
310  AssertDimension(fe1.get_fe().n_components(), dim);
315 
316  const double fi = int_factor;
317  const double fe = (ext_factor < 0) ? int_factor : ext_factor;
318  const double penalty = .5 * pen * (fi + fe);
319 
320 
321  for (unsigned int k=0; k<fe1.n_quadrature_points; ++k)
322  {
323  const double dx = fe1.JxW(k);
324  const Tensor<1,dim> n = fe1.normal_vector(k);
325  double uni = 0.;
326  double une = 0.;
327  double divui = 0.;
328  double divue = 0.;
329  for (unsigned int d=0; d<dim; ++d)
330  {
331  uni += input1[d][k]*n[d];
332  une += input2[d][k]*n[d];
333  divui += Dinput1[d][k][d];
334  divue += Dinput2[d][k][d];
335  }
336 
337  for (unsigned int i=0; i<n1; ++i)
338  {
339  double vni = 0.;
340  double vne = 0.;
341  const double divvi = fe1[FEValuesExtractors::Vector(0)].divergence(i,k);
342  const double divve = fe2[FEValuesExtractors::Vector(0)].divergence(i,k);
343  for (unsigned int d=0; d<dim; ++d)
344  {
345  vni += fe1.shape_value_component(i,k,d)*n[d];
346  vne += fe2.shape_value_component(i,k,d)*n[d];
347  }
348 
349  result1(i) += dx*(-.5*fi*divvi*uni-.5*fi*divui*vni+penalty*uni*vni);
350  result1(i) += dx*( .5*fi*divvi*une-.5*fe*divue*vni-penalty*vni*une);
351  result2(i) += dx*(-.5*fe*divve*uni+.5*fi*divui*vne-penalty*uni*vne);
352  result2(i) += dx*( .5*fe*divve*une+.5*fe*divue*vne+penalty*une*vne);
353  }
354  }
355  }
356  }
357 }
358 
359 DEAL_II_NAMESPACE_CLOSE
360 
361 
362 #endif
size_type m() const
void nitsche_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
Definition: grad_div.h:121
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
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:295
const unsigned int dofs_per_cell
Definition: fe_values.h:1459
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double factor=1.)
Definition: grad_div.h:52
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:223
#define AssertVectorVectorDimension(vec, dim1, dim2)
Definition: exceptions.h:1157
const FiniteElement< dim, spacedim > & get_fe() const
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
size_type n() const
const unsigned int n_quadrature_points
Definition: fe_values.h:1452
double JxW(const unsigned int quadrature_point) 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: grad_div.h:87
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:174
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const