Amandus: Simulations based on multilevel Schwarz methods
readiff/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 __readiff_matrix_h
8 #define __readiff_matrix_h
9 
10 #include <amandus/integrator.h>
11 #include <deal.II/integrators/divergence.h>
12 #include <deal.II/integrators/l2.h>
13 #include <deal.II/integrators/laplace.h>
14 #include <deal.II/meshworker/integration_info.h>
15 #include <readiff/parameters.h>
16 
17 using namespace dealii;
18 using namespace LocalIntegrators;
19 
32 {
38 template <int dim>
39 class Matrix : public AmandusIntegrator<dim>
40 {
41 public:
42  Matrix(const Parameters& par);
43 
44  virtual void cell(MeshWorker::DoFInfo<dim>& dinfo, MeshWorker::IntegrationInfo<dim>& info) const;
45  virtual void boundary(MeshWorker::DoFInfo<dim>& dinfo,
46  MeshWorker::IntegrationInfo<dim>& info) const;
47  virtual void face(MeshWorker::DoFInfo<dim>& dinfo1, MeshWorker::DoFInfo<dim>& dinfo2,
48  MeshWorker::IntegrationInfo<dim>& info1,
49  MeshWorker::IntegrationInfo<dim>& info2) const;
50 
51 private:
52  SmartPointer<const Parameters, class Matrix<dim>> parameters;
53 };
54 
55 template <int dim>
56 Matrix<dim>::Matrix(const Parameters& par)
57  : parameters(&par)
58 {
59  this->use_boundary = false;
60 }
61 
62 template <int dim>
63 void
64 Matrix<dim>::cell(MeshWorker::DoFInfo<dim>& dinfo, MeshWorker::IntegrationInfo<dim>& info) const
65 {
66  AssertDimension(dinfo.n_matrices(), 4);
67  // Assert (info.values.size() >0, ExcInternalError());
68 
69  const Parameters& p = *parameters;
70  Laplace::cell_matrix(dinfo.matrix(0, false).matrix, info.fe_values(0), p.alpha1);
71  Laplace::cell_matrix(dinfo.matrix(3, false).matrix, info.fe_values(0), p.alpha2);
72  if (info.values.size() > 0)
73  {
74  AssertDimension(info.values[0].size(), 2);
75  AssertDimension(info.values[0][0].size(), info.fe_values(0).n_quadrature_points);
76  AssertDimension(info.values[0][1].size(), info.fe_values(0).n_quadrature_points);
77  std::vector<double> Du_ru(info.fe_values(0).n_quadrature_points);
78  std::vector<double> Dv_ru(info.fe_values(0).n_quadrature_points);
79  std::vector<double> Dv_rv(info.fe_values(0).n_quadrature_points);
80  std::vector<double> Du_rv(info.fe_values(0).n_quadrature_points);
81  for (unsigned int k = 0; k < Du_ru.size(); ++k)
82  {
83  const double u = info.values[0][0][k];
84  const double v = info.values[0][1][k];
85  Du_ru[k] = p.B1 + p.D1 * v + 2. * p.E1 * u + 2. * p.G1 * u * v + p.H1 * v * v;
86  Dv_ru[k] = p.C1 + p.D1 * u + 2. * p.F1 * v + 2. * p.H1 * u * v + p.G1 * u * u;
87  Dv_rv[k] = p.C2 + p.D2 * u + 2. * p.F2 * v + 2. * p.H2 * u * v + p.G2 * u * u;
88  Du_rv[k] = p.B2 + p.D2 * v + 2. * p.E2 * u + 2. * p.G2 * u * v + p.H2 * v * v;
89  }
90  L2::weighted_mass_matrix(dinfo.matrix(0, false).matrix, info.fe_values(0), Du_ru);
91  L2::weighted_mass_matrix(dinfo.matrix(1, false).matrix, info.fe_values(0), Dv_ru);
92  L2::weighted_mass_matrix(dinfo.matrix(2, false).matrix, info.fe_values(0), Du_rv);
93  L2::weighted_mass_matrix(dinfo.matrix(3, false).matrix, info.fe_values(0), Dv_rv);
94  }
95 }
96 
97 template <int dim>
98 void
99 Matrix<dim>::boundary(MeshWorker::DoFInfo<dim>&, typename MeshWorker::IntegrationInfo<dim>&) const
100 {
101 }
102 
103 template <int dim>
104 void
105 Matrix<dim>::face(MeshWorker::DoFInfo<dim>& dinfo1, MeshWorker::DoFInfo<dim>& dinfo2,
106  MeshWorker::IntegrationInfo<dim>& info1,
107  MeshWorker::IntegrationInfo<dim>& info2) const
108 {
109  const unsigned int deg = info1.fe_values(0).get_fe().tensor_degree();
110  const Parameters& p = *parameters;
111  Laplace::ip_matrix(dinfo1.matrix(0, false).matrix,
112  dinfo1.matrix(0, true).matrix,
113  dinfo2.matrix(0, true).matrix,
114  dinfo2.matrix(0, false).matrix,
115  info1.fe_values(0),
116  info2.fe_values(0),
117  Laplace::compute_penalty(dinfo1, dinfo2, deg, deg),
118  p.alpha1);
119  Laplace::ip_matrix(dinfo1.matrix(3, false).matrix,
120  dinfo1.matrix(3, true).matrix,
121  dinfo2.matrix(3, true).matrix,
122  dinfo2.matrix(3, false).matrix,
123  info1.fe_values(0),
124  info2.fe_values(0),
125  Laplace::compute_penalty(dinfo1, dinfo2, deg, deg),
126  p.alpha2);
127 }
128 }
129 
130 #endif
void weighted_mass_matrix(dealii::FullMatrix< double > &M, const dealii::FEValuesBase< dim > &fe, const dealii::TensorFunction< 2, dim > &weight)
Definition: darcy/integrators.h:41
Definition: readiff/matrix.h:39
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(MeshWorker::DoFInfo< dim > &dinfo, MeshWorker::IntegrationInfo< dim > &info) const
Definition: readiff/matrix.h:64
virtual void face(MeshWorker::DoFInfo< dim > &dinfo1, MeshWorker::DoFInfo< dim > &dinfo2, MeshWorker::IntegrationInfo< dim > &info1, MeshWorker::IntegrationInfo< dim > &info2) const
Definition: readiff/matrix.h:105
Definition: readiff/matrix.h:31
virtual void boundary(MeshWorker::DoFInfo< dim > &dinfo, MeshWorker::IntegrationInfo< dim > &info) const
Definition: readiff/matrix.h:99
SmartPointer< const Parameters, class Matrix< dim > > parameters
Definition: readiff/matrix.h:52
Definition: integrator.h:29