Amandus: Simulations based on multilevel Schwarz methods
laplace/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 __matrix_laplace_h
8 #define __matrix_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 
30 namespace LaplaceIntegrators
31 {
35 template <int dim>
36 class Matrix : public AmandusIntegrator<dim>
37 {
38 public:
39  virtual void cell(MeshWorker::DoFInfo<dim>& dinfo, MeshWorker::IntegrationInfo<dim>& info) const;
40  virtual void boundary(MeshWorker::DoFInfo<dim>& dinfo,
41  MeshWorker::IntegrationInfo<dim>& info) const;
42  virtual void face(MeshWorker::DoFInfo<dim>& dinfo1, MeshWorker::DoFInfo<dim>& dinfo2,
43  MeshWorker::IntegrationInfo<dim>& info1,
44  MeshWorker::IntegrationInfo<dim>& info2) const;
45 };
46 
47 template <int dim>
48 void
49 Matrix<dim>::cell(MeshWorker::DoFInfo<dim>& dinfo, MeshWorker::IntegrationInfo<dim>& info) const
50 {
51  AssertDimension(dinfo.n_matrices(), 1);
52  Laplace::cell_matrix(dinfo.matrix(0, false).matrix, info.fe_values(0));
53 }
54 
55 template <int dim>
56 void
57 Matrix<dim>::boundary(MeshWorker::DoFInfo<dim>& dinfo,
58  typename MeshWorker::IntegrationInfo<dim>& info) const
59 {
60  if (info.fe_values(0).get_fe().conforms(FiniteElementData<dim>::H1))
61  return;
62 
63  const unsigned int deg = info.fe_values(0).get_fe().tensor_degree();
64  Laplace::nitsche_matrix(dinfo.matrix(0, false).matrix,
65  info.fe_values(0),
66  Laplace::compute_penalty(dinfo, dinfo, deg, deg));
67 }
68 
69 template <int dim>
70 void
71 Matrix<dim>::face(MeshWorker::DoFInfo<dim>& dinfo1, MeshWorker::DoFInfo<dim>& dinfo2,
72  MeshWorker::IntegrationInfo<dim>& info1,
73  MeshWorker::IntegrationInfo<dim>& info2) const
74 {
75  if (info1.fe_values(0).get_fe().conforms(FiniteElementData<dim>::H1))
76  return;
77 
78  const unsigned int deg = info1.fe_values(0).get_fe().tensor_degree();
79  Laplace::ip_matrix(dinfo1.matrix(0, false).matrix,
80  dinfo1.matrix(0, true).matrix,
81  dinfo2.matrix(0, true).matrix,
82  dinfo2.matrix(0, false).matrix,
83  info1.fe_values(0),
84  info2.fe_values(0),
85  Laplace::compute_penalty(dinfo1, dinfo2, deg, deg));
86 }
87 }
88 
89 #endif
Integrator for the matrix of the Laplace operator.
Definition: laplace/matrix.h:36
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
Definition: laplace/eigen.h:19
Definition: integrator.h:29