45template <
int dim,
typename VectorType,
typename DoFHandlerType>
47 const DoFHandlerType &dof)
48 : dof_handler(&dof, typeid(*this).name())
53 DoFHandlerType::dimension,
54 DoFHandlerType::space_dimension
> *>(
56 ExcMessage(
"You are calling the ::SolutionTransfer class "
57 "with a DoF handler that is built on a "
58 "parallel::distributed::Triangulation. This will not "
59 "work for parallel computations. You probably want to "
60 "use the parallel::distributed::SolutionTransfer class."));
65template <
int dim,
typename VectorType,
typename DoFHandlerType>
73template <
int dim,
typename VectorType,
typename DoFHandlerType>
77 indices_on_cell.clear();
78 dof_values_on_cell.clear();
86template <
int dim,
typename VectorType,
typename DoFHandlerType>
90 Assert(prepared_for != pure_refinement, ExcAlreadyPrepForRef());
91 Assert(prepared_for != coarsening_and_refinement,
92 ExcAlreadyPrepForCoarseAndRef());
104 const internal::parallel::shared::
105 TemporarilyRestoreSubdomainIds<dim, DoFHandlerType::space_dimension>
106 subdomain_modifier(dof_handler->get_triangulation());
109 dof_handler->get_triangulation().n_active_cells();
110 n_dofs_old = dof_handler->n_dofs();
114 .
swap(indices_on_cell);
116 for (
const auto &cell : dof_handler->active_cell_iterators())
118 const unsigned int i = cell->active_cell_index();
119 indices_on_cell[i].resize(cell->get_fe().n_dofs_per_cell());
126 cell->get_dof_indices(indices_on_cell[i]);
127 cell_map[std::make_pair(cell->level(), cell->index())] =
130 prepared_for = pure_refinement;
135template <
int dim,
typename VectorType,
typename DoFHandlerType>
138 const VectorType &in,
139 VectorType & out)
const
141 Assert(prepared_for == pure_refinement, ExcNotPrepared());
143 Assert(out.size() == dof_handler->n_dofs(),
146 ExcMessage(
"Vectors cannot be used as input and output"
147 " at the same time!"));
157 const internal::parallel::shared::
158 TemporarilyRestoreSubdomainIds<dim, DoFHandlerType::space_dimension>
159 subdomain_modifier(dof_handler->get_triangulation());
163 typename std::map<std::pair<unsigned int, unsigned int>,
165 cell_map_end = cell_map.end();
167 for (
const auto &cell : dof_handler->cell_iterators())
170 cell_map.find(std::make_pair(cell->level(), cell->index()));
172 if (pointerstruct != cell_map_end)
180 const unsigned int this_fe_index =
182 const unsigned int dofs_per_cell =
183 cell->get_dof_handler().get_fe(this_fe_index).n_dofs_per_cell();
184 local_values.
reinit(dofs_per_cell,
true);
189 Assert(dofs_per_cell == (*pointerstruct->second.indices_ptr).size(),
191 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
193 in, (*pointerstruct->second.indices_ptr)[i]);
194 cell->set_dof_values_by_interpolation(local_values,
218 template <
int dim,
int spacedim>
223 if (dof.has_hp_capabilities() ==
false)
226 const ::hp::FECollection<dim, spacedim> &fe = dof.get_fe_collection();
227 matrices.reinit(fe.size(), fe.size());
228 for (
unsigned int i = 0; i < fe.size(); ++i)
229 for (
unsigned int j = 0; j < fe.size(); ++j)
232 matrices(i, j).reinit(fe[i].n_dofs_per_cell(),
233 fe[j].n_dofs_per_cell());
244 fe[i].get_interpolation_matrix(fe[j], matrices(i, j));
247 ExcInterpolationNotImplemented &)
249 matrices(i, j).reinit(0, 0);
255 template <
int dim,
int spacedim>
258 std::vector<std::vector<bool>> &)
261 template <
int dim,
int spacedim>
264 std::vector<std::vector<bool>> &restriction_is_additive)
266 restriction_is_additive.resize(fe.size());
267 for (
unsigned int f = 0; f < fe.size(); ++f)
269 restriction_is_additive[f].resize(fe[f].n_dofs_per_cell());
270 for (
unsigned int i = 0; i < fe[f].n_dofs_per_cell(); ++i)
271 restriction_is_additive[f][i] = fe[f].restriction_is_additive(i);
278template <
int dim,
typename VectorType,
typename DoFHandlerType>
283 Assert(prepared_for != pure_refinement, ExcAlreadyPrepForRef());
284 Assert(prepared_for != coarsening_and_refinement,
285 ExcAlreadyPrepForCoarseAndRef());
288 n_dofs_old = dof_handler->n_dofs();
289 const unsigned int in_size = all_in.size();
293 ExcMessage(
"The array of input vectors you pass to this "
294 "function has no elements. This is not useful."));
295 for (
unsigned int i = 0; i < in_size; ++i)
297 Assert(all_in[i].size() == n_dofs_old,
310 const internal::parallel::shared::
311 TemporarilyRestoreSubdomainIds<dim, DoFHandlerType::space_dimension>
312 subdomain_modifier(dof_handler->get_triangulation());
317 unsigned int n_cells_to_coarsen = 0;
318 unsigned int n_cells_to_stay_or_refine = 0;
319 for (
const auto &act_cell : dof_handler->active_cell_iterators())
321 if (act_cell->coarsen_flag_set())
322 ++n_cells_to_coarsen;
324 ++n_cells_to_stay_or_refine;
326 Assert((n_cells_to_coarsen + n_cells_to_stay_or_refine) ==
327 dof_handler->get_triangulation().n_active_cells(),
330 unsigned int n_coarsen_fathers = 0;
331 for (
const auto &cell : dof_handler->cell_iterators())
332 if (!cell->is_active() && cell->child(0)->coarsen_flag_set())
339 std::vector<std::vector<types::global_dof_index>>(n_cells_to_stay_or_refine)
340 .
swap(indices_on_cell);
342 std::vector<std::vector<Vector<typename VectorType::value_type>>>(
344 std::vector<Vector<typename VectorType::value_type>>(in_size))
345 .swap(dof_values_on_cell);
348 std::vector<std::vector<bool>> restriction_is_additive;
352 restriction_is_additive);
357 unsigned int n_sr = 0, n_cf = 0;
358 for (
const auto &cell : dof_handler->cell_iterators())
361 if (cell->is_active() && !cell->coarsen_flag_set())
363 const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
364 indices_on_cell[n_sr].resize(dofs_per_cell);
369 cell->get_dof_indices(indices_on_cell[n_sr]);
370 cell_map[std::make_pair(cell->level(), cell->index())] =
371 Pointerstruct(&indices_on_cell[n_sr], cell->active_fe_index());
376 else if (cell->has_children() && cell->child(0)->coarsen_flag_set())
384 std::set<unsigned int> fe_indices_children;
385 for (
const auto &child : cell->child_iterators())
387 Assert(child->is_active() && child->coarsen_flag_set(),
389 dim>::ExcInconsistentCoarseningFlags());
391 fe_indices_children.insert(child->active_fe_index());
395 const unsigned int target_fe_index =
396 dof_handler->get_fe_collection().find_dominated_fe_extended(
397 fe_indices_children, 0);
403 const unsigned int dofs_per_cell =
404 dof_handler->get_fe(target_fe_index).n_dofs_per_cell();
406 std::vector<Vector<typename VectorType::value_type>>(
408 .
swap(dof_values_on_cell[n_cf]);
416 for (
unsigned int j = 0; j < in_size; ++j)
417 cell->get_interpolated_dof_values(all_in[j],
418 dof_values_on_cell[n_cf][j],
420 cell_map[std::make_pair(cell->level(), cell->index())] =
428 prepared_for = coarsening_and_refinement;
433template <
int dim,
typename VectorType,
typename DoFHandlerType>
438 std::vector<VectorType> all_in(1, in);
439 prepare_for_coarsening_and_refinement(all_in);
444template <
int dim,
typename VectorType,
typename DoFHandlerType>
447 const std::vector<VectorType> &all_in,
448 std::vector<VectorType> & all_out)
const
450 const unsigned int size = all_in.size();
452 Assert(prepared_for == coarsening_and_refinement, ExcNotPrepared());
454 for (
unsigned int i = 0; i < size; ++i)
455 Assert(all_in[i].size() == n_dofs_old,
457 for (
unsigned int i = 0; i < all_out.size(); ++i)
458 Assert(all_out[i].size() == dof_handler->n_dofs(),
460 for (
unsigned int i = 0; i < size; ++i)
461 for (
unsigned int j = 0; j < size; ++j)
462 Assert(&all_in[i] != &all_out[j],
463 ExcMessage(
"Vectors cannot be used as input and output"
464 " at the same time!"));
475 const internal::parallel::shared::
476 TemporarilyRestoreSubdomainIds<dim, DoFHandlerType::space_dimension>
477 subdomain_modifier(dof_handler->get_triangulation());
480 std::vector<types::global_dof_index> dofs;
482 typename std::map<std::pair<unsigned int, unsigned int>,
484 cell_map_end = cell_map.end();
490 for (
const auto &cell : dof_handler->cell_iterators())
493 cell_map.find(std::make_pair(cell->level(), cell->index()));
495 if (pointerstruct != cell_map_end)
497 const std::vector<types::global_dof_index> *
const indexptr =
498 pointerstruct->second.indices_ptr;
500 const std::vector<Vector<typename VectorType::value_type>>
501 *
const valuesptr = pointerstruct->second.dof_values_ptr;
508 const unsigned int old_fe_index =
509 pointerstruct->second.active_fe_index;
513 unsigned int in_size = indexptr->size();
514 for (
unsigned int j = 0; j < size; ++j)
516 tmp.
reinit(in_size,
true);
517 for (
unsigned int i = 0; i < in_size; ++i)
522 cell->set_dof_values_by_interpolation(tmp,
533 const unsigned int dofs_per_cell =
534 cell->get_fe().n_dofs_per_cell();
535 dofs.resize(dofs_per_cell);
538 cell->get_dof_indices(dofs);
541 for (
unsigned int j = 0; j < size; ++j)
548 const unsigned int active_fe_index = cell->active_fe_index();
549 if (active_fe_index != pointerstruct->second.active_fe_index)
551 const unsigned int old_index =
552 pointerstruct->second.active_fe_index;
554 interpolation_hp(active_fe_index, old_index);
557 if (interpolation_matrix.empty())
558 tmp.
reinit(dofs_per_cell,
false);
561 tmp.
reinit(dofs_per_cell,
true);
563 interpolation_matrix.
n());
565 interpolation_matrix.
vmult(tmp, (*valuesptr)[j]);
570 data = &(*valuesptr)[j];
573 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
588template <
int dim,
typename VectorType,
typename DoFHandlerType>
591 const VectorType &in,
592 VectorType & out)
const
595 Assert(out.size() == dof_handler->n_dofs(),
598 std::vector<VectorType> all_in(1);
600 std::vector<VectorType> all_out(1);
608template <
int dim,
typename VectorType,
typename DoFHandlerType>
618 sizeof(prepared_for) +
625template <
int dim,
typename VectorType,
typename DoFHandlerType>
630 return sizeof(*this);
635#define SPLIT_INSTANTIATIONS_COUNT 4
636#ifndef SPLIT_INSTANTIATIONS_INDEX
637# define SPLIT_INSTANTIATIONS_INDEX 0
639#include "solution_transfer.inst"
void swap(BlockIndices &u, BlockIndices &v)
void vmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
void prepare_for_pure_refinement()
SolutionTransfer(const DoFHandlerType &dof)
std::size_t memory_consumption() const
SmartPointer< const DoFHandlerType, SolutionTransfer< dim, VectorType, DoFHandlerType > > dof_handler
void refine_interpolate(const VectorType &in, VectorType &out) const
void interpolate(const std::vector< VectorType > &all_in, std::vector< VectorType > &all_out) const
void prepare_for_coarsening_and_refinement(const std::vector< VectorType > &all_in)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcNoDominatedFiniteElementOnChildren()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual void reinit(const size_type N, const bool omit_zeroing_entries=false)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
unsigned int n_active_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
void extract_interpolation_matrices(const ::DoFHandler< dim, spacedim > &dof, ::Table< 2, FullMatrix< double > > &matrices)
void restriction_additive(const FiniteElement< dim, spacedim > &, std::vector< std::vector< bool > > &)
static const unsigned int invalid_unsigned_int
TriangulationBase< dim, spacedim > Triangulation
unsigned int active_fe_index
std::size_t memory_consumption() const
static VectorType::value_type get(const VectorType &V, const types::global_dof_index i)