Amandus: Simulations based on multilevel Schwarz methods
advection/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 #ifndef __advection_polynomial_h
8 #define __advection_polynomial_h
9 
10 #include <advection/parameters.h>
11 #include <amandus/integrator.h>
12 #include <deal.II/base/polynomial.h>
13 #include <deal.II/base/tensor.h>
14 #include <deal.II/integrators/advection.h>
15 #include <deal.II/integrators/l2.h>
16 
17 using namespace dealii;
18 using namespace LocalIntegrators;
19 using namespace MeshWorker;
20 
21 namespace Advection
22 {
35 template <int dim>
36 class PolynomialRHS : public AmandusIntegrator<dim>
37 {
38 public:
39  PolynomialRHS(const Parameters& par,
40  const std::vector<Polynomials::Polynomial<double>> potentials_1d);
41 
42  virtual void cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
43  virtual void boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
44 
45 private:
46  dealii::SmartPointer<const Parameters, class PolynomialRHS<dim>> parameters;
47  std::vector<Polynomials::Polynomial<double>> potentials_1d;
48 
49  std::vector<std::vector<double>> velocity;
50 };
51 
77 template <int dim>
79 {
80 public:
81  PolynomialError(const Parameters& par,
82  const std::vector<Polynomials::Polynomial<double>> potentials_1d);
83 
84  virtual void cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
85  virtual void boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
86  virtual void face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2, IntegrationInfo<dim>& info1,
87  IntegrationInfo<dim>& info2) const;
88 
89 private:
90  dealii::SmartPointer<const Parameters, class PolynomialError<dim>> parameters;
91  std::vector<Polynomials::Polynomial<double>> potentials_1d;
92 };
93 
94 //----------------------------------------------------------------------//
95 
96 template <int dim>
97 PolynomialRHS<dim>::PolynomialRHS(const Parameters& par,
98  const std::vector<Polynomials::Polynomial<double>> potentials_1d)
99  : parameters(&par)
100  , potentials_1d(potentials_1d)
101  , velocity(dim, std::vector<double>(1))
102 {
103  this->use_boundary = true;
104  this->use_face = false;
105  velocity[0][0] = 1.;
106  velocity[1][0] = 2.;
107 }
108 
109 template <int dim>
110 void
111 PolynomialRHS<dim>::cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
112 {
113  std::vector<double> rhs(info.fe_values(0).n_quadrature_points, 0.);
114 
115  std::vector<double> px(2);
116  std::vector<double> py(2);
117  for (unsigned int k = 0; k < info.fe_values(0).n_quadrature_points; ++k)
118  {
119  const double x = info.fe_values(0).quadrature_point(k)(0);
120  const double y = info.fe_values(0).quadrature_point(k)(1);
121  potentials_1d[0].value(x, px);
122  potentials_1d[0].value(y, py);
123 
124  rhs[k] = velocity[0][0] * px[1] * py[0] + velocity[1][0] * px[0] * py[1];
125  }
126 
127  L2::L2(dinfo.vector(0).block(0), info.fe_values(0), rhs);
128 }
129 
130 template <int dim>
131 void
132 PolynomialRHS<dim>::boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
133 {
134  std::vector<double> rhs(info.fe_values(0).n_quadrature_points, 0.);
135  std::vector<double> null(info.fe_values(0).n_quadrature_points, 0.);
136 
137  std::vector<double> px(2);
138  std::vector<double> py(2);
139  for (unsigned int k = 0; k < info.fe_values(0).n_quadrature_points; ++k)
140  {
141  const double x = info.fe_values(0).quadrature_point(k)(0);
142  const double y = info.fe_values(0).quadrature_point(k)(1);
143  potentials_1d[0].value(x, px);
144  potentials_1d[0].value(y, py);
145 
146  rhs[k] = px[0] * py[0];
147  }
148 
149  dealii::LocalIntegrators::Advection::upwind_value_residual(
150  dinfo.vector(0).block(0), info.fe_values(0), null, rhs, velocity);
151 }
152 
153 //----------------------------------------------------------------------//
154 
155 template <int dim>
157  const Parameters& par, const std::vector<Polynomials::Polynomial<double>> potentials_1d)
158  : parameters(&par)
160 {
161  this->use_boundary = false;
162  this->use_face = false;
163  this->error_types.push_back(2);
164  this->error_types.push_back(2);
165 }
166 
167 template <int dim>
168 void
169 PolynomialError<dim>::cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
170 {
171  Assert(dinfo.n_values() >= 2, ExcDimensionMismatch(dinfo.n_values(), 4));
172 
173  std::vector<double> px(3);
174  std::vector<double> py(3);
175  for (unsigned int k = 0; k < info.fe_values(0).n_quadrature_points; ++k)
176  {
177  const double x = info.fe_values(0).quadrature_point(k)(0);
178  const double y = info.fe_values(0).quadrature_point(k)(1);
179  potentials_1d[0].value(x, px);
180  potentials_1d[0].value(y, py);
181  const double dx = info.fe_values(0).JxW(k);
182 
183  Tensor<1, dim> Du = info.gradients[0][0][k];
184  Du[0] -= px[1] * py[0];
185  Du[1] -= px[0] * py[1];
186  double u = info.values[0][0][k];
187  u -= px[0] * py[0];
188 
189  // 0. L^2(u)
190  dinfo.value(0) += (u * u) * dx;
191  // 1. H^1(u)
192  dinfo.value(1) += (Du * Du) * dx;
193  }
194 
195  for (unsigned int i = 0; i < 2; ++i)
196  dinfo.value(i) = std::sqrt(dinfo.value(i));
197 }
198 
199 template <int dim>
200 void
201 PolynomialError<dim>::boundary(DoFInfo<dim>&, IntegrationInfo<dim>&) const
202 {
203 }
204 
205 template <int dim>
206 void
207 PolynomialError<dim>::face(DoFInfo<dim>&, DoFInfo<dim>&, IntegrationInfo<dim>&,
208  IntegrationInfo<dim>&) const
209 {
210 }
211 }
212 
213 #endif
virtual void boundary(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: advection/polynomial.h:132
STL namespace.
std::vector< Polynomials::Polynomial< double > > potentials_1d
Definition: advection/polynomial.h:91
Definition: advection/polynomial.h:78
virtual void face(DoFInfo< dim > &dinfo1, DoFInfo< dim > &dinfo2, IntegrationInfo< dim > &info1, IntegrationInfo< dim > &info2) const
Definition: advection/polynomial.h:207
Definition: advection/matrix.h:17
dealii::SmartPointer< const Parameters, class PolynomialError< dim > > parameters
Definition: advection/polynomial.h:90
virtual void cell(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: advection/polynomial.h:111
Definition: advection/polynomial.h:36
dealii::SmartPointer< const Parameters, class PolynomialRHS< dim > > parameters
Definition: advection/polynomial.h:46
std::vector< unsigned int > error_types
Definition: integrator.h:113
PolynomialError(const Parameters &par, const std::vector< Polynomials::Polynomial< double >> potentials_1d)
Definition: advection/polynomial.h:156
std::vector< std::vector< double > > velocity
Definition: advection/polynomial.h:49
std::vector< Polynomials::Polynomial< double > > potentials_1d
Definition: advection/polynomial.h:47
virtual void cell(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: advection/polynomial.h:169
virtual void boundary(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: advection/polynomial.h:201
Definition: integrator.h:29