<html>
<head>
<meta http-equiv="Content-Type" content="text/html; charset=us-ascii">
<style type="text/css" style="display:none;"> P {margin-top:0;margin-bottom:0;} </style>
</head>
<body dir="ltr">
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
Hi <font size="2"><span style="font-size:11pt">Berend</span></font>,</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
I was not able to reproduce the issue you are having, but the following 1D example (and similar 2D examples) worked fine for me using the latest PETSc. Please note that DMPlexCoordinatesLoad() now takes a PetscSF object as the third argument, but the default
 behavior is unchanged.<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
/* test_periodic_io.c */</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
#include <petscdmplex.h>
<div>#include <petscsf.h></div>
<div>#include <petsclayouthdf5.h></div>
<div><br>
</div>
<div>int main(int argc, char **argv)</div>
<div>{</div>
<div>  DM                 dm;</div>
<div>  Vec                coordinates;</div>
<div>  PetscViewer        viewer;</div>
<div>  PetscViewerFormat  format = PETSC_VIEWER_HDF5_PETSC;</div>
<div>  PetscSF            sfO;</div>
<div>  PetscErrorCode     ierr;</div>
<div><br>
</div>
<div>  ierr = PetscInitialize(&argc, &argv, NULL, NULL); if (ierr) return ierr;</div>
<div>  /* Save */</div>
<div>  ierr = PetscViewerHDF5Open(PETSC_COMM_WORLD, "periodic_example.h5", FILE_MODE_WRITE, &viewer);CHKERRQ(ierr);</div>
<div>  {</div>
<div>    DM              pdm;</div>
<div>    PetscInt        dim = 1;</div>
<div>    const PetscInt  faces[1] = {4};</div>
<div>    DMBoundaryType  periodicity[] = {DM_BOUNDARY_PERIODIC};</div>
<div>    PetscInt        overlap = 1;</div>
<div><br>
</div>
<div>    ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, faces, NULL, NULL, periodicity, PETSC_TRUE, &dm);CHKERRQ(ierr);</div>
<div>    ierr = DMPlexDistribute(dm, overlap, NULL, &pdm);CHKERRQ(ierr);</div>
<div>    ierr = DMDestroy(&dm);CHKERRQ(ierr);</div>
<div>    dm = pdm;</div>
<div>    ierr = PetscObjectSetName((PetscObject)dm, "periodicDM");CHKERRQ(ierr);</div>
<div>  }</div>
<div>  ierr = DMGetCoordinates(dm, &coordinates);CHKERRQ(ierr);</div>
<div>  ierr = PetscPrintf(PETSC_COMM_WORLD, "Coordinates before saving:\n");CHKERRQ(ierr);</div>
<div>  ierr = VecView(coordinates, NULL);CHKERRQ(ierr);</div>
<div>  ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr);</div>
<div>  ierr = DMPlexTopologyView(dm, viewer);CHKERRQ(ierr);</div>
<div>  ierr = DMPlexCoordinatesView(dm, viewer);CHKERRQ(ierr);</div>
<div>  ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);</div>
<div>  ierr = DMDestroy(&dm);CHKERRQ(ierr);</div>
  ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
  /* Load */
<div>  ierr = PetscViewerHDF5Open(PETSC_COMM_WORLD, "periodic_example.h5", FILE_MODE_READ, &viewer);CHKERRQ(ierr);</div>
<div>  ierr = DMCreate(PETSC_COMM_WORLD, &dm);CHKERRQ(ierr);</div>
<div>  ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);</div>
<div>  ierr = PetscObjectSetName((PetscObject)dm, "periodicDM");CHKERRQ(ierr);</div>
<div>  ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr);</div>
<div>  ierr = DMPlexTopologyLoad(dm, viewer, &sfO);CHKERRQ(ierr);</div>
<div>  ierr = DMPlexCoordinatesLoad(dm, viewer, sfO);CHKERRQ(ierr);</div>
<div>  ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);</div>
<div>  ierr = DMGetCoordinates(dm, &coordinates);CHKERRQ(ierr);</div>
<div>  ierr = PetscPrintf(PETSC_COMM_WORLD, "Coordinates after loading:\n");CHKERRQ(ierr);</div>
<div>  ierr = VecView(coordinates, NULL);CHKERRQ(ierr);</div>
<div>  ierr = PetscSFDestroy(&sfO);CHKERRQ(ierr);</div>
<div>  ierr = DMDestroy(&dm);CHKERRQ(ierr);</div>
<div>  ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);</div>
<div>  ierr = PetscFinalize();</div>
<div>  return ierr;</div>
}<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
mpiexec -n 2 ./test_periodic_io</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
Coordinates before saving:
<div>Vec Object: coordinates 2 MPI processes</div>
<div>  type: mpi</div>
<div>Process [0]</div>
<div>0.</div>
<div>Process [1]</div>
<div>0.25</div>
<div>0.5</div>
<div>0.75</div>
<div>Coordinates after loading:</div>
<div>Vec Object: vertices 2 MPI processes</div>
<div>  type: mpi</div>
<div>Process [0]</div>
<div>0.</div>
<div>0.25</div>
<div>0.5</div>
<div>0.75</div>
<span>Process [1]</span><br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
I would also like to note that, with the latest update, we can optionally load coordinates directly on the distributed dm as (using your notation):</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
  /* Distribute dm */<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
  ...<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<font size="2"><span style="font-size:11pt">  PetscSFCompose(sfO, sfDist, &sf);</span></font><br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
  DMPlexCoordinatesLoad(dm, viewer, sf);</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
To use this feature, we need to pass "-dm_plex_view_hdf5_storage_version 2.0.0" option when saving topology/coordinates.<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<br>
</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
Thanks,</div>
<div style="font-family: Calibri, Arial, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
Koki<br>
</div>
<div id="appendonsend"></div>
<hr style="display:inline-block;width:98%" tabindex="-1">
<div id="divRplyFwdMsg" dir="ltr"><font face="Calibri, sans-serif" style="font-size:11pt" color="#000000"><b>From:</b> Berend van Wachem <berend.vanwachem@ovgu.de><br>
<b>Sent:</b> Wednesday, November 17, 2021 3:16 PM<br>
<b>To:</b> Hapla Vaclav <vaclav.hapla@erdw.ethz.ch>; PETSc users list <petsc-users@mcs.anl.gov>; Lawrence Mitchell <wence@gmx.li>; Sagiyama, Koki <k.sagiyama@imperial.ac.uk><br>
<b>Subject:</b> Re: [petsc-users] DMView and DMLoad</font>
<div> </div>
</div>
<div class="BodyFragment"><font size="2"><span style="font-size:11pt;">
<div class="PlainText"><br>
*******************<br>
This email originates from outside Imperial. Do not click on links and attachments unless you recognise the sender.
<br>
If you trust the sender, add them to your safe senders list <a href="https://spam.ic.ac.uk/SpamConsole/Senders.aspx">
https://spam.ic.ac.uk/SpamConsole/Senders.aspx</a> to disable email stamping for this address.<br>
*******************<br>
Dear Vaclav, Lawrence, Koki,<br>
<br>
Thanks for your help! Following your advice and following your example <br>
(<a href="https://petsc.org/main/docs/manual/dmplex/#saving-and-loading-data-with-hdf5">https://petsc.org/main/docs/manual/dmplex/#saving-and-loading-data-with-hdf5</a>)
<br>
we are able to save and load the DM with a wrapped Vector in h5 format <br>
(PETSC_VIEWER_HDF5_PETSC) successfully.<br>
<br>
For saving, we use something similar to:<br>
<br>
     DMPlexTopologyView(dm, viewer);<br>
     DMClone(dm, &sdm);<br>
     ...<br>
     DMPlexSectionView(dm, viewer, sdm);<br>
     DMGetLocalVector(sdm, &vec);<br>
     ...<br>
     DMPlexLocalVectorView(dm, viewer, sdm, vec);<br>
<br>
and for loading:<br>
<br>
     DMCreate(PETSC_COMM_WORLD, &dm);<br>
     DMSetType(dm, DMPLEX);<br>
         ...<br>
       PetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_PETSC);<br>
     DMPlexTopologyLoad(dm, viewer, &sfO);<br>
     DMPlexLabelsLoad(dm, viewer);<br>
     DMPlexCoordinatesLoad(dm, viewer);<br>
     PetscViewerPopFormat(viewer);<br>
     ...<br>
     PetscSFCompose(sfO, sfDist, &sf);<br>
     ...<br>
     DMClone(dm, &sdm);<br>
     DMPlexSectionLoad(dm, viewer, sdm, sf, &globalDataSF, &localDataSF);<br>
     DMGetLocalVector(sdm, &vec);<br>
     ...<br>
     DMPlexLocalVectorLoad(dm, viewer, sdm, localDataSF, vec);<br>
<br>
<br>
This works fine for non-periodic DMs but for periodic cases the line:<br>
<br>
     DMPlexCoordinatesLoad(dm, H5Viewer);<br>
<br>
delivers the error message: invalid argument and the number of loaded <br>
coordinates does not match the number of vertices.<br>
<br>
Is this a known shortcoming, or have we forgotten something to load <br>
periodic DMs?<br>
<br>
Best regards,<br>
<br>
Berend.<br>
<br>
<br>
<br>
On 9/22/21 20:59, Hapla Vaclav wrote:<br>
> To avoid confusions here, Berend seems to be specifically demanding XDMF <br>
> (PETSC_VIEWER_HDF5_XDMF). The stuff we are now working on is parallel <br>
> checkpointing in our own HDF5 format (PETSC_VIEWER_HDF5_PETSC), I will <br>
> make a series of MRs on this topic in the following days.<br>
> <br>
> For XDMF, we are specifically missing the ability to write/load DMLabels <br>
> properly. XDMF uses specific cell-local numbering for faces for <br>
> specification of face sets, and face-local numbering for specification <br>
> of edge sets, which is not great wrt DMPlex design. And ParaView doesn't <br>
> show any of these properly so it's hard to debug. Matt, we should talk <br>
> about this soon.<br>
> <br>
> Berend, for now, could you just load the mesh initially from XDMF and <br>
> then use our PETSC_VIEWER_HDF5_PETSC format for subsequent saving/loading?<br>
> <br>
> Thanks,<br>
> <br>
> Vaclav<br>
> <br>
>> On 17 Sep 2021, at 15:46, Lawrence Mitchell <wence@gmx.li <br>
>> <<a href="mailto:wence@gmx.li">mailto:wence@gmx.li</a>>> wrote:<br>
>><br>
>> Hi Berend,<br>
>><br>
>>> On 14 Sep 2021, at 12:23, Matthew Knepley <knepley@gmail.com <br>
>>> <<a href="mailto:knepley@gmail.com">mailto:knepley@gmail.com</a>>> wrote:<br>
>>><br>
>>> On Tue, Sep 14, 2021 at 5:15 AM Berend van Wachem <br>
>>> <berend.vanwachem@ovgu.de <<a href="mailto:berend.vanwachem@ovgu.de">mailto:berend.vanwachem@ovgu.de</a>>> wrote:<br>
>>> Dear PETSc-team,<br>
>>><br>
>>> We are trying to save and load distributed DMPlex and its associated<br>
>>> physical fields (created with DMCreateGlobalVector)  (Uvelocity,<br>
>>> VVelocity,  ...) in HDF5_XDMF format. To achieve this, we do the <br>
>>> following:<br>
>>><br>
>>> 1) save in the same xdmf.h5 file:<br>
>>> DMView( DM         , H5_XDMF_Viewer );<br>
>>> VecView( UVelocity, H5_XDMF_Viewer );<br>
>>><br>
>>> 2) load the dm:<br>
>>> DMPlexCreateFromfile(PETSC_COMM_WORLD, Filename, PETSC_TRUE, DM);<br>
>>><br>
>>> 3) load the physical field:<br>
>>> VecLoad( UVelocity, H5_XDMF_Viewer );<br>
>>><br>
>>> There are no errors in the execution, but the loaded DM is distributed<br>
>>> differently to the original one, which results in the incorrect<br>
>>> placement of the values of the physical fields (UVelocity etc.) in the<br>
>>> domain.<br>
>>><br>
>>> This approach is used to restart the simulation with the last saved DM.<br>
>>> Is there something we are missing, or there exists alternative routes to<br>
>>> this goal? Can we somehow get the IS of the redistribution, so we can<br>
>>> re-distribute the vector data as well?<br>
>>><br>
>>> Many thanks, best regards,<br>
>>><br>
>>> Hi Berend,<br>
>>><br>
>>> We are in the midst of rewriting this. We want to support saving <br>
>>> multiple meshes, with fields attached to each,<br>
>>> and preserving the discretization (section) information, and allowing <br>
>>> us to load up on a different number of<br>
>>> processes. We plan to be done by October. Vaclav and I are doing this <br>
>>> in collaboration with Koki Sagiyama,<br>
>>> David Ham, and Lawrence Mitchell from the Firedrake team.<br>
>><br>
>> The core load/save cycle functionality is now in PETSc main. So if <br>
>> you're using main rather than a release, you can get access to it now. <br>
>> This section of the manual shows an example of how to do <br>
>> thingshttps://petsc.org/main/docs/manual/dmplex/#saving-and-loading-data-with-hdf5
<br>
>> <<a href="https://petsc.org/main/docs/manual/dmplex/#saving-and-loading-data-with-hdf5">https://petsc.org/main/docs/manual/dmplex/#saving-and-loading-data-with-hdf5</a>><br>
>><br>
>> Let us know if things aren't clear!<br>
>><br>
>> Thanks,<br>
>><br>
>> Lawrence<br>
> <br>
</div>
</span></font></div>
</body>
</html>