Error in writing and reading cell based data for restart

63 views
Skip to first unread message

Sambit Das

unread,
Apr 9, 2018, 10:26:19 AM4/9/18
to deal.II User Group
Hi all,

I have written the following failing minimal example (also attached) where I create a single element parallel distributed triangulation and try to write and read a double. 

using namespace dealii;
int main (int argc, char *argv[])
{
  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv);

  const double L=20;
  parallel::distributed::Triangulation<3> triangulation(MPI_COMM_WORLD);
  GridGenerator::hyper_cube (triangulation, -L, L);

  unsigned int offset = triangulation.register_data_attach
          (sizeof(double),
   [&](const typename dealii::parallel::distributed::Triangulation<3>::cell_iterator &cell,
       const typename dealii::parallel::distributed::Triangulation<3>::CellStatus status,
       void * data) -> void
       {

     double* dataStore = reinterpret_cast<double*>(data);
     *dataStore=0.0;
       }
   );

  std::cout<< "offset=" << offset << std::endl;
  std::string filename="triangulationChk";
  triangulation.save(filename.c_str());

  triangulation.load(filename.c_str());
  triangulation.notify_ready_to_unpack
      (offset,[&](const typename dealii::parallel::distributed::Triangulation<3>::cell_iterator &cell,
        const typename dealii::parallel::distributed::Triangulation<3>::CellStatus status,
        const void * data) -> void
        {
        }
       );
}


I get the following error:

An error occurred in line <3236> of file </home/vikramg/DFT-FE-softwares/softwareCentos/dealiiDev/dealii/source/distributed/tria.cc> in function
    void dealii::parallel::distributed::Triangulation<dim, spacedim>::notify_ready_to_unpack(unsigned int, const std::function<void (const dealii::Triangulation<dim, spacedim>::cell_iterator &, dealii::Triangulation<dim, spacedim>::CellStatus, const void *)> &) [with int dim = 3, int spacedim = 3]
The violated condition was: 
    offset < sizeof(CellStatus)+attached_data_size
Additional information: 
    invalid offset in notify_ready_to_unpack()

I am wondering if I am missing something in my implementation.

Thanks,
Sambit

minimalExampleCellQuadDataWriteRead.cc

Sambit Das

unread,
Apr 9, 2018, 12:49:40 PM4/9/18
to deal.II User Group

An error occurred in line <3236> of file </home/vikramg/DFT-FE-softwares/softwareCentos/dealiiDev/dealii/source/distributed/tria.cc> in function
    void dealii::parallel::distributed::Triangulation<dim, spacedim>::notify_ready_to_unpack(unsigned int, const std::function<void (const dealii::Triangulation<dim, spacedim>::cell_iterator &, dealii::Triangulation<dim, spacedim>::CellStatus, const void *)> &) [with int dim = 3, int spacedim = 3]
The violated condition was: 
    offset < sizeof(CellStatus)+attached_data_size
Additional information: 
    invalid offset in notify_ready_to_unpack()

I tried to debug this using ddt debugger- the reason the above condition is violated is due to attached_data_size=0, which is strange as the triangulationChk.info file mentions one object to be attached:

 triangulationChk.info 
version nproc attached_bytes n_attached_objs n_coarse_cells
2 1 12 1 1 

Wolfgang Bangerth

unread,
Apr 11, 2018, 12:56:55 AM4/11/18
to dea...@googlegroups.com

> I get the following error:
>
> /An error occurred in line <3236> of file
> </home/vikramg/DFT-FE-softwares/softwareCentos/dealiiDev/dealii/source/distributed/tria.cc>
> in function/
> / void dealii::parallel::distributed::Triangulation<dim,
> spacedim>::notify_ready_to_unpack(unsigned int, const std::function<void
> (const dealii::Triangulation<dim, spacedim>::cell_iterator &,
> dealii::Triangulation<dim, spacedim>::CellStatus, const void *)> &) [with int
> dim = 3, int spacedim = 3]/
> /The violated condition was: /
> / offset < sizeof(CellStatus)+attached_data_size/
> /Additional information: /
> / invalid offset in notify_ready_to_unpack()/
> /
> /
> I am wondering if I am missing something in my implementation.

Nothing obvious at least. I copied the issue here:
https://github.com/dealii/dealii/issues/6223

I don't recall why we reset the size to zero. The intention was probably that
we don't ever call `load()` in the same execution as we called `save()`
before. We may be calling `save()` many times, of course, in the course of
long-running programs, but call `load()` only at the beginning of a restart.
Consequently, this may simply be a bug that we have just never encountered, or
there may have been a reason to do it that way -- I can't recall.

I imagine you can work around this is you simply call `register_data_attach()`
again before calling `load()`.

Best
W.

--
------------------------------------------------------------------------
Wolfgang Bangerth email: bang...@colostate.edu
www: http://www.math.colostate.edu/~bangerth/

Sambit Das

unread,
Apr 11, 2018, 3:05:31 AM4/11/18
to deal.II User Group
Prof. Bangerth, Thank you for your reply and creating the issue.


I imagine you can work around this is you simply call `register_data_attach()`
again before calling `load()`.

I checked that calling  `register_data_attach()`again before calling `load()` 
still triggers the same assert, but calling `register_data_attach()` again just after 
calling `load()`  and before calling 'notify_ready_to_unpack` works- I 
am able to read and print the correct cell data from the file.

Best,
Sambit
Reply all
Reply to author
Forward
0 new messages