* Save the refinement information from the coarse mesh into the given
* file. This file needs to be reachable from all nodes in the computation
* on a shared network file system. See the SolutionTransfer class
- * on how to store solution vectors into this file.
+ * on how to store solution vectors into this file. Additional cell-based data can be saved
+ * using register_data_attach().
*/
void save(const char *filename) const;
/**
* Load the refinement information saved with save() back in. The mesh
- * must contain the same coarse mesh that was used in save(). You do not
+ * must contain the same coarse mesh that was used in save() before calling
+ * this function. You do not
* need to load with the same number of MPI processes that you saved
* with. Rather, if a mesh is loaded with a different number of MPI
* processes than used at the time of saving, the mesh is repartitioned
- * appropriately.
+ * appropriately. Cell-based data that was saved with register_data_attach()
+ * can be read in with notify_ready_to_unpack() after calling load().
*/
void load(const char *filename);
if (attached_data_size>0)
real_data_size = attached_data_size+sizeof(CellStatus);
+ Assert(this->n_cells()>0, ExcMessage("Can not save() an empty Triangulation."));
+
if (my_subdomain==0)
{
std::string fname=std::string(filename)+".info";
std::ofstream f(fname.c_str());
- f << Utilities::MPI::n_mpi_processes (mpi_communicator) << " "
+ f << "version nproc attached_bytes n_attached_objs n_coarse_cells" << std::endl
+ << 2 << " "
+ << Utilities::MPI::n_mpi_processes (mpi_communicator) << " "
<< real_data_size << " "
- << attached_data_pack_callbacks.size() << std::endl;
+ << attached_data_pack_callbacks.size() << " "
+ << this->n_cells(0)
+ << std::endl;
}
if (attached_data_size>0)
Triangulation<dim,spacedim>::
load(const char *filename)
{
+ Assert(this->n_cells()>0, ExcMessage("load() only works if Triangulation already contains the same coarse mesh!"));
+ Assert(this->n_levels()==1, ExcMessage("Triangulation may only contain coarse cells when calling load()."));
+
+
if (parallel_ghost != 0)
{
dealii::internal::p4est::functions<dim>::ghost_destroy (parallel_ghost);
dealii::internal::p4est::functions<dim>::destroy (parallel_forest);
parallel_forest = 0;
dealii::internal::p4est::functions<dim>::connectivity_destroy (connectivity);
- connectivity=0;
+ connectivity = 0;
- unsigned int numcpus, attached_size, attached_count;
+ unsigned int version, numcpus, attached_size, attached_count, n_coarse_cells;
{
std::string fname=std::string(filename)+".info";
std::ifstream f(fname.c_str());
- f >> numcpus >> attached_size >> attached_count;
+ std::string firstline;
+ getline(f, firstline); //skip first line
+ f >> version >> numcpus >> attached_size >> attached_count >> n_coarse_cells;
}
+ Assert(version == 2, ExcMessage("Incompatible version found in .info file."));
+ Assert(this->n_cells(0) == n_coarse_cells, ExcMessage("Number of coarse cells differ!"));
+
attached_data_size = 0;
n_attached_datas = 0;
n_attached_deserialize = attached_count;