-
Notifications
You must be signed in to change notification settings - Fork 33
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Relative file paths & field data with parallel file formats #139
Comments
Thank you for your interest!
Definitely! That sounds like a useful feature which should be easy to add. I'll look into it.
This sounds like a bug, but may also be an issue in ParaView. Could you provide a minimal working example for testing? (Maybe based on this example.) |
Yes, it is relatively easy to reproduce: using WriteVTK
all_data = [
# Process 1
(
points = rand(3, 5), # 5 points on process 1
cells = [ # 2 cells on process 1
MeshCell(VTKCellTypes.VTK_TRIANGLE, [1, 4, 2]),
MeshCell(VTKCellTypes.VTK_QUAD, [2, 4, 3, 5]),
],
),
# Process 2
(
points = rand(3, 4), # 4 points on process 2
cells = [ # 1 cell on process 2
MeshCell(VTKCellTypes.VTK_QUAD, [1, 2, 3, 4]),
]
),
]
saved_files = Vector{Vector{String}}(undef, 2) # files saved by each "process"
for part = 1:2
data = all_data[part]
saved_files[part] = pvtk_grid(
"simulation", data.points, data.cells;
part = part, nparts = 2,
) do pvtk
pvtk["Pressure"] = sum(data.points; dims = 1)
pvtk["Time"] = 0.1234 # <-- just added this
end
end Now the file <?xml version="1.0" encoding="utf-8"?>
<VTKFile type="PUnstructuredGrid" version="1.0" byte_order="LittleEndian">
<PUnstructuredGrid GhostLevel="0">
<Piece Source="simulation/simulation_1.vtu"/>
<Piece Source="simulation/simulation_2.vtu"/>
<PPoints>
<PDataArray type="Float64" Name="Points" NumberOfComponents="3"/>
</PPoints>
<PPointData>
<PDataArray type="Float64" Name="Pressure" NumberOfComponents="1"/>
</PPointData>
</PUnstructuredGrid>
</VTKFile> However, the <?xml version="1.0" encoding="utf-8"?>
<VTKFile type="UnstructuredGrid" version="1.0" byte_order="LittleEndian" header_type="UInt64" compressor="vtkZLibDataCompressor">
<UnstructuredGrid>
<Piece NumberOfPoints="5" NumberOfCells="2">
<Points>
<DataArray type="Float64" Name="Points" NumberOfComponents="3" format="appended" offset="0"/>
</Points>
<Cells>
<DataArray type="Int64" Name="connectivity" NumberOfComponents="1" format="appended" offset="163"/>
<DataArray type="Int64" Name="offsets" NumberOfComponents="1" format="appended" offset="219"/>
<DataArray type="UInt8" Name="types" NumberOfComponents="1" format="appended" offset="265"/>
</Cells>
<PointData>
<DataArray type="Float64" Name="Pressure" NumberOfComponents="1" format="appended" offset="307"/>
</PointData>
</Piece>
<FieldData>
<DataArray type="Float64" Name="Time" NumberOfComponents="1" NumberOfTuples="1" format="appended" offset="390"/>
</FieldData>
</UnstructuredGrid>
<AppendedData encoding="raw">
...
</AppendedData>
</VTKFile> |
Thanks! I think I just managed to make the time information visible by ParaView. In fact the Can you try opening this modified <?xml version="1.0" encoding="utf-8"?>
<VTKFile type="PUnstructuredGrid" version="1.0" byte_order="LittleEndian">
<PUnstructuredGrid GhostLevel="0">
<Piece Source="simulation/simulation_1.vtu"/>
<Piece Source="simulation/simulation_2.vtu"/>
<PPoints>
<PDataArray type="Float64" Name="Points" NumberOfComponents="3"/>
</PPoints>
<FieldData>
<DataArray type="Float64" Name="Time" NumberOfComponents="1">
0.1234
</DataArray>
</FieldData>
<PPointData>
<PDataArray type="Float64" Name="Pressure" NumberOfComponents="1"/>
</PPointData>
</PUnstructuredGrid>
</VTKFile> |
Yes, it also works on my side! |
As for your point 1, I decided not to create a The changes are in #141 and soon in WriteVTK 1.19.0. |
First, thank you very much for writing such a great package! I am the author of the
Peridynamics.jl
package and have usedWriteVTK.jl
from the beginning to export simulation results.Recently, we extended the package to support MPI and now use
pvtk_grid
:https://github.com/kaipartmann/Peridynamics.jl/blob/7af8831891c4ae02fe27e3c26c932e2a4a16186f/src/auxiliary/io.jl#L83-L88
However, with this change, we stumbled across two minor issues.
1. Relative paths to the
vtu
files inside thepvtu
fileWe export files into a subdirectory
vtk
. To have relocatable files, we mustcd
into this directory where we want to write all the vtk files and then give only a filename (without the path to the subdirectory) topvtk_grid
. Otherwise, a path relative to the pwd of the simulation would be written inside thepvtu
file, leading to problems with ParaView.It would be nice to have a keyword within
pvtk_grid
, likerelative_path=true
, so that the path to thevtu
files inside thepvtu
file would be relative to the latter. Would it be possible to add such a feature?2. Field data
We export the current time as field data, as specified in this example:
https://juliavtk.github.io/WriteVTK.jl/stable/grids/datasets/#Writing-datasets
However, the field data fields are missing in the written
.pvtu
file. Previously, when using just.vtu
files, ParaView was able to find thetime
array, but now it doesn't.What is the recommended workflow for writing time meta data in parallel vtk files?
The text was updated successfully, but these errors were encountered: