Amandus: Simulations based on multilevel Schwarz methods
darcy/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 _amandus_matrix_darcy_h
8 #define _amandus_matrix_darcy_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/meshworker/integration_info.h>
14 
15 using namespace dealii;
16 using namespace LocalIntegrators;
17 
21 namespace DarcyIntegrators
22 {
23 template <int dim>
24 class Matrix : public AmandusIntegrator<dim>
25 {
26 public:
27  Matrix();
28  virtual void cell(dealii::MeshWorker::DoFInfo<dim>& dinfo,
29  dealii::MeshWorker::IntegrationInfo<dim>& info) const;
30  virtual void boundary(dealii::MeshWorker::DoFInfo<dim>& dinfo,
31  dealii::MeshWorker::IntegrationInfo<dim>& info) const;
32  virtual void face(dealii::MeshWorker::DoFInfo<dim>& dinfo1,
33  dealii::MeshWorker::DoFInfo<dim>& dinfo2,
34  dealii::MeshWorker::IntegrationInfo<dim>& info1,
35  dealii::MeshWorker::IntegrationInfo<dim>& info2) const;
36 
37  std::vector<double> resistance;
38 };
39 
40 template <int dim>
42  : resistance(1, 1.)
43 {
44  this->use_boundary = false;
45  this->use_face = false;
46 }
47 
48 template <int dim>
49 void
50 Matrix<dim>::cell(dealii::MeshWorker::DoFInfo<dim>& dinfo,
51  dealii::MeshWorker::IntegrationInfo<dim>& info) const
52 {
53  AssertDimension(dinfo.n_matrices(), 4);
54  const unsigned int id = dinfo.cell->material_id();
55  AssertIndexRange(id, resistance.size());
56  const double R = resistance[id];
57 
58  L2::mass_matrix(dinfo.matrix(0, false).matrix, info.fe_values(0), R);
59  Divergence::cell_matrix(dinfo.matrix(2, false).matrix, info.fe_values(0), info.fe_values(1));
60  dinfo.matrix(1, false).matrix.copy_transposed(dinfo.matrix(2, false).matrix);
61 }
62 
63 template <int dim>
64 void
65 Matrix<dim>::boundary(dealii::MeshWorker::DoFInfo<dim>& /*dinfo*/,
66  typename dealii::MeshWorker::IntegrationInfo<dim>& /*info*/) const
67 {
68 }
69 
70 template <int dim>
71 void
72 Matrix<dim>::face(dealii::MeshWorker::DoFInfo<dim>& /*dinfo1*/,
73  dealii::MeshWorker::DoFInfo<dim>& /*dinfo2*/,
74  dealii::MeshWorker::IntegrationInfo<dim>& /*info1*/,
75  dealii::MeshWorker::IntegrationInfo<dim>& /*info2*/) const
76 {
77 }
78 }
79 
80 #endif
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 face(dealii::MeshWorker::DoFInfo< dim > &dinfo1, dealii::MeshWorker::DoFInfo< dim > &dinfo2, dealii::MeshWorker::IntegrationInfo< dim > &info1, dealii::MeshWorker::IntegrationInfo< dim > &info2) const
Definition: darcy/matrix.h:72
std::vector< double > resistance
Definition: darcy/matrix.h:37
virtual void cell(dealii::MeshWorker::DoFInfo< dim > &dinfo, dealii::MeshWorker::IntegrationInfo< dim > &info) const
Definition: darcy/matrix.h:50
virtual void boundary(dealii::MeshWorker::DoFInfo< dim > &dinfo, dealii::MeshWorker::IntegrationInfo< dim > &info) const
Definition: darcy/matrix.h:65
Definition: estimator.h:33
Definition: darcy/matrix.h:24
Definition: integrator.h:29