Amandus: Simulations based on multilevel Schwarz methods
explicit.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 __brusselator_explicit_h
10 #define __brusselator_explicit_h
11 
12 #include <amandus/brusselator/implicit.h>
13 #include <amandus/brusselator/parameters.h>
14 #include <amandus/integrator.h>
15 #include <deal.II/base/polynomial.h>
16 #include <deal.II/base/tensor.h>
17 #include <deal.II/integrators/divergence.h>
18 #include <deal.II/integrators/l2.h>
19 #include <deal.II/integrators/laplace.h>
20 
21 using namespace dealii;
22 using namespace LocalIntegrators;
23 using namespace MeshWorker;
24 
25 namespace Brusselator
26 {
34 template <int dim>
36 {
37 public:
38  ExplicitResidual(const Parameters& par);
39 
40  virtual void cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
41  virtual void boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const;
42  virtual void face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2, IntegrationInfo<dim>& info1,
43  IntegrationInfo<dim>& info2) const;
44 
45 private:
46  SmartPointer<const Parameters, class ImplicitResidual<dim>> parameters;
47 };
48 
49 //----------------------------------------------------------------------//
50 
51 template <int dim>
53  : parameters(&par)
54 {
55  this->use_boundary = false;
56  this->use_face = true;
57 }
58 
59 template <int dim>
60 void
61 ExplicitResidual<dim>::cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
62 {
63  Assert(this->timestep != 0, ExcMessage("Only for transient problems"));
64  Assert(info.values.size() >= 1, ExcDimensionMismatch(info.values.size(), 1));
65  Assert(info.gradients.size() >= 1, ExcDimensionMismatch(info.values.size(), 1));
66 
67  std::vector<double> rhs0(info.fe_values(0).n_quadrature_points, 0.);
68  std::vector<double> rhs1(info.fe_values(0).n_quadrature_points, 0.);
69 
70  for (unsigned int k = 0; k < info.fe_values(0).n_quadrature_points; ++k)
71  {
72  const double u = info.values[0][0][k];
73  const double v = info.values[0][1][k];
74  rhs0[k] = u - this->timestep * (-parameters->B - u * u * v + (parameters->A + 1.) * u);
75  rhs1[k] = v - this->timestep * (-parameters->A * u + u * u * v);
76  }
77 
78  L2::L2(dinfo.vector(0).block(0), info.fe_values(0), rhs0);
79  L2::L2(dinfo.vector(0).block(1), info.fe_values(0), rhs1);
80  Laplace::cell_residual(dinfo.vector(0).block(0),
81  info.fe_values(0),
82  info.gradients[0][0],
83  -parameters->alpha0 * this->timestep);
84  Laplace::cell_residual(dinfo.vector(0).block(1),
85  info.fe_values(0),
86  info.gradients[0][1],
87  -parameters->alpha1 * this->timestep);
88 }
89 
90 template <int dim>
91 void
92 ExplicitResidual<dim>::boundary(DoFInfo<dim>&, IntegrationInfo<dim>&) const
93 {
94 }
95 
96 template <int dim>
97 void
98 ExplicitResidual<dim>::face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2, IntegrationInfo<dim>& info1,
99  IntegrationInfo<dim>& info2) const
100 {
101  const unsigned int deg = info1.fe_values(0).get_fe().tensor_degree();
102  Laplace::ip_residual(dinfo1.vector(0).block(0),
103  dinfo2.vector(0).block(0),
104  info1.fe_values(0),
105  info2.fe_values(0),
106  info1.values[0][0],
107  info1.gradients[0][0],
108  info2.values[0][0],
109  info2.gradients[0][0],
110  Laplace::compute_penalty(dinfo1, dinfo2, deg, deg),
111  -parameters->alpha0 * this->timestep);
112  Laplace::ip_residual(dinfo1.vector(0).block(1),
113  dinfo2.vector(0).block(1),
114  info1.fe_values(0),
115  info2.fe_values(0),
116  info1.values[0][1],
117  info1.gradients[0][1],
118  info2.values[0][1],
119  info2.gradients[0][1],
120  Laplace::compute_penalty(dinfo1, dinfo2, deg, deg),
121  -parameters->alpha1 * this->timestep);
122 }
123 }
124 
125 #endif
virtual void cell(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: explicit.h:61
Definition: explicit.h:25
Definition: explicit.h:35
double timestep
Current timestep if applicable.
Definition: integrator.h:48
Definition: brusselator/parameters.h:14
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 boundary(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: explicit.h:92
virtual void face(DoFInfo< dim > &dinfo1, DoFInfo< dim > &dinfo2, IntegrationInfo< dim > &info1, IntegrationInfo< dim > &info2) const
Definition: explicit.h:98
Definition: integrator.h:29
SmartPointer< const Parameters, class ImplicitResidual< dim > > parameters
Definition: explicit.h:46