19#ifdef DEAL_II_WITH_P4EST
58 template <
typename value_type>
61 const unsigned int dofs_per_cell)
63 for (
const auto &
values : dof_values)
69 const std::size_t bytes_per_entry =
sizeof(value_type) * dofs_per_cell;
71 std::vector<char> buffer(dof_values.size() * bytes_per_entry);
72 for (
unsigned int i = 0; i < dof_values.size(); ++i)
73 std::memcpy(&buffer[i * bytes_per_entry],
85 template <
typename value_type>
86 std::vector<Vector<value_type>>
88 const boost::iterator_range<std::vector<char>::const_iterator> &data_range,
89 const unsigned int dofs_per_cell)
91 const std::size_t bytes_per_entry =
sizeof(value_type) * dofs_per_cell;
92 const unsigned int n_elements = data_range.size() / bytes_per_entry;
96 std::vector<Vector<value_type>> unpacked_data;
97 unpacked_data.reserve(n_elements);
98 for (
unsigned int i = 0; i < n_elements; ++i)
101 std::memcpy(&dof_values(0),
102 &(*std::next(data_range.begin(), i * bytes_per_entry)),
104 unpacked_data.emplace_back(std::move(dof_values));
107 return unpacked_data;
115 namespace distributed
117 template <
int dim,
typename VectorType,
typename DoFHandlerType>
119 const DoFHandlerType &dof)
120 : dof_handler(&dof, typeid(*this).name())
126 DoFHandlerType::space_dimension
> *>(
129 "parallel::distributed::SolutionTransfer requires a parallel::distributed::Triangulation object."));
134 template <
int dim,
typename VectorType,
typename DoFHandlerType>
138 const std::vector<const VectorType *> &all_in)
140 for (
unsigned int i = 0; i < all_in.size(); ++i)
141 Assert(all_in[i]->size() == dof_handler->n_dofs(),
144 input_vectors = all_in;
145 register_data_attach();
150 template <
int dim,
typename VectorType,
typename DoFHandlerType>
157 DoFHandlerType::space_dimension
> *>(
159 *
>(&dof_handler->get_triangulation())));
162 handle = tria->register_data_attach(
165 cell_iterator &cell_,
167 CellStatus status) {
return this->pack_callback(cell_, status); },
168 dof_handler->has_hp_capabilities());
173 template <
int dim,
typename VectorType,
typename DoFHandlerType>
178 std::vector<const VectorType *> all_in(1, &in);
179 prepare_for_coarsening_and_refinement(all_in);
184 template <
int dim,
typename VectorType,
typename DoFHandlerType>
189 std::vector<const VectorType *> all_in(1, &in);
190 prepare_for_serialization(all_in);
195 template <
int dim,
typename VectorType,
typename DoFHandlerType>
200 prepare_for_coarsening_and_refinement(all_in);
205 template <
int dim,
typename VectorType,
typename DoFHandlerType>
210 std::vector<VectorType *> all_in(1, &in);
216 template <
int dim,
typename VectorType,
typename DoFHandlerType>
219 std::vector<VectorType *> &all_in)
221 register_data_attach();
224 input_vectors.resize(all_in.size());
230 template <
int dim,
typename VectorType,
typename DoFHandlerType>
233 std::vector<VectorType *> &all_out)
235 Assert(input_vectors.size() == all_out.size(),
237 for (
unsigned int i = 0; i < all_out.size(); ++i)
238 Assert(all_out[i]->size() == dof_handler->n_dofs(),
244 DoFHandlerType::space_dimension
> *>(
246 *
>(&dof_handler->get_triangulation())));
249 tria->notify_ready_to_unpack(
253 cell_iterator &cell_,
256 const boost::iterator_range<std::vector<char>::const_iterator>
258 this->unpack_callback(cell_, status, data_range, all_out);
261 for (
typename std::vector<VectorType *>::iterator it = all_out.begin();
266 input_vectors.clear();
271 template <
int dim,
typename VectorType,
typename DoFHandlerType>
276 std::vector<VectorType *> all_out(1, &out);
282 template <
int dim,
typename VectorType,
typename DoFHandlerType>
286 cell_iterator &cell_,
288 DoFHandlerType::space_dimension>::CellStatus
291 typename DoFHandlerType::cell_iterator cell(*cell_, dof_handler);
294 std::vector<::Vector<typename VectorType::value_type>> dof_values(
295 input_vectors.size());
297 unsigned int fe_index = 0;
298 if (dof_handler->has_hp_capabilities())
304 DoFHandlerType::space_dimension>::CELL_PERSIST:
307 DoFHandlerType::space_dimension>::CELL_REFINE:
309 fe_index = cell->future_fe_index();
315 DoFHandlerType::space_dimension>::CELL_COARSEN:
321 for (
const auto &child : cell->child_iterators())
322 Assert(child->is_active() && child->coarsen_flag_set(),
324 dim>::ExcInconsistentCoarseningFlags());
327 fe_index = ::internal::hp::DoFHandlerImplementation::
328 dominated_future_fe_on_children<
330 DoFHandlerType::space_dimension>(cell);
340 const unsigned int dofs_per_cell =
341 dof_handler->get_fe(fe_index).n_dofs_per_cell();
343 if (dofs_per_cell == 0)
344 return std::vector<char>();
346 auto it_input = input_vectors.cbegin();
347 auto it_output = dof_values.begin();
348 for (; it_input != input_vectors.cend(); ++it_input, ++it_output)
350 it_output->reinit(dofs_per_cell);
351 cell->get_interpolated_dof_values(*(*it_input), *it_output, fe_index);
354 return pack_dof_values<typename VectorType::value_type>(dof_values,
360 template <
int dim,
typename VectorType,
typename DoFHandlerType>
364 cell_iterator &cell_,
366 DoFHandlerType::space_dimension>::CellStatus
368 const boost::iterator_range<std::vector<char>::const_iterator>
370 std::vector<VectorType *> &all_out)
372 typename DoFHandlerType::cell_iterator cell(*cell_, dof_handler);
374 unsigned int fe_index = 0;
375 if (dof_handler->has_hp_capabilities())
381 DoFHandlerType::space_dimension>::CELL_PERSIST:
384 DoFHandlerType::space_dimension>::CELL_COARSEN:
386 fe_index = cell->active_fe_index();
392 DoFHandlerType::space_dimension>::CELL_REFINE:
399 fe_index = cell->child(0)->active_fe_index();
400 for (
unsigned int child_index = 1;
401 child_index < cell->n_children();
403 Assert(cell->child(child_index)->active_fe_index() ==
415 const unsigned int dofs_per_cell =
416 dof_handler->get_fe(fe_index).n_dofs_per_cell();
418 if (dofs_per_cell == 0)
421 const std::vector<::Vector<typename VectorType::value_type>>
423 unpack_dof_values<typename VectorType::value_type>(data_range,
431 for (
auto it_dof_values = dof_values.begin();
432 it_dof_values != dof_values.end();
435 dofs_per_cell == it_dof_values->size(),
437 "The transferred data was packed with a different number of dofs than the "
438 "currently registered FE object assigned to the DoFHandler has."));
441 auto it_input = dof_values.cbegin();
442 auto it_output = all_out.begin();
443 for (; it_input != dof_values.cend(); ++it_input, ++it_output)
444 cell->set_dof_values_by_interpolation(*it_input,
453# include "solution_transfer.inst"
std::vector< char > pack_callback(const typename Triangulation< dim, DoFHandlerType::space_dimension >::cell_iterator &cell, const typename Triangulation< dim, DoFHandlerType::space_dimension >::CellStatus status)
SmartPointer< const DoFHandlerType, SolutionTransfer< dim, VectorType, DoFHandlerType > > dof_handler
void unpack_callback(const typename Triangulation< dim, DoFHandlerType::space_dimension >::cell_iterator &cell, const typename Triangulation< dim, DoFHandlerType::space_dimension >::CellStatus status, const boost::iterator_range< std::vector< char >::const_iterator > &data_range, std::vector< VectorType * > &all_out)
void prepare_for_serialization(const VectorType &in)
void register_data_attach()
void deserialize(VectorType &in)
void prepare_for_coarsening_and_refinement(const std::vector< const VectorType * > &all_in)
SolutionTransfer(const DoFHandlerType &dof)
void interpolate(std::vector< VectorType * > &all_out)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
static const unsigned int invalid_unsigned_int
TriangulationBase< dim, spacedim > Triangulation