Amandus: Simulations based on multilevel Schwarz methods
matrix_factor.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 __matrixfaktor_laplace_h
8 #define __matrixfaktor_laplace_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 
16 using namespace dealii;
17 using namespace LocalIntegrators;
18 
19 namespace LaplaceIntegrators
20 {
32 template <int dim>
33 class MatrixFaktor : public AmandusIntegrator<dim>
34 {
35 public:
36  MatrixFaktor(double faktor);
37 
38  virtual void cell(MeshWorker::DoFInfo<dim>& dinfo, MeshWorker::IntegrationInfo<dim>& info) const;
39  virtual void boundary(MeshWorker::DoFInfo<dim>& dinfo,
40  MeshWorker::IntegrationInfo<dim>& info) const;
41  virtual void face(MeshWorker::DoFInfo<dim>& dinfo1, MeshWorker::DoFInfo<dim>& dinfo2,
42  MeshWorker::IntegrationInfo<dim>& info1,
43  MeshWorker::IntegrationInfo<dim>& info2) const;
44 
45 private:
46  double faktor;
47 };
48 
49 template <int dim>
51  : faktor(faktor)
52 {
53 }
54 
55 template <int dim>
56 void
57 MatrixFaktor<dim>::cell(MeshWorker::DoFInfo<dim>& dinfo,
58  MeshWorker::IntegrationInfo<dim>& info) const
59 {
60  AssertDimension(dinfo.n_matrices(), 1);
61  Laplace::cell_matrix(dinfo.matrix(0, false).matrix, info.fe_values(0), faktor);
62 }
63 
64 template <int dim>
65 void
66 MatrixFaktor<dim>::boundary(MeshWorker::DoFInfo<dim>& dinfo,
67  typename MeshWorker::IntegrationInfo<dim>& info) const
68 {
69  const unsigned int deg = info.fe_values(0).get_fe().tensor_degree();
70  Laplace::nitsche_matrix(dinfo.matrix(0, false).matrix,
71  info.fe_values(0),
72  Laplace::compute_penalty(dinfo, dinfo, deg, deg),
73  faktor);
74 }
75 
76 template <int dim>
77 void
78 MatrixFaktor<dim>::face(MeshWorker::DoFInfo<dim>& dinfo1, MeshWorker::DoFInfo<dim>& dinfo2,
79  MeshWorker::IntegrationInfo<dim>& info1,
80  MeshWorker::IntegrationInfo<dim>& info2) const
81 {
82  const unsigned int deg = info1.fe_values(0).get_fe().tensor_degree();
83  Laplace::ip_matrix(dinfo1.matrix(0, false).matrix,
84  dinfo1.matrix(0, true).matrix,
85  dinfo2.matrix(0, true).matrix,
86  dinfo2.matrix(0, false).matrix,
87  info1.fe_values(0),
88  info2.fe_values(0),
89  Laplace::compute_penalty(dinfo1, dinfo2, deg, deg),
90  faktor);
91 }
92 }
93 
94 #endif
double faktor
Definition: matrix_factor.h:46
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: matrix_factor.h:57
virtual void face(MeshWorker::DoFInfo< dim > &dinfo1, MeshWorker::DoFInfo< dim > &dinfo2, MeshWorker::IntegrationInfo< dim > &info1, MeshWorker::IntegrationInfo< dim > &info2) const
Definition: matrix_factor.h:78
Definition: laplace/eigen.h:19
virtual void boundary(MeshWorker::DoFInfo< dim > &dinfo, MeshWorker::IntegrationInfo< dim > &info) const
Definition: matrix_factor.h:66
Definition: matrix_factor.h:33
Definition: integrator.h:29