7 #ifndef __advection_matrix_h 8 #define __advection_matrix_h 10 #include <amandus/integrator.h> 11 #include <deal.II/integrators/advection.h> 12 #include <deal.II/integrators/l2.h> 13 #include <deal.II/meshworker/integration_info.h> 34 Matrix(
const Parameters& par);
36 virtual void cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info)
const;
37 virtual void boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info)
const;
38 virtual void face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2, IntegrationInfo<dim>& info1,
39 IntegrationInfo<dim>& info2)
const;
42 dealii::SmartPointer<const Parameters, class Matrix<dim>>
parameters;
49 , direction(dim,
std::vector<double>(1))
60 AssertDimension(dinfo.n_matrices(), 1);
62 dinfo.matrix(0,
false).matrix, info.fe_values(0), info.fe_values(0),
direction);
69 dealii::LocalIntegrators::Advection::upwind_value_matrix(
70 dinfo.matrix(0,
false).matrix, info.fe_values(0), info.fe_values(0),
direction);
76 IntegrationInfo<dim>& info2)
const 78 dealii::LocalIntegrators::Advection::upwind_value_matrix(dinfo1.matrix(0,
false).matrix,
79 dinfo1.matrix(0,
true).matrix,
80 dinfo2.matrix(0,
true).matrix,
81 dinfo2.matrix(0,
false).matrix,
std::vector< std::vector< double > > direction
Definition: advection/matrix.h:43
Definition: advection/matrix.h:31
void cell_matrix(dealii::FullMatrix< double > &M, 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: matrix_integrators.h:23
virtual void cell(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: advection/matrix.h:58
Definition: advection/matrix.h:17
virtual void face(DoFInfo< dim > &dinfo1, DoFInfo< dim > &dinfo2, IntegrationInfo< dim > &info1, IntegrationInfo< dim > &info2) const
Definition: advection/matrix.h:75
virtual void boundary(DoFInfo< dim > &dinfo, IntegrationInfo< dim > &info) const
Definition: advection/matrix.h:67
Definition: integrator.h:29
dealii::SmartPointer< const Parameters, class Matrix< dim > > parameters
Definition: advection/matrix.h:42