49 template <
int dim,
int spacedim>
52 const ::FE_Q_Bubbles<dim, spacedim> & fe,
54 const bool isotropic_only)
56 const unsigned int dpc = fe.dofs_per_cell;
57 const unsigned int degree = fe.degree;
60 std::unique_ptr<Quadrature<dim>> q_fine;
62 std::vector<double>(1, 1.));
67 q_fine = std_cxx14::make_unique<QGauss<dim>>(degree + 1);
68 else if (spacedim == 2)
69 q_fine = std_cxx14::make_unique<QAnisotropic<dim>>(
72 q_fine = std_cxx14::make_unique<QAnisotropic<dim>>(
77 q_fine = std_cxx14::make_unique<QGauss<dim>>(degree + 1);
79 q_fine = std_cxx14::make_unique<QAnisotropic<dim>>(
83 q_fine = std_cxx14::make_unique<QGauss<dim>>(degree + 1);
90 const unsigned int nq = q_fine->size();
93 unsigned int ref_case = (isotropic_only) ?
96 for (; ref_case <= RefinementCase<dim>::isotropic_refinement;
99 const unsigned int nc =
102 for (
unsigned int i = 0; i < nc; ++i)
104 Assert(matrices[ref_case - 1][i].n() == dpc,
107 Assert(matrices[ref_case - 1][i].m() == dpc,
127 const unsigned int n_dofs = dh.
n_dofs();
132 std::vector<std::vector<types::global_dof_index>> child_ldi(
133 nc, std::vector<types::global_dof_index>(fe.dofs_per_cell));
136 unsigned int child_no = 0;
137 typename ::DoFHandler<dim>::active_cell_iterator cell =
139 for (; cell != dh.
end(); ++cell, ++child_no)
142 cell->get_dof_indices(child_ldi[child_no]);
144 for (
unsigned int q = 0; q < nq; ++q)
145 for (
unsigned int i = 0; i < dpc; ++i)
146 for (
unsigned int j = 0; j < dpc; ++j)
148 const unsigned int gdi = child_ldi[child_no][i];
149 const unsigned int gdj = child_ldi[child_no][j];
150 fine_mass(gdi, gdj) += fine.shape_value(i, q) *
151 fine.shape_value(j, q) *
154 for (
unsigned int k = 0; k < dim; ++k)
155 quad_tmp(k) = fine.quadrature_point(q)(k);
156 coarse_rhs_matrix(gdi, j) +=
157 fine.shape_value(i, q) * fe.shape_value(j, quad_tmp) *
165 fine_mass.
mmult(solution, coarse_rhs_matrix);
168 for (
unsigned int child_no = 0; child_no < nc; ++child_no)
169 for (
unsigned int i = 0; i < dpc; ++i)
170 for (
unsigned int j = 0; j < dpc; ++j)
172 const unsigned int gdi = child_ldi[child_no][i];
175 matrices[ref_case - 1][child_no](i, j) = solution(gdi, j);
184 template <
int dim,
int spacedim>
194 get_riaf_vector(q_degree))
195 , n_bubbles((q_degree <= 1) ? 1 : dim)
198 ExcMessage(
"This element can only be used for polynomial degrees "
199 "greater than zero"));
205 for (
unsigned int d = 0;
d < dim; ++
d)
207 for (
unsigned int i = 0; i <
n_bubbles; ++i)
224 template <
int dim,
int spacedim>
233 get_riaf_vector(points.size() - 1))
234 , n_bubbles((points.size() - 1 <= 1) ? 1 : dim)
237 ExcMessage(
"This element can only be used for polynomial degrees "
244 for (
unsigned int d = 0;
d < dim; ++
d)
246 for (
unsigned int i = 0; i <
n_bubbles; ++i)
263 template <
int dim,
int spacedim>
271 std::ostringstream namebuf;
273 const unsigned int n_points = this->degree;
274 std::vector<double> points(n_points);
275 const unsigned int dofs_per_cell = this->dofs_per_cell;
276 const std::vector<Point<dim>> &unit_support_points =
277 this->unit_support_points;
278 unsigned int index = 0;
281 for (
unsigned int j = 0; j < dofs_per_cell; j++)
283 if ((dim > 1) ? (unit_support_points[j](1) == 0 &&
284 ((dim > 2) ? unit_support_points[j](2) == 0 :
true)) :
288 points[index] = unit_support_points[j](0);
290 points[n_points - 1] = unit_support_points[j](0);
292 points[index - 1] = unit_support_points[j](0);
298 Assert(index == n_points || (dim == 1 && index == n_points + n_bubbles),
300 "Could not decode support points in one coordinate direction."));
303 for (
unsigned int j = 0; j < n_points; j++)
304 if (
std::fabs(points[j] -
static_cast<double>(j) / (this->degree - 1)) >
313 if (this->degree > 3)
315 <<
">(QIterated(QTrapez()," << this->degree - 1 <<
"))";
318 <<
">(" << this->degree - 1 <<
")";
325 for (
unsigned int j = 0; j < n_points; j++)
326 if (points[j] != points_gl.
point(j)(0))
332 namebuf <<
"FE_Q_Bubbles<" << dim <<
">(" << this->degree - 1 <<
")";
334 namebuf <<
"FE_Q_Bubbles<" << dim <<
">(QUnknownNodes(" << this->degree
337 return namebuf.str();
342 template <
int dim,
int spacedim>
343 std::unique_ptr<FiniteElement<dim, spacedim>>
346 return std_cxx14::make_unique<FE_Q_Bubbles<dim, spacedim>>(*this);
351 template <
int dim,
int spacedim>
356 std::vector<double> & nodal_values)
const
358 Assert(support_point_values.size() == this->unit_support_points.size(),
360 this->unit_support_points.size()));
361 Assert(nodal_values.size() == this->dofs_per_cell,
367 for (
unsigned int i = 0; i < this->dofs_per_cell - 1; ++i)
369 const std::pair<unsigned int, unsigned int> index =
370 this->system_to_component_index(i);
371 nodal_values[i] = support_point_values[i](index.first);
375 for (
unsigned int i = 0; i < n_bubbles; ++i)
376 nodal_values[nodal_values.size() - i - 1] = 0.;
381 template <
int dim,
int spacedim>
393 (x_source_fe.
get_name().find(
"FE_Q_Bubbles<") == 0) ||
394 (
dynamic_cast<const FEQBUBBLES *
>(&x_source_fe) !=
nullptr),
396 Assert(interpolation_matrix.
m() == this->dofs_per_cell,
403 auto casted_fe =
dynamic_cast<const FEQBUBBLES *
>(&x_source_fe);
404 if (casted_fe !=
nullptr && casted_fe->degree == this->degree)
405 for (
unsigned int i = 0; i < interpolation_matrix.
m(); ++i)
406 interpolation_matrix.
set(i, i, 1.);
412 spacedim>::ExcInterpolationNotImplemented()));
417 template <
int dim,
int spacedim>
421 unsigned int n_cont_dofs = Utilities::fixed_power<dim>(q_deg + 1);
422 const unsigned int n_bubbles = (q_deg <= 1 ? 1 : dim);
423 return std::vector<bool>(n_cont_dofs + n_bubbles,
true);
428 template <
int dim,
int spacedim>
429 std::vector<unsigned int>
432 std::vector<unsigned int> dpo(dim + 1, 1
U);
433 for (
unsigned int i = 1; i < dpo.size(); ++i)
434 dpo[i] = dpo[i - 1] * (q_deg - 1);
437 (q_deg <= 1 ? 1 : dim);
443 template <
int dim,
int spacedim>
446 const unsigned int shape_index,
447 const unsigned int face_index)
const
450 if (shape_index >= this->n_dofs_per_cell() - n_bubbles)
454 has_support_on_face(shape_index, face_index);
459 template <
int dim,
int spacedim>
462 const unsigned int child,
469 "Prolongation matrices are only available for refined cells!"));
472 Assert(this->prolongation[refinement_case - 1][child].n() != 0,
473 ExcMessage(
"This prolongation matrix has not been computed yet!"));
475 return this->prolongation[refinement_case - 1][child];
480 template <
int dim,
int spacedim>
483 const unsigned int child,
490 "Restriction matrices are only available for refined cells!"));
493 Assert(this->restriction[refinement_case - 1][child].n() != 0,
494 ExcMessage(
"This restriction matrix has not been computed yet!"));
497 return this->restriction[refinement_case - 1][child];
502 template <
int dim,
int spacedim>
506 const unsigned int codim)
const
525 if (this->degree < fe_bubbles_other->degree)
527 else if (this->degree == fe_bubbles_other->degree)
535 if (fe_nothing->is_dominating())
550 #include "fe_q_bubbles.inst"