[petsc-users] Generating xdmf from h5 file.

Matthew Knepley knepley at gmail.com
Fri Sep 26 14:52:52 CDT 2014


On Thu, Sep 25, 2014 at 10:29 PM, subramanya sadasiva <potaman at outlook.com>
wrote:

> Hi Matt,
> Sorry about that,
> I changed
> if 'time' in h5:
>     time        = np.array(h5['time']).flatten()
>   else:
>     time        = np.empty(1)
>
> The code now fails in the writeSpaceGridHeader  function.  with the error,
>
>
> Traceback (most recent call last):
>
>   File
> "/Users/ssadasiv/software/petsc/bin/pythonscripts/petsc_gen_xdmf.py", line
> 232, in <module>
>
>     generateXdmf(sys.argv[1])
>
>   File
> "/Users/ssadasiv/software/petsc/bin/pythonscripts/petsc_gen_xdmf.py", line
> 227, in generateXdmf
>
>     Xdmf(xdmfFilename).write(hdfFilename, topoPath, numCells, numCorners,
> cellDim, geomPath, numVertices, spaceDim, time, vfields, cfields)
>
>   File
> "/Users/ssadasiv/software/petsc/bin/pythonscripts/petsc_gen_xdmf.py", line
> 180, in write
>
>     self.writeSpaceGridHeader(fp, numCells, numCorners, cellDim, spaceDim)
>
>   File
> "/Users/ssadasiv/software/petsc/bin/pythonscripts/petsc_gen_xdmf.py", line
> 64, in writeSpaceGridHeader
>     print self.cellMap[cellDim][numCorners]
>
> The error is due to the fact that numCorners is set to be 1 , while
> celldim=2. cellMap has the following elements.
>
> {1: {1: 'Polyvertex', 2: 'Polyline'}, 2: {3: 'Triangle', 4:
> 'Quadrilateral'}, 3: {8: 'Hexahedron', 4: 'Tetrahedron'}}
>
>
> I  also  tried
>
> ./ex12 -dm_view vtk:my.vtk:vtk_vtu . This doesn't seem to do anything.  Is
> there any specific option I need to build petsc with to get vtk output? My
> current build has hdf5 and netcdf enabled.
>
I made the change in a branch, and merged to 'next', so it will work if you
pull. Here is a test from ex12 (testnum 39 in builder.py):

  ./ex12 -run_type full -refinement_limit 0.015625 -interpolate 1
-petscspace_order 2 -pc_type gamg -ksp_rtol 1.0e-10 -ksp_monitor_short
-ksp_converged_reason -snes_monitor_short -snes_converged_reason -dm_view
hdf5:sol.h5 -snes_view_solution hdf5:sol.h5::append

which makes sol.h5. Then I run

  ./bin/pythonscripts/petsc_gen_xdmf.py sol.h5

which makes sol.xmf. I load it up in Paraview and it makes the attached
picture.

  Thanks,

     Matt


> Thanks,
>
> Subramanya
>
> Date: Wed, 24 Sep 2014 17:36:52 -0500
> Subject: Re: [petsc-users] Generating xdmf from h5 file.
> From: knepley at gmail.com
> To: potaman at outlook.com; petsc-maint at mcs.anl.gov; petsc-users at mcs.anl.gov
>
> On Wed, Sep 24, 2014 at 5:29 PM, subramanya sadasiva <potaman at outlook.com>
> wrote:
>
> Hi Matt,
> That did not help.
>
>
> That's not enough description to fix anything, and fixing it will require
> programming.
>
>
> Is there any other way to output the mesh to something that paraview can
> view?  I tried outputting the file to a vtk file using
> ex12 -dm_view vtk:my.vtk:ascii_vtk
>
> which, I saw in another post on the forums, but that did not give me any
> output.
>
>
> This is mixing two different things. PETSc has a diagnostic ASCII vtk
> output, so the type would be ascii, not vtk,
> and format ascii_vtk . It also has a production VTU output, which is type
> vtk with format vtk_vtu.
>
>   Thanks,
>
>      Matt
>
>
>
> Subramanya
>
> ------------------------------
> Date: Wed, 24 Sep 2014 17:19:51 -0500
> Subject: Re: [petsc-users] Generating xdmf from h5 file.
> From: knepley at gmail.com
> To: potaman at outlook.com
> CC: petsc-users at mcs.anl.gov
>
> On Wed, Sep 24, 2014 at 5:08 PM, subramanya sadasiva <potaman at outlook.com>
> wrote:
>
> Hi,
> i was trying to use petsc_gen_xdmf.py to convert a h5 file to a xdmf file.
> The h5 file was generated by snes/ex12 which was run as,
>
> ex12 -dm_view hdf5:my.h5
>
> When I do,
> petsc_gen_xdmf.py my.h5
>
> I get the following error,
>
>  File "/home/ssadasiv/software/petsc/bin/pythonscripts/petsc_gen_xdmf.py",
> line 220, in <module>
>     generateXdmf(sys.argv[1])
>   File
> "/home/ssadasiv/software/petsc/bin/pythonscripts/petsc_gen_xdmf.py", line
> 208, in generateXdmf
>     time        = np.array(h5['time']).flatten()
>   File "/usr/lib/python2.7/dist-packages/h5py/_hl/group.py", line 153, in
> __getitem__
>     oid = h5o.open(self.id, self._e(name), lapl=self._lapl)
>   File "h5o.pyx", line 173, in h5py.h5o.open (h5py/h5o.c:3403)
> KeyError: "unable to open object (Symbol table: Can't open object)"
>
> I am not sure if the error is on my end. This is on Ubuntu 14.04 with the
> serial version of hdf5. I built petsc with --download-hdf5, is it necessary
> to use the same version of hdf5 to generate the xdmf file?
>
>
> That code is alpha, and mainly built for me to experiment with an
> application here, so it is not user-friendly. In your
> HDF5 file, there is no 'time' since you are not running a TS. This access
> to h5['time'] should just be protected, and
> an empty array should be put in if its not there.
>
>   Matt
>
>
> Thanks
> Subramanya
>
>
>
>
> --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> -- Norbert Wiener
>
>
>
>
> --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> -- Norbert Wiener
>



-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140926/ebcc4632/attachment-0001.html>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: sol.pdf
Type: application/pdf
Size: 1331865 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140926/ebcc4632/attachment-0001.pdf>


More information about the petsc-users mailing list