7 #ifndef __matrix_curl_curl_h 8 #define __matrix_curl_curl_h 10 #include <amandus/integrator.h> 11 #include <deal.II/integrators/divergence.h> 12 #include <deal.II/integrators/l2.h> 13 #include <deal.II/integrators/maxwell.h> 14 #include <deal.II/meshworker/integration_info.h> 69 virtual void cell(dealii::MeshWorker::DoFInfo<dim>& dinfo,
70 dealii::MeshWorker::IntegrationInfo<dim>& info)
const;
71 virtual void boundary(dealii::MeshWorker::DoFInfo<dim>& dinfo,
72 dealii::MeshWorker::IntegrationInfo<dim>& info)
const;
73 virtual void face(dealii::MeshWorker::DoFInfo<dim>& dinfo1,
74 dealii::MeshWorker::DoFInfo<dim>& dinfo2,
75 dealii::MeshWorker::IntegrationInfo<dim>& info1,
76 dealii::MeshWorker::IntegrationInfo<dim>& info2)
const;
84 : curl_coefficient(1, 1.)
85 , mass_coefficient(1, 0.)
87 this->use_boundary =
false;
88 this->use_face =
false;
94 dealii::MeshWorker::IntegrationInfo<dim>& info)
const 96 AssertDimension(dinfo.n_matrices(), 4);
97 const unsigned int id = dinfo.cell->material_id();
103 Maxwell::curl_curl_matrix(dinfo.matrix(0,
false).matrix, info.fe_values(0), mu);
105 L2::mass_matrix(dinfo.matrix(0,
false).matrix, info.fe_values(0), sigma);
106 Divergence::gradient_matrix(dinfo.matrix(1,
false).matrix, info.fe_values(1), info.fe_values(0));
107 dinfo.matrix(2,
false).matrix.copy_transposed(dinfo.matrix(1,
false).matrix);
113 typename dealii::MeshWorker::IntegrationInfo<dim>& )
const 120 dealii::MeshWorker::DoFInfo<dim>& ,
121 dealii::MeshWorker::IntegrationInfo<dim>& ,
122 dealii::MeshWorker::IntegrationInfo<dim>& )
const Definition: maxwell/eigen.h:22
virtual void cell(dealii::MeshWorker::DoFInfo< dim > &dinfo, dealii::MeshWorker::IntegrationInfo< dim > &info) const
Definition: maxwell/matrix.h:93
virtual void boundary(dealii::MeshWorker::DoFInfo< dim > &dinfo, dealii::MeshWorker::IntegrationInfo< dim > &info) const
Definition: maxwell/matrix.h:112
std::vector< double > curl_coefficient
Definition: maxwell/matrix.h:78
Definition: maxwell/matrix.h:61
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: maxwell/matrix.h:119
std::vector< double > mass_coefficient
Definition: maxwell/matrix.h:79
Definition: integrator.h:29