9 #ifndef __brusselator_explicit_h 10 #define __brusselator_explicit_h 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> 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;
46 SmartPointer<const Parameters, class ImplicitResidual<dim>>
parameters;
55 this->use_boundary =
false;
56 this->use_face =
true;
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));
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.);
70 for (
unsigned int k = 0; k < info.fe_values(0).n_quadrature_points; ++k)
72 const double u = info.values[0][0][k];
73 const double v = info.values[0][1][k];
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);
99 IntegrationInfo<dim>& info2)
const 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),
107 info1.gradients[0][0],
109 info2.gradients[0][0],
110 Laplace::compute_penalty(dinfo1, dinfo2, deg, deg),
112 Laplace::ip_residual(dinfo1.vector(0).block(1),
113 dinfo2.vector(0).block(1),
117 info1.gradients[0][1],
119 info2.gradients[0][1],
120 Laplace::compute_penalty(dinfo1, dinfo2, deg, deg),
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