Amandus: Simulations based on multilevel Schwarz methods
allen_cahn/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_allen_cahn_h
8 #define __matrix_allen_cahn_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 
19 namespace AllenCahn
20 {
21 using namespace dealii;
22 using namespace LocalIntegrators;
23 
24 template <int dim>
25 class Matrix : public AmandusIntegrator<dim>
26 {
27 public:
28  Matrix(double diffusion);
29 
30  virtual void cell(MeshWorker::DoFInfo<dim>& dinfo, MeshWorker::IntegrationInfo<dim>& info) const;
31  virtual void boundary(MeshWorker::DoFInfo<dim>& dinfo,
32  MeshWorker::IntegrationInfo<dim>& info) const;
33  virtual void face(MeshWorker::DoFInfo<dim>& dinfo1, MeshWorker::DoFInfo<dim>& dinfo2,
34  MeshWorker::IntegrationInfo<dim>& info1,
35  MeshWorker::IntegrationInfo<dim>& info2) const;
36 
37 private:
38  double D;
39 };
40 
41 template <int dim>
42 Matrix<dim>::Matrix(double diffusion)
43  : D(diffusion)
44 {
45  this->use_boundary = false;
46 }
47 
48 template <int dim>
49 void
50 Matrix<dim>::cell(MeshWorker::DoFInfo<dim>& dinfo, MeshWorker::IntegrationInfo<dim>& info) const
51 {
52  AssertDimension(dinfo.n_matrices(), 1);
53  // Assert (info.values.size() >0, ExcInternalError());
54  Laplace::cell_matrix(dinfo.matrix(0, false).matrix, info.fe_values(0), D);
55  if (info.values.size() > 0)
56  {
57  AssertDimension(info.values[0][0].size(), info.fe_values(0).n_quadrature_points);
58  std::vector<double> fu(info.fe_values(0).n_quadrature_points);
59  for (unsigned int k = 0; k < fu.size(); ++k)
60  {
61  const double u = info.values[0][0][k];
62  fu[k] = (3. * u * u - 1.);
63  }
64  L2::weighted_mass_matrix(dinfo.matrix(0, false).matrix, info.fe_values(0), fu);
65  }
66 }
67 
68 template <int dim>
69 void
70 Matrix<dim>::boundary(MeshWorker::DoFInfo<dim>&, typename MeshWorker::IntegrationInfo<dim>&) const
71 {
72 }
73 
74 template <int dim>
75 void
76 Matrix<dim>::face(MeshWorker::DoFInfo<dim>& dinfo1, MeshWorker::DoFInfo<dim>& dinfo2,
77  MeshWorker::IntegrationInfo<dim>& info1,
78  MeshWorker::IntegrationInfo<dim>& info2) const
79 {
80  const unsigned int deg = info1.fe_values(0).get_fe().tensor_degree();
81  Laplace::ip_matrix(dinfo1.matrix(0, false).matrix,
82  dinfo1.matrix(0, true).matrix,
83  dinfo2.matrix(0, true).matrix,
84  dinfo2.matrix(0, false).matrix,
85  info1.fe_values(0),
86  info2.fe_values(0),
87  Laplace::compute_penalty(dinfo1, dinfo2, deg, deg),
88  D);
89 }
90 }
91 
92 #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
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: allen_cahn/matrix.h:50
virtual void face(MeshWorker::DoFInfo< dim > &dinfo1, MeshWorker::DoFInfo< dim > &dinfo2, MeshWorker::IntegrationInfo< dim > &info1, MeshWorker::IntegrationInfo< dim > &info2) const
Definition: allen_cahn/matrix.h:76
virtual void boundary(MeshWorker::DoFInfo< dim > &dinfo, MeshWorker::IntegrationInfo< dim > &info) const
Definition: allen_cahn/matrix.h:70
Matrix(double diffusion)
Definition: allen_cahn/matrix.h:42
double D
Definition: allen_cahn/matrix.h:38
Definition: allen_cahn/matrix.h:25
Definition: allen_cahn/matrix.h:19
Definition: integrator.h:29