--- /dev/null
+New: The method parallel::distributed::Triagnulation::load()
+can now also accept a p4est/p8est forest, which can be queried
+from an existing triangulation via parallel::distributed::Triagnulation::get_p4est().
+<br>
+(Marc Fehling, Peter Munch, 2021/02/09)
+
p4est_init_t init_fn,
void * user_pointer);
+ static types<2>::forest *(©_forest)(types<2>::forest *input,
+ int copy_data);
+
static void (&destroy)(types<2>::forest *p4est);
static void (&refine)(types<2>::forest *p4est,
p8est_init_t init_fn,
void * user_pointer);
+ static types<3>::forest *(©_forest)(types<3>::forest *input,
+ int copy_data);
+
static void (&destroy)(types<3>::forest *p8est);
static void (&refine)(types<3>::forest *p8est,
load(const std::string &filename,
const bool autopartition = true) override;
+ /**
+ * Load the refinement information from a given parallel forest. This
+ * forest might be obtained from the function call to
+ * parallel::distributed::Triangulation::get_p4est();
+ */
+ void
+ load(const typename dealii::internal::p4est::types<dim>::forest *forest);
+
/**
* Return a permutation vector for the order the coarse cells are handed
* off to p4est. For example the value of the $i$th element in this
}
+ namespace internal
+ {
+ template <int dim, int spacedim>
+ void
+ load_forest(const parallel::distributed::Triangulation<dim, spacedim> &in,
+ parallel::distributed::Triangulation<dim, spacedim> & out)
+ {
+#ifndef DEAL_II_WITH_P4EST
+ (void)in;
+ (void)out;
+#else
+ out.load(in.get_p4est());
+#endif
+ }
+
+ template <int spacedim>
+ void
+ load_forest(const parallel::distributed::Triangulation<1, spacedim> &in,
+ parallel::distributed::Triangulation<1, spacedim> & out)
+ {
+ (void)in;
+ (void)out;
+ }
+ } // namespace internal
template <int dim, int spacedim>
std::vector<std::shared_ptr<Triangulation<dim, spacedim>>>
// create coarse meshes
for (unsigned int l = max_level; l > 0; --l)
{
- // store finer mesh to file
- dynamic_cast<
- const parallel::distributed::Triangulation<dim, spacedim> *>(
- coarse_grid_triangulations[l].get())
- ->save("mesh");
- MPI_Barrier(fine_triangulation->get_communicator());
-
// create empty triangulation
auto new_tria =
std::make_shared<parallel::distributed::Triangulation<dim, spacedim>>(
new_tria->set_manifold(i, fine_triangulation->get_manifold(i));
// create refinement hierarchy (by loading stored mesh)
- new_tria->load("mesh", false);
- MPI_Barrier(fine_triangulation->get_communicator());
+ internal::load_forest(
+ *dynamic_cast<
+ const parallel::distributed::Triangulation<dim, spacedim> *>(
+ coarse_grid_triangulations[l].get()),
+ *new_tria);
// coarsen mesh
new_tria->coarsen_global();
p4est_init_t init_fn,
void * user_pointer) = p4est_new_ext;
+ types<2>::forest *(&functions<2>::copy_forest)(types<2>::forest *input,
+ int copy_data) = p4est_copy;
+
void (&functions<2>::destroy)(types<2>::forest *p4est) = p4est_destroy;
void (&functions<2>::refine)(types<2>::forest *p4est,
p8est_init_t init_fn,
void * user_pointer) = p8est_new_ext;
+ types<3>::forest *(&functions<3>::copy_forest)(types<3>::forest *input,
+ int copy_data) = p8est_copy;
+
void (&functions<3>::destroy)(types<3>::forest *p8est) = p8est_destroy;
void (&functions<3>::refine)(types<3>::forest *p8est,
}
+ template <int dim, int spacedim>
+ void
+ Triangulation<dim, spacedim>::load(
+ const typename dealii::internal::p4est::types<dim>::forest *forest)
+ {
+ if (parallel_ghost != nullptr)
+ {
+ dealii::internal::p4est::functions<dim>::ghost_destroy(
+ parallel_ghost);
+ parallel_ghost = nullptr;
+ }
+ dealii::internal::p4est::functions<dim>::destroy(parallel_forest);
+ parallel_forest = nullptr;
+
+ // note: we can keep the connectivity, since the coarse grid does not
+ // change
+
+ // create deep copy
+ typename dealii::internal::p4est::types<dim>::forest *temp =
+ const_cast<typename dealii::internal::p4est::types<dim>::forest *>(
+ forest);
+ parallel_forest =
+ dealii::internal::p4est::functions<dim>::copy_forest(temp, false);
+ parallel_forest->user_pointer = this;
+
+ try
+ {
+ copy_local_forest_to_triangulation();
+ }
+ catch (const typename Triangulation<dim>::DistortedCellList &)
+ {
+ // the underlying
+ // triangulation should not
+ // be checking for
+ // distorted cells
+ Assert(false, ExcInternalError());
+ }
+
+ this->update_periodic_face_map();
+ this->update_number_cache();
+ }
+
+
template <int dim, int spacedim>
unsigned int