Amandus: Simulations based on multilevel Schwarz methods
advection/matrix.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 #ifndef __advection_matrix_h
8 #define __advection_matrix_h
9 
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>
14 
15 using namespace dealii::MeshWorker;
16 
17 namespace Advection
18 {
30 template <int dim>
31 class Matrix : public AmandusIntegrator<dim>
32 {
33 public:
34  Matrix(const Parameters& par);
35 
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;
40 
41 private:
42  dealii::SmartPointer<const Parameters, class Matrix<dim>> parameters;
43  std::vector<std::vector<double>> direction;
44 };
45 
46 template <int dim>
47 Matrix<dim>::Matrix(const Parameters& par)
48  : parameters(&par)
49  , direction(dim, std::vector<double>(1))
50 {
51  // this->input_vector_names.push_back("Newton iterate");
52  direction[0][0] = 1.;
53  direction[1][0] = 2.;
54 }
55 
56 template <int dim>
57 void
58 Matrix<dim>::cell(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
59 {
60  AssertDimension(dinfo.n_matrices(), 1);
62  dinfo.matrix(0, false).matrix, info.fe_values(0), info.fe_values(0), direction);
63 }
64 
65 template <int dim>
66 void
67 Matrix<dim>::boundary(DoFInfo<dim>& dinfo, IntegrationInfo<dim>& info) const
68 {
69  dealii::LocalIntegrators::Advection::upwind_value_matrix(
70  dinfo.matrix(0, false).matrix, info.fe_values(0), info.fe_values(0), direction);
71 }
72 
73 template <int dim>
74 void
75 Matrix<dim>::face(DoFInfo<dim>& dinfo1, DoFInfo<dim>& dinfo2, IntegrationInfo<dim>& info1,
76  IntegrationInfo<dim>& info2) const
77 {
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,
82  info1.fe_values(0),
83  info2.fe_values(0),
84  info1.fe_values(0),
85  info2.fe_values(0),
86  direction);
87 }
88 }
89 
90 #endif
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
STL namespace.
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