44 template <
int dim,
typename VectorType,
typename DoFHandlerType>
46 const DoFHandlerType &dof)
47 : dof_handler(&dof, typeid(*this).name())
52 DoFHandlerType::dimension,
53 DoFHandlerType::space_dimension
> *>(
55 ExcMessage(
"You are calling the ::SolutionTransfer class "
56 "with a DoF handler that is built on a "
57 "parallel::distributed::Triangulation. This will not "
58 "work for parallel computations. You probably want to "
59 "use the parallel::distributed::SolutionTransfer class."));
64 template <
int dim,
typename VectorType,
typename DoFHandlerType>
72 template <
int dim,
typename VectorType,
typename DoFHandlerType>
76 indices_on_cell.clear();
77 dof_values_on_cell.clear();
85 template <
int dim,
typename VectorType,
typename DoFHandlerType>
89 Assert(prepared_for != pure_refinement, ExcAlreadyPrepForRef());
90 Assert(prepared_for != coarsening_and_refinement,
91 ExcAlreadyPrepForCoarseAndRef());
96 dof_handler->get_triangulation().n_active_cells();
97 n_dofs_old = dof_handler->n_dofs();
101 .
swap(indices_on_cell);
103 typename DoFHandlerType::active_cell_iterator cell =
104 dof_handler->begin_active(),
105 endc = dof_handler->end();
107 for (
unsigned int i = 0; cell != endc; ++cell, ++i)
109 indices_on_cell[i].resize(cell->get_fe().dofs_per_cell);
116 cell->get_dof_indices(indices_on_cell[i]);
117 cell_map[std::make_pair(cell->level(), cell->index())] =
120 prepared_for = pure_refinement;
125 template <
int dim,
typename VectorType,
typename DoFHandlerType>
131 Assert(prepared_for == pure_refinement, ExcNotPrepared());
133 Assert(out.size() == dof_handler->n_dofs(),
136 ExcMessage(
"Vectors cannot be used as input and output"
137 " at the same time!"));
141 typename DoFHandlerType::cell_iterator cell = dof_handler->begin(),
142 endc = dof_handler->end();
144 typename std::map<std::pair<unsigned int, unsigned int>,
146 cell_map_end = cell_map.end();
148 for (; cell != endc; ++cell)
151 cell_map.find(std::make_pair(cell->level(), cell->index()));
153 if (pointerstruct != cell_map_end)
161 const unsigned int this_fe_index =
162 pointerstruct->second.active_fe_index;
163 const unsigned int dofs_per_cell =
164 cell->get_dof_handler().get_fe(this_fe_index).dofs_per_cell;
165 local_values.
reinit(dofs_per_cell,
true);
170 Assert(dofs_per_cell == (*pointerstruct->second.indices_ptr).size(),
172 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
174 in, (*pointerstruct->second.indices_ptr)[i]);
175 cell->set_dof_values_by_interpolation(local_values,
199 template <
typename DoFHandlerType>
205 template <
int dim,
int spacedim>
208 const ::hp::DoFHandler<dim, spacedim> &dof,
211 const ::hp::FECollection<dim, spacedim> &fe = dof.get_fe_collection();
212 matrices.reinit(fe.size(), fe.size());
213 for (
unsigned int i = 0; i < fe.size(); ++i)
214 for (
unsigned int j = 0; j < fe.size(); ++j)
217 matrices(i, j).reinit(fe[i].dofs_per_cell, fe[j].dofs_per_cell);
228 fe[i].get_interpolation_matrix(fe[j], matrices(i, j));
231 ExcInterpolationNotImplemented &)
233 matrices(i, j).reinit(0, 0);
239 template <
int dim,
int spacedim>
242 std::vector<std::vector<bool>> &)
245 template <
int dim,
int spacedim>
248 std::vector<std::vector<bool>> &restriction_is_additive)
250 restriction_is_additive.resize(fe.size());
251 for (
unsigned int f = 0; f < fe.size(); ++f)
253 restriction_is_additive[f].resize(fe[f].dofs_per_cell);
254 for (
unsigned int i = 0; i < fe[f].dofs_per_cell; ++i)
255 restriction_is_additive[f][i] = fe[f].restriction_is_additive(i);
262 template <
int dim,
typename VectorType,
typename DoFHandlerType>
267 Assert(prepared_for != pure_refinement, ExcAlreadyPrepForRef());
268 Assert(prepared_for != coarsening_and_refinement,
269 ExcAlreadyPrepForCoarseAndRef());
271 const unsigned int in_size = all_in.size();
273 ExcMessage(
"The array of input vectors you pass to this "
274 "function has no elements. This is not useful."));
279 dof_handler->get_triangulation().n_active_cells();
281 n_dofs_old = dof_handler->n_dofs();
283 for (
unsigned int i = 0; i < in_size; ++i)
285 Assert(all_in[i].size() == n_dofs_old,
292 unsigned int n_cells_to_coarsen = 0;
293 unsigned int n_cells_to_stay_or_refine = 0;
294 for (
const auto &act_cell : dof_handler->active_cell_iterators())
296 if (act_cell->coarsen_flag_set())
297 ++n_cells_to_coarsen;
299 ++n_cells_to_stay_or_refine;
304 unsigned int n_coarsen_fathers = 0;
305 for (
typename DoFHandlerType::cell_iterator cell = dof_handler->begin();
306 cell != dof_handler->end();
308 if (!cell->is_active() && cell->child(0)->coarsen_flag_set())
315 std::vector<std::vector<types::global_dof_index>>(n_cells_to_stay_or_refine)
316 .
swap(indices_on_cell);
318 std::vector<std::vector<Vector<typename VectorType::value_type>>>(
320 std::vector<Vector<typename VectorType::value_type>>(in_size))
321 .swap(dof_values_on_cell);
324 std::vector<std::vector<bool>> restriction_is_additive;
328 restriction_is_additive);
333 unsigned int n_sr = 0, n_cf = 0;
334 for (
typename DoFHandlerType::cell_iterator cell = dof_handler->begin();
335 cell != dof_handler->end();
339 if (cell->is_active() && !cell->coarsen_flag_set())
341 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
342 indices_on_cell[n_sr].resize(dofs_per_cell);
347 cell->get_dof_indices(indices_on_cell[n_sr]);
348 cell_map[std::make_pair(cell->level(), cell->index())] =
349 Pointerstruct(&indices_on_cell[n_sr], cell->active_fe_index());
354 else if (cell->has_children() && cell->child(0)->coarsen_flag_set())
362 std::set<unsigned int> fe_indices_children;
363 for (
unsigned int child_index = 0; child_index < cell->n_children();
366 const auto &child = cell->child(child_index);
367 Assert(child->is_active() && child->coarsen_flag_set(),
369 dim>::ExcInconsistentCoarseningFlags());
371 fe_indices_children.insert(child->active_fe_index());
375 const unsigned int target_fe_index =
376 dof_handler->get_fe_collection().find_dominated_fe_extended(
377 fe_indices_children, 0);
380 typename ::hp::FECollection<
381 dim>::ExcNoDominatedFiniteElementAmongstChildren());
383 const unsigned int dofs_per_cell =
384 dof_handler->get_fe(target_fe_index).dofs_per_cell;
386 std::vector<Vector<typename VectorType::value_type>>(
388 .
swap(dof_values_on_cell[n_cf]);
396 for (
unsigned int j = 0; j < in_size; ++j)
397 cell->get_interpolated_dof_values(all_in[j],
398 dof_values_on_cell[n_cf][j],
400 cell_map[std::make_pair(cell->level(), cell->index())] =
408 prepared_for = coarsening_and_refinement;
413 template <
int dim,
typename VectorType,
typename DoFHandlerType>
418 std::vector<VectorType> all_in(1, in);
419 prepare_for_coarsening_and_refinement(all_in);
424 template <
int dim,
typename VectorType,
typename DoFHandlerType>
427 const std::vector<VectorType> &all_in,
428 std::vector<VectorType> & all_out)
const
430 Assert(prepared_for == coarsening_and_refinement, ExcNotPrepared());
431 const unsigned int size = all_in.size();
433 for (
unsigned int i = 0; i < size; ++i)
434 Assert(all_in[i].size() == n_dofs_old,
436 for (
unsigned int i = 0; i < all_out.size(); ++i)
437 Assert(all_out[i].size() == dof_handler->n_dofs(),
439 for (
unsigned int i = 0; i < size; ++i)
440 for (
unsigned int j = 0; j < size; ++j)
441 Assert(&all_in[i] != &all_out[j],
442 ExcMessage(
"Vectors cannot be used as input and output"
443 " at the same time!"));
446 std::vector<types::global_dof_index> dofs;
448 typename std::map<std::pair<unsigned int, unsigned int>,
450 cell_map_end = cell_map.end();
456 for (
const auto &cell : dof_handler->cell_iterators())
459 cell_map.find(std::make_pair(cell->level(), cell->index()));
461 if (pointerstruct != cell_map_end)
463 const std::vector<types::global_dof_index> *
const indexptr =
464 pointerstruct->second.indices_ptr;
466 const std::vector<Vector<typename VectorType::value_type>>
467 *
const valuesptr = pointerstruct->second.dof_values_ptr;
474 const unsigned int old_fe_index =
475 pointerstruct->second.active_fe_index;
479 unsigned int in_size = indexptr->size();
480 for (
unsigned int j = 0; j < size; ++j)
482 tmp.
reinit(in_size,
true);
483 for (
unsigned int i = 0; i < in_size; ++i)
488 cell->set_dof_values_by_interpolation(tmp,
499 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
500 dofs.resize(dofs_per_cell);
503 cell->get_dof_indices(dofs);
506 for (
unsigned int j = 0; j < size; ++j)
513 const unsigned int active_fe_index = cell->active_fe_index();
514 if (active_fe_index != pointerstruct->second.active_fe_index)
516 const unsigned int old_index =
517 pointerstruct->second.active_fe_index;
519 interpolation_hp(active_fe_index, old_index);
522 if (interpolation_matrix.
empty())
523 tmp.
reinit(dofs_per_cell,
false);
526 tmp.
reinit(dofs_per_cell,
true);
528 interpolation_matrix.
n());
530 interpolation_matrix.
vmult(tmp, (*valuesptr)[j]);
535 data = &(*valuesptr)[j];
538 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
553 template <
int dim,
typename VectorType,
typename DoFHandlerType>
560 Assert(out.size() == dof_handler->n_dofs(),
563 std::vector<VectorType> all_in(1);
565 std::vector<VectorType> all_out(1);
573 template <
int dim,
typename VectorType,
typename DoFHandlerType>
583 sizeof(prepared_for) +
590 template <
int dim,
typename VectorType,
typename DoFHandlerType>
595 return sizeof(*this);
600 #define SPLIT_INSTANTIATIONS_COUNT 4
601 #ifndef SPLIT_INSTANTIATIONS_INDEX
602 # define SPLIT_INSTANTIATIONS_INDEX 0
604 #include "solution_transfer.inst"