Amandus: Simulations based on multilevel Schwarz methods
allen_cahn/residual.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_residual_h
10 #define __allen_cahn_residual_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>
33 class Residual : public AmandusIntegrator<dim>
34 {
35 public:
36  Residual(double diffusion);
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 };
46 
47 //----------------------------------------------------------------------//
48 
49 template <int dim>
50 Residual<dim>::Residual(double diffusion)
51  : D(diffusion)
52 {
53  this->use_boundary = false;
54  this->use_face = true;
55 }
56 
57 template <int dim>
58 void
59 Residual<dim>::cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
60 {
61  Assert(info.values.size() >= 1, ExcDimensionMismatch(info.values.size(), 1));
62  Assert(info.gradients.size() >= 1, ExcDimensionMismatch(info.values.size(), 1));
63 
64  std::vector<double> rhs(info.fe_values(0).n_quadrature_points, 0.);
65 
66  for (unsigned int k = 0; k < info.fe_values(0).n_quadrature_points; ++k)
67  {
68  const double u = info.values[0][0][k];
69  rhs[k] += u * (u * u - 1.);
70  }
71 
72  L2::L2(dinfo.vector(0).block(0), info.fe_values(0), rhs);
73  Laplace::cell_residual(dinfo.vector(0).block(0), info.fe_values(0), info.gradients[0][0], D);
74 }
75 
76 template <int dim>
77 void
78 Residual<dim>::boundary(DoFInfo<dim>&, IntegrationInfo<dim>&) const
79 {
80 }
81 
82 template <int dim>
83 void
84 Residual<dim>::face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2, IntegrationInfo<dim>& info1,
85  IntegrationInfo<dim>& info2) const
86 {
87  const unsigned int deg = info1.fe_values(0).get_fe().tensor_degree();
88  Laplace::ip_residual(dinfo1.vector(0).block(0),
89  dinfo2.vector(0).block(0),
90  info1.fe_values(0),
91  info2.fe_values(0),
92  info1.values[0][0],
93  info1.gradients[0][0],
94  info2.values[0][0],
95  info2.gradients[0][0],
96  Laplace::compute_penalty(dinfo1, dinfo2, deg, deg),
97  D);
98 }
99 }
100 
101 #endif
virtual void cell(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: allen_cahn/residual.h:59
Residual(double diffusion)
Definition: allen_cahn/residual.h:50
Definition: allen_cahn/residual.h:33
double D
Definition: allen_cahn/residual.h:44
virtual void boundary(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: allen_cahn/residual.h:78
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
virtual void face(DoFInfo< dim > &dinfo1, DoFInfo< dim > &dinfo2, IntegrationInfo< dim > &info1, IntegrationInfo< dim > &info2) const
Definition: allen_cahn/residual.h:84
Definition: allen_cahn/matrix.h:19
Definition: integrator.h:29