17 #include <deal.II/base/quadrature.h> 18 #include <deal.II/base/qprojector.h> 19 #include <deal.II/base/template_constraints.h> 20 #include <deal.II/base/std_cxx11/unique_ptr.h> 21 #include <deal.II/fe/fe_q_bubbles.h> 22 #include <deal.II/fe/fe_nothing.h> 23 #include <deal.II/fe/fe_tools.h> 24 #include <deal.II/fe/fe_values.h> 25 #include <deal.II/fe/mapping_q1.h> 26 #include <deal.II/base/quadrature_lib.h> 27 #include <deal.II/dofs/dof_accessor.h> 28 #include <deal.II/dofs/dof_handler.h> 29 #include <deal.II/grid/tria.h> 31 #include <deal.II/grid/grid_generator.h> 37 DEAL_II_NAMESPACE_OPEN
44 template <
int dim,
int spacedim>
49 const bool isotropic_only)
52 const unsigned int degree = fe.
degree;
55 std_cxx11::unique_ptr<Quadrature<dim> > q_fine;
81 const unsigned int nq = q_fine->size();
84 unsigned int ref_case = (isotropic_only)
87 for (; ref_case <= RefinementCase<dim>::isotropic_refinement; ++ref_case)
92 for (
unsigned int i=0; i<nc; ++i)
94 Assert(matrices[ref_case-1][i].n() == dpc,
96 Assert(matrices[ref_case-1][i].m() == dpc,
107 dh.distribute_dofs(fe);
113 const unsigned int n_dofs = dh.n_dofs();
118 std::vector<std::vector<types::global_dof_index> > child_ldi
119 (nc, std::vector<types::global_dof_index>(fe.
dofs_per_cell));
122 unsigned int child_no = 0;
123 typename ::DoFHandler<dim>::active_cell_iterator cell
125 for (; cell!=dh.end(); ++cell, ++child_no)
128 cell->get_dof_indices(child_ldi[child_no]);
130 for (
unsigned int q=0; q<nq; ++q)
131 for (
unsigned int i=0; i<dpc; ++i)
132 for (
unsigned int j=0; j<dpc; ++j)
134 const unsigned int gdi=child_ldi[child_no][i];
135 const unsigned int gdj=child_ldi[child_no][j];
136 fine_mass(gdi, gdj)+=fine.shape_value(i,q)
137 *fine.shape_value(j,q)
140 for (
unsigned int k=0; k<dim; ++k)
141 quad_tmp(k) = fine.quadrature_point(q)(k);
142 coarse_rhs_matrix(gdi, j)
143 +=fine.shape_value(i,q)
151 fine_mass.gauss_jordan();
152 fine_mass.mmult(solution, coarse_rhs_matrix);
155 for (
unsigned int child_no=0; child_no<nc; ++child_no)
156 for (
unsigned int i=0; i<dpc; ++i)
157 for (
unsigned int j=0; j<dpc; ++j)
159 const unsigned int gdi=child_ldi[child_no][i];
161 if (std::fabs(solution(gdi, j)) > 1.e-12)
162 matrices[ref_case-1][child_no](i,j)=solution(gdi, j);
170 template <
int dim,
int spacedim>
178 get_riaf_vector(q_degree)),
179 n_bubbles((q_degree<=1)?1:dim)
182 ExcMessage (
"This element can only be used for polynomial degrees " 183 "greater than zero"));
189 for (
unsigned int d=0; d<dim; ++d)
198 FE_Q_Bubbles_Helper::compute_embedding_matrices
208 template <
int dim,
int spacedim>
216 get_riaf_vector(points.size()-1)),
217 n_bubbles((points.size()-1<=1)?1:dim)
220 ExcMessage (
"This element can only be used for polynomial degrees " 227 for (
unsigned int d=0; d<dim; ++d)
229 for (
unsigned int i=0; i<
n_bubbles; ++i)
236 FE_Q_Bubbles_Helper::compute_embedding_matrices
245 template <
int dim,
int spacedim>
253 std::ostringstream namebuf;
255 const unsigned int n_points = this->degree;
256 std::vector<double> points(n_points);
257 const unsigned int dofs_per_cell = this->dofs_per_cell;
258 const std::vector<Point<dim> > &unit_support_points = this->unit_support_points;
259 unsigned int index = 0;
262 for (
unsigned int j=0; j<dofs_per_cell; j++)
264 if ((dim>1) ? (unit_support_points[j](1)==0 &&
265 ((dim>2) ? unit_support_points[j](2)==0:
true)) :
true)
268 points[index] = unit_support_points[j](0);
270 points[n_points-1] = unit_support_points[j](0);
272 points[index-1] = unit_support_points[j](0);
278 Assert (index == n_points || (dim==1 && index == n_points+n_bubbles),
279 ExcMessage (
"Could not decode support points in one coordinate direction."));
282 for (
unsigned int j=0; j<n_points; j++)
283 if (std::fabs(points[j] - (
double)j/(this->degree-1)) > 1e-15)
291 if (this->degree > 3)
292 namebuf <<
"FE_Q_Bubbles<" 294 <<
">(QIterated(QTrapez()," << this->degree-1 <<
"))";
296 namebuf <<
"FE_Q_Bubbles<" 298 <<
">(" << this->degree-1 <<
")";
306 for (
unsigned int j=0; j<n_points; j++)
307 if (points[j] != points_gl.
point(j)(0))
313 namebuf <<
"FE_Q_Bubbles<" << dim <<
">(" << this->degree-1 <<
")";
315 namebuf <<
"FE_Q_Bubbles<" << dim <<
">(QUnknownNodes(" << this->degree <<
"))";
317 return namebuf.str();
322 template <
int dim,
int spacedim>
331 template <
int dim,
int spacedim>
335 std::vector<double> &nodal_values)
const 337 Assert (support_point_values.size() == this->unit_support_points.size(),
339 this->unit_support_points.size()));
340 Assert (nodal_values.size() == this->dofs_per_cell,
345 for (
unsigned int i=0; i<this->dofs_per_cell-1; ++i)
347 const std::pair<unsigned int, unsigned int> index
348 = this->system_to_component_index(i);
349 nodal_values[i] = support_point_values[i](index.first);
353 for (
unsigned int i = 0; i<n_bubbles; ++i)
354 nodal_values[nodal_values.size()-i-1] = 0.;
359 template <
int dim,
int spacedim>
362 const std::vector<double> &values)
const 364 Assert (values.size() == this->unit_support_points.size(),
366 this->unit_support_points.size()));
367 Assert (local_dofs.size() == this->dofs_per_cell,
372 std::copy(values.begin(), values.end(), local_dofs.begin());
375 for (
unsigned int i = 0; i<n_bubbles; ++i)
376 local_dofs[local_dofs.size()-i-1] = 0.;
381 template <
int dim,
int spacedim>
385 const unsigned int offset)
const 387 Assert (values.size() == this->unit_support_points.size(),
389 this->unit_support_points.size()));
390 Assert (local_dofs.size() == this->dofs_per_cell,
395 for (
unsigned int i=0; i<this->dofs_per_cell-1; ++i)
397 const std::pair<unsigned int, unsigned int> index
398 = this->system_to_component_index(i);
399 local_dofs[i] = values[i](offset+index.first);
403 for (
unsigned int i = 0; i<n_bubbles; ++i)
404 local_dofs[local_dofs.size()-i-1] = 0.;
409 template <
int dim,
int spacedim>
412 std::vector<double> &local_dofs,
413 const VectorSlice<
const std::vector<std::vector<double> > > &values)
const 415 Assert (values[0].size() == this->unit_support_points.size(),
417 this->unit_support_points.size()));
418 Assert (local_dofs.size() == this->dofs_per_cell,
423 for (
unsigned int i=0; i<this->dofs_per_cell-1; ++i)
425 const std::pair<unsigned int, unsigned int> index
426 = this->system_to_component_index(i);
427 local_dofs[i] = values[index.first][i];
431 for (
unsigned int i = 0; i<n_bubbles; ++i)
432 local_dofs[local_dofs.size()-i-1] = 0.;
437 template <
int dim,
int spacedim>
451 (dynamic_cast<const FEQBUBBLES *>(&x_source_fe) != 0)
454 ExcInterpolationNotImplemented());
455 Assert (interpolation_matrix.
m() == this->dofs_per_cell,
457 this->dofs_per_cell));
463 if (dynamic_cast<const FEQBUBBLES *>(&x_source_fe)->degree == this->degree)
464 for (
unsigned int i=0; i<interpolation_matrix.
m(); ++i)
465 interpolation_matrix.
set(i,i,1.);
468 Assert(
false,
typename FEL::ExcInterpolationNotImplemented())
473 template <
int dim,
int spacedim>
477 unsigned int n_cont_dofs = Utilities::fixed_power<dim>(q_deg+1);
478 const unsigned int n_bubbles = (q_deg<=1?1:dim);
479 return std::vector<bool> (n_cont_dofs+n_bubbles,
true);
484 template <
int dim,
int spacedim>
485 std::vector<unsigned int>
488 std::vector<unsigned int> dpo(dim+1, 1U);
489 for (
unsigned int i=1; i<dpo.size(); ++i)
490 dpo[i]=dpo[i-1]*(q_deg-1);
492 dpo[dim]+=(q_deg<=1?1:dim);
498 template <
int dim,
int spacedim>
501 (
const unsigned int shape_index,
502 const unsigned int face_index)
const 505 if (shape_index >= this->n_dofs_per_cell()-n_bubbles)
513 template <
int dim,
int spacedim>
516 (
const unsigned int child,
522 ExcMessage(
"Prolongation matrices are only available for refined cells!"));
526 Assert (this->prolongation[refinement_case-1][child].n() != 0,
527 ExcMessage(
"This prolongation matrix has not been computed yet!"));
529 return this->prolongation[refinement_case-1][child];
534 template <
int dim,
int spacedim>
537 (
const unsigned int child,
543 ExcMessage(
"Restriction matrices are only available for refined cells!"));
547 Assert(this->restriction[refinement_case-1][child].n() != 0,
548 ExcMessage(
"This restriction matrix has not been computed yet!"));
551 return this->restriction[refinement_case-1][child];
555 #include "fe_q_bubbles.inst" 557 DEAL_II_NAMESPACE_CLOSE
Transformed quadrature weights.
const unsigned int n_bubbles
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const 1
std::vector< std::vector< FullMatrix< double > > > restriction
#define AssertDimension(dim1, dim2)
const std::vector< Point< dim > > & get_points() const
const unsigned int degree
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
const Point< dim > & point(const unsigned int i) const
active_cell_iterator begin_active(const unsigned int level=0) const
Transformed quadrature points.
#define AssertThrow(cond, exc)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
std::vector< Point< dim > > unit_support_points
void set(const size_type i, const size_type j, const number value)
virtual void execute_coarsening_and_refinement()
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case) const
static ::ExceptionBase & ExcMessage(std::string arg1)
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
std::vector< std::vector< FullMatrix< double > > > prolongation
#define Assert(cond, exc)
virtual FiniteElement< dim, spacedim > * clone() const
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case) const
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual std::string get_name() const =0
virtual std::string get_name() const
std::string dim_string(const int dim, const int spacedim)
unsigned int size() const
const unsigned int dofs_per_cell
void initialize(const std::vector< Point< 1 > > &support_points_1d)
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
static std::vector< bool > get_riaf_vector(const unsigned int degree)
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
FE_Q_Bubbles(const unsigned int p)
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
virtual void convert_generalized_support_point_values_to_nodal_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
static ::ExceptionBase & ExcInternalError()