Amandus: Simulations based on multilevel Schwarz methods
stokes/polynomial.h
Go to the documentation of this file.
1 /**********************************************************************
2  * Copyright (C) 2011 - 2014 by the authors
3  * Distributed under the MIT License
4  *
5  * See the files AUTHORS and LICENSE in the project root directory
6  *
7  **********************************************************************/
8 
9 #ifndef __stokes_polynomial_h
10 #define __stokes_polynomial_h
11 
12 #include <amandus/integrator.h>
13 #include <deal.II/base/polynomial.h>
14 #include <deal.II/base/tensor.h>
15 #include <deal.II/integrators/divergence.h>
16 #include <deal.II/integrators/l2.h>
17 #include <deal.II/integrators/laplace.h>
18 
19 using namespace dealii;
20 using namespace LocalIntegrators;
21 using namespace MeshWorker;
22 
23 namespace StokesIntegrators
24 {
39 template <int dim>
40 class PolynomialRHS : public AmandusIntegrator<dim>
41 {
42 public:
43  PolynomialRHS(const Polynomials::Polynomial<double> curl_potential_1d,
44  const Polynomials::Polynomial<double> grad_potential_1d);
45 
46  virtual void cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
47  virtual void boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
48  virtual void face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2, IntegrationInfo<dim>& info1,
49  IntegrationInfo<dim>& info2) const;
50 
51 private:
52  Polynomials::Polynomial<double> curl_potential_1d;
53  Polynomials::Polynomial<double> grad_potential_1d;
54 };
55 
75 template <int dim>
77 {
78 public:
79  PolynomialResidual(const Polynomials::Polynomial<double> curl_potential_1d,
80  const Polynomials::Polynomial<double> grad_potential_1d);
81 
82  virtual void cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
83  virtual void boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
84  virtual void face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2, IntegrationInfo<dim>& info1,
85  IntegrationInfo<dim>& info2) const;
86 
87 private:
88  Polynomials::Polynomial<double> curl_potential_1d;
89  Polynomials::Polynomial<double> grad_potential_1d;
90 };
91 
125 template <int dim>
127 {
128 public:
129  PolynomialError(const Polynomials::Polynomial<double> curl_potential_1d,
130  const Polynomials::Polynomial<double> grad_potential_1d);
131 
132  virtual void cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
133  virtual void boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
134  virtual void face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2, IntegrationInfo<dim>& info1,
135  IntegrationInfo<dim>& info2) const;
136 
137 private:
138  Polynomials::Polynomial<double> curl_potential_1d;
139  Polynomials::Polynomial<double> grad_potential_1d;
140 };
141 
142 //----------------------------------------------------------------------//
143 
144 template <int dim>
145 PolynomialRHS<dim>::PolynomialRHS(const Polynomials::Polynomial<double> curl_potential_1d,
146  const Polynomials::Polynomial<double> grad_potential_1d)
147  : curl_potential_1d(curl_potential_1d)
148  , grad_potential_1d(grad_potential_1d)
149 {
150  this->use_boundary = false;
151  this->use_face = false;
152 }
153 
154 template <int dim>
155 void
156 PolynomialRHS<dim>::cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
157 {
158  std::vector<std::vector<double>> rhs(dim,
159  std::vector<double>(info.fe_values(0).n_quadrature_points));
160 
161  std::vector<double> px(4);
162  std::vector<double> py(4);
163  for (unsigned int k = 0; k < info.fe_values(0).n_quadrature_points; ++k)
164  {
165  const double x = info.fe_values(0).quadrature_point(k)(0);
166  const double y = info.fe_values(0).quadrature_point(k)(1);
167  curl_potential_1d.value(x, px);
168  curl_potential_1d.value(y, py);
169 
170  rhs[0][k] = -px[2] * py[1] - px[0] * py[3];
171  rhs[1][k] = px[3] * py[0] + px[1] * py[2];
172 
173  // Add a gradient part to the
174  // right hand side to test for
175  // pressure
176  grad_potential_1d.value(x, px);
177  grad_potential_1d.value(y, py);
178  rhs[0][k] += px[1] * py[0];
179  rhs[1][k] += px[0] * py[1];
180  }
181 
182  L2::L2(dinfo.vector(0).block(0), info.fe_values(0), rhs);
183 }
184 
185 template <int dim>
186 void
187 PolynomialRHS<dim>::boundary(DoFInfo<dim>&, IntegrationInfo<dim>&) const
188 {
189 }
190 
191 template <int dim>
192 void
193 PolynomialRHS<dim>::face(DoFInfo<dim>&, DoFInfo<dim>&, IntegrationInfo<dim>&,
194  IntegrationInfo<dim>&) const
195 {
196 }
197 
198 //----------------------------------------------------------------------//
199 
200 template <int dim>
202  const Polynomials::Polynomial<double> grad_potential_1d)
203  : curl_potential_1d(curl_potential_1d)
204  , grad_potential_1d(grad_potential_1d)
205 {
206  this->use_boundary = true;
207  this->use_face = true;
208  this->input_vector_names.push_back("Newton iterate");
209 }
210 
211 template <int dim>
212 void
213 PolynomialResidual<dim>::cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
214 {
215  Assert(info.values.size() >= 1, ExcDimensionMismatch(info.values.size(), 1));
216  Assert(info.gradients.size() >= 1, ExcDimensionMismatch(info.values.size(), 1));
217 
218  std::vector<std::vector<double>> rhs(dim,
219  std::vector<double>(info.fe_values(0).n_quadrature_points));
220 
221  std::vector<double> px(4);
222  std::vector<double> py(4);
223  for (unsigned int k = 0; k < info.fe_values(0).n_quadrature_points; ++k)
224  {
225  const double x = info.fe_values(0).quadrature_point(k)(0);
226  const double y = info.fe_values(0).quadrature_point(k)(1);
227  curl_potential_1d.value(x, px);
228  curl_potential_1d.value(y, py);
229 
230  rhs[0][k] = -px[2] * py[1] - px[0] * py[3];
231  rhs[1][k] = px[3] * py[0] + px[1] * py[2];
232 
233  // Add a gradient part to the
234  // right hand side to test for
235  // pressure
236  grad_potential_1d.value(x, px);
237  grad_potential_1d.value(y, py);
238  rhs[0][k] += px[1] * py[0];
239  rhs[1][k] += px[0] * py[1];
240  }
241 
242  L2::L2(dinfo.vector(0).block(0), info.fe_values(0), rhs, -1.);
244  dinfo.vector(0).block(0), info.fe_values(0), make_slice(info.gradients[0], 0, dim));
245  Divergence::gradient_residual(
246  dinfo.vector(0).block(0), info.fe_values(0), info.values[0][dim], -1.);
248  dinfo.vector(0).block(1), info.fe_values(1), make_slice(info.gradients[0], 0, dim), 1.);
249 }
250 
251 template <int dim>
252 void
253 PolynomialResidual<dim>::boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
254 {
255  std::vector<std::vector<double>> null(
256  dim, std::vector<double>(info.fe_values(0).n_quadrature_points, 0.));
257 
258  const unsigned int deg = info.fe_values(0).get_fe().tensor_degree();
259  Laplace::nitsche_residual(dinfo.vector(0).block(0),
260  info.fe_values(0),
261  make_slice(info.values[0], 0, dim),
262  make_slice(info.gradients[0], 0, dim),
263  null,
264  Laplace::compute_penalty(dinfo, dinfo, deg, deg));
265 }
266 
267 template <int dim>
268 void
269 PolynomialResidual<dim>::face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2,
270  IntegrationInfo<dim>& info1, IntegrationInfo<dim>& info2) const
271 {
272  const unsigned int deg = info1.fe_values(0).get_fe().tensor_degree();
273  Laplace::ip_residual(dinfo1.vector(0).block(0),
274  dinfo2.vector(0).block(0),
275  info1.fe_values(0),
276  info2.fe_values(0),
277  make_slice(info1.values[0], 0, dim),
278  make_slice(info1.gradients[0], 0, dim),
279  make_slice(info2.values[0], 0, dim),
280  make_slice(info2.gradients[0], 0, dim),
281  Laplace::compute_penalty(dinfo1, dinfo2, deg, deg));
282 }
283 
284 //----------------------------------------------------------------------//
285 
286 template <int dim>
287 PolynomialError<dim>::PolynomialError(const Polynomials::Polynomial<double> curl_potential_1d,
288  const Polynomials::Polynomial<double> grad_potential_1d)
289  : curl_potential_1d(curl_potential_1d)
290  , grad_potential_1d(grad_potential_1d)
291 {
292  this->error_types.push_back(2);
293  this->error_names.push_back("L2u");
294  this->error_types.push_back(2);
295  this->error_names.push_back("H1u");
296  this->error_types.push_back(2);
297  this->error_names.push_back("L2divu");
298  this->error_types.push_back(2);
299  this->error_names.push_back("L2p");
300  this->error_types.push_back(0);
301  this->error_names.push_back("Linftyp");
302  this->error_types.push_back(2);
303  this->error_names.push_back("H1p");
304  this->error_types.push_back(1);
305  this->error_names.push_back("meanp");
306  this->use_boundary = false;
307  this->use_face = false;
308 }
309 
310 template <int dim>
311 void
312 PolynomialError<dim>::cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
313 {
314  // Assert(dinfo.n_values() >= 4, ExcDimensionMismatch(dinfo.n_values(), 4));
315 
316  std::vector<double> px(3);
317  std::vector<double> py(3);
318  for (unsigned int k = 0; k < info.fe_values(0).n_quadrature_points; ++k)
319  {
320  const double x = info.fe_values(0).quadrature_point(k)(0);
321  const double y = info.fe_values(0).quadrature_point(k)(1);
322  curl_potential_1d.value(x, px);
323  curl_potential_1d.value(y, py);
324  const double dx = info.fe_values(0).JxW(k);
325 
326  Tensor<1, dim> Du0 = info.gradients[0][0][k];
327  Du0[0] -= px[1] * py[1];
328  Du0[1] -= px[0] * py[2];
329  Tensor<1, dim> Du1 = info.gradients[0][1][k];
330  Du1[0] += px[2] * py[0];
331  Du1[1] += px[1] * py[1];
332  double u0 = info.values[0][0][k];
333  u0 -= px[0] * py[1];
334  double u1 = info.values[0][1][k];
335  u1 += px[1] * py[0];
336 
337  grad_potential_1d.value(x, px);
338  grad_potential_1d.value(y, py);
339  double p = info.values[0][dim][k];
340  p += px[0] * py[0];
341  Tensor<1, dim> Dp = info.gradients[0][dim][k];
342  Dp[0] += px[1] * py[0];
343  Dp[1] += px[0] * py[1];
344  // 0. L^2(u)
345  dinfo.value(0) += (u0 * u0 + u1 * u1) * dx;
346  // 1. H^1(u)
347  dinfo.value(1) += ((Du0 * Du0) + (Du1 * Du1)) * dx;
348  // 2. div u
349  dinfo.value(2) = Divergence::norm(info.fe_values(0), make_slice(info.gradients[0], 0, dim));
350  // 3. L^2(p) up to mean value
351  dinfo.value(3) += p * p * dx;
352  // 3. L^\infty(p) up to mean value
353  dinfo.value(4) = std::max(dinfo.value(4), std::abs(p));
354  // 4. H^1(p)
355  dinfo.value(5) += Dp * Dp * dx;
356  // 5. Mean value of p
357  dinfo.value(6) += p * dx;
358  }
359 
360  for (unsigned int i = 0; i <= 4; ++i)
361  if (this->error_type(i) == 2)
362  dinfo.value(i) = std::sqrt(dinfo.value(i));
363 }
364 
365 template <int dim>
366 void
367 PolynomialError<dim>::boundary(DoFInfo<dim>&, IntegrationInfo<dim>&) const
368 {
369 }
370 
371 template <int dim>
372 void
373 PolynomialError<dim>::face(DoFInfo<dim>&, DoFInfo<dim>&, IntegrationInfo<dim>&,
374  IntegrationInfo<dim>&) const
375 {
376 }
377 }
378 
379 #endif
virtual void boundary(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: stokes/polynomial.h:187
virtual void boundary(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: stokes/polynomial.h:367
std::vector< std::string > error_names
Definition: integrator.h:120
PolynomialResidual(const Polynomials::Polynomial< double > curl_potential_1d, const Polynomials::Polynomial< double > grad_potential_1d)
Definition: stokes/polynomial.h:201
Polynomials::Polynomial< double > grad_potential_1d
Definition: stokes/polynomial.h:89
Polynomials::Polynomial< double > grad_potential_1d
Definition: stokes/polynomial.h:139
Polynomials::Polynomial< double > curl_potential_1d
Definition: stokes/polynomial.h:88
Definition: stokes/polynomial.h:126
virtual void boundary(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: stokes/polynomial.h:253
virtual void face(DoFInfo< dim > &dinfo1, DoFInfo< dim > &dinfo2, IntegrationInfo< dim > &info1, IntegrationInfo< dim > &info2) const
Definition: stokes/polynomial.h:193
Definition: stokes/polynomial.h:40
virtual void cell(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: stokes/polynomial.h:213
Definition: stokes/eigen.h:37
virtual void cell(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: stokes/polynomial.h:312
virtual void face(DoFInfo< dim > &dinfo1, DoFInfo< dim > &dinfo2, IntegrationInfo< dim > &info1, IntegrationInfo< dim > &info2) const
Definition: stokes/polynomial.h:269
virtual void cell(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: stokes/polynomial.h:156
Polynomials::Polynomial< double > grad_potential_1d
Definition: stokes/polynomial.h:53
std::vector< unsigned int > error_types
Definition: integrator.h:113
PolynomialError(const Polynomials::Polynomial< double > curl_potential_1d, const Polynomials::Polynomial< double > grad_potential_1d)
Definition: stokes/polynomial.h:287
void cell_residual(dealii::Vector< number > &result, const dealii::FEValuesBase< dim > &fe, const dealii::VectorSlice< const std::vector< std::vector< dealii::Tensor< 1, dim >>>> &input, double lambda=0., double mu=1.)
Definition: elasticity/integrators.h:23
Polynomials::Polynomial< double > curl_potential_1d
Definition: stokes/polynomial.h:138
Definition: stokes/polynomial.h:76
Polynomials::Polynomial< double > curl_potential_1d
Definition: stokes/polynomial.h:52
virtual void face(DoFInfo< dim > &dinfo1, DoFInfo< dim > &dinfo2, IntegrationInfo< dim > &info1, IntegrationInfo< dim > &info2) const
Definition: stokes/polynomial.h:373
unsigned int error_type(unsigned int i) const
Definition: integrator.h:437
Definition: integrator.h:29