Amandus: Simulations based on multilevel Schwarz methods
allen_cahn/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 __allen_cahn_polynomial_h
10 #define __allen_cahn_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 namespace AllenCahn
20 {
21 using namespace dealii;
22 using namespace LocalIntegrators;
23 using namespace MeshWorker;
24 
32 template <int dim>
34 {
35 public:
36  PolynomialResidual(double diffusion, const Polynomials::Polynomial<double> solution_1d);
37 
38  virtual void cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
39  virtual void boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
40  virtual void face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2, IntegrationInfo<dim>& info1,
41  IntegrationInfo<dim>& info2) const;
42 
43 private:
44  double D;
45  Polynomials::Polynomial<double> solution_1d;
46 };
47 
48 template <int dim>
50 {
51 public:
52  PolynomialError(const Polynomials::Polynomial<double> solution_1d);
53 
54  virtual void cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
55  virtual void boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
56  virtual void face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2, IntegrationInfo<dim>& info1,
57  IntegrationInfo<dim>& info2) const;
58 
59 private:
60  Polynomials::Polynomial<double> solution_1d;
61 };
62 
63 //----------------------------------------------------------------------//
64 
65 template <int dim>
67  const Polynomials::Polynomial<double> solution_1d)
68  : D(diffusion)
69  , solution_1d(solution_1d)
70 {
71  this->use_boundary = false;
72  this->use_face = true;
73  this->input_vector_names.push_back("Newton iterate");
74 }
75 
76 template <int dim>
77 void
78 PolynomialResidual<dim>::cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
79 {
80  Assert(info.values.size() >= 1, ExcDimensionMismatch(info.values.size(), 1));
81  Assert(info.gradients.size() >= 1, ExcDimensionMismatch(info.values.size(), 1));
82 
83  std::vector<double> rhs(info.fe_values(0).n_quadrature_points, 0.);
84 
85  std::vector<double> px(3);
86  std::vector<double> py(3);
87  for (unsigned int k = 0; k < info.fe_values(0).n_quadrature_points; ++k)
88  {
89  const double x = info.fe_values(0).quadrature_point(k)(0);
90  const double y = info.fe_values(0).quadrature_point(k)(1);
91  solution_1d.value(x, px);
92  solution_1d.value(y, py);
93 
94  // negative Laplacian
95  rhs[k] = D * px[2] * py[0] + D * px[0] * py[2];
96  // nonlinearity of true solution
97  rhs[k] -= px[0] * py[0] * (px[0] * py[0] * px[0] * py[0] - 1.);
98 
99  const double u = info.values[0][0][k];
100  rhs[k] += u * (u * u - 1.);
101  }
102 
103  L2::L2(dinfo.vector(0).block(0), info.fe_values(0), rhs);
104  Laplace::cell_residual(dinfo.vector(0).block(0), info.fe_values(0), info.gradients[0][0], D);
105 }
106 
107 template <int dim>
108 void
109 PolynomialResidual<dim>::boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
110 {
111  std::vector<double> null(info.fe_values(0).n_quadrature_points, 0.);
112 
113  const unsigned int deg = info.fe_values(0).get_fe().tensor_degree();
114  Laplace::nitsche_residual(dinfo.vector(0).block(0),
115  info.fe_values(0),
116  info.values[0][0],
117  info.gradients[0][0],
118  null,
119  Laplace::compute_penalty(dinfo, dinfo, deg, deg),
120  D);
121 }
122 
123 template <int dim>
124 void
125 PolynomialResidual<dim>::face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2,
126  IntegrationInfo<dim>& info1, IntegrationInfo<dim>& info2) const
127 {
128  const unsigned int deg = info1.fe_values(0).get_fe().tensor_degree();
129  Laplace::ip_residual(dinfo1.vector(0).block(0),
130  dinfo2.vector(0).block(0),
131  info1.fe_values(0),
132  info2.fe_values(0),
133  info1.values[0][0],
134  info1.gradients[0][0],
135  info2.values[0][0],
136  info2.gradients[0][0],
137  Laplace::compute_penalty(dinfo1, dinfo2, deg, deg),
138  D);
139 }
140 
141 //----------------------------------------------------------------------//
142 
143 template <int dim>
144 PolynomialError<dim>::PolynomialError(const Polynomials::Polynomial<double> solution_1d)
145  : solution_1d(solution_1d)
146 {
147  this->use_boundary = false;
148  this->use_face = false;
149 }
150 
151 template <int dim>
152 void
153 PolynomialError<dim>::cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
154 {
155  // Assert(dinfo.n_values() >= 4, ExcDimensionMismatch(dinfo.n_values(), 4));
156 
157  std::vector<double> px(2);
158  std::vector<double> py(2);
159  for (unsigned int k = 0; k < info.fe_values(0).n_quadrature_points; ++k)
160  {
161  const double x = info.fe_values(0).quadrature_point(k)(0);
162  const double y = info.fe_values(0).quadrature_point(k)(1);
163  solution_1d.value(x, px);
164  solution_1d.value(y, py);
165  const double dx = info.fe_values(0).JxW(k);
166 
167  Tensor<1, dim> Du = info.gradients[0][0][k];
168  Du[0] -= px[1] * py[0];
169  Du[1] -= px[0] * py[1];
170  double u = info.values[0][0][k];
171  u -= px[0] * py[0];
172  unsigned int i = 0;
173  // 0. L^2(u)
174  dinfo.value(i++) += u * u * dx;
175  // 1. H^1(u)
176  dinfo.value(i++) += (Du * Du) * dx;
177  }
178 }
179 
180 template <int dim>
181 void
182 PolynomialError<dim>::boundary(DoFInfo<dim>&, IntegrationInfo<dim>&) const
183 {
184 }
185 
186 template <int dim>
187 void
188 PolynomialError<dim>::face(DoFInfo<dim>&, DoFInfo<dim>&, IntegrationInfo<dim>&,
189  IntegrationInfo<dim>&) const
190 {
191 }
192 }
193 
194 #endif
virtual void face(DoFInfo< dim > &dinfo1, DoFInfo< dim > &dinfo2, IntegrationInfo< dim > &info1, IntegrationInfo< dim > &info2) const
Definition: allen_cahn/polynomial.h:188
virtual void boundary(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: allen_cahn/polynomial.h:182
Polynomials::Polynomial< double > solution_1d
Definition: allen_cahn/polynomial.h:45
virtual void cell(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: allen_cahn/polynomial.h:78
virtual void face(DoFInfo< dim > &dinfo1, DoFInfo< dim > &dinfo2, IntegrationInfo< dim > &info1, IntegrationInfo< dim > &info2) const
Definition: allen_cahn/polynomial.h:125
PolynomialResidual(double diffusion, const Polynomials::Polynomial< double > solution_1d)
Definition: allen_cahn/polynomial.h:66
Definition: allen_cahn/polynomial.h:49
Polynomials::Polynomial< double > solution_1d
Definition: allen_cahn/polynomial.h:60
virtual void boundary(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: allen_cahn/polynomial.h:109
double D
Definition: allen_cahn/polynomial.h:44
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
Definition: allen_cahn/polynomial.h:33
PolynomialError(const Polynomials::Polynomial< double > solution_1d)
Definition: allen_cahn/polynomial.h:144
virtual void cell(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: allen_cahn/polynomial.h:153
Definition: allen_cahn/matrix.h:19
Definition: integrator.h:29