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