Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
D
DG1d.jl
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Snippets
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Model registry
Operate
Environments
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
DG
DG1d.jl
Commits
9a635bc8
Commit
9a635bc8
authored
1 year ago
by
Florian Atteneder
Browse files
Options
Downloads
Patches
Plain Diff
clean up xdm
parent
f75fa9db
Branches
fa/amr-1d
No related tags found
No related merge requests found
Changes
1
Hide whitespace changes
Inline
Side-by-side
Showing
1 changed file
src/xdmf.jl
+54
-45
54 additions, 45 deletions
src/xdmf.jl
with
54 additions
and
45 deletions
src/xdmf.jl
+
54
−
45
View file @
9a635bc8
...
@@ -4,6 +4,7 @@ const XDMF_END = """
...
@@ -4,6 +4,7 @@ const XDMF_END = """
"""
"""
mutable struct
XDMFOutput
mutable struct
XDMFOutput
buffer
::
IOBuffer
dir
::
String
dir
::
String
h5_filename
::
String
h5_filename
::
String
mesh_filename
::
String
mesh_filename
::
String
...
@@ -25,7 +26,9 @@ function XDMFOutput(dir)
...
@@ -25,7 +26,9 @@ function XDMFOutput(dir)
open
(
mesh_filename
,
"w"
)
do
file
open
(
mesh_filename
,
"w"
)
do
file
end
end
return
XDMFOutput
(
dir
,
h5_filename
,
mesh_filename
,
0
,
0
)
buffer
=
IOBuffer
()
return
XDMFOutput
(
buffer
,
dir
,
h5_filename
,
mesh_filename
,
0
,
0
)
end
end
...
@@ -35,21 +38,20 @@ function xdmfwrite_mesh(xdmf::XDMFOutput, mesh::Mesh)
...
@@ -35,21 +38,20 @@ function xdmfwrite_mesh(xdmf::XDMFOutput, mesh::Mesh)
throw
(
ArgumentError
(
"Mesh has already been write"
))
throw
(
ArgumentError
(
"Mesh has already been write"
))
end
end
buffer
=
IOBuffer
()
#
buffer = IOBuffer()
write_mesh_header
(
buffer
)
write_mesh_header
(
xdmf
.
buffer
)
# Using "w" here to overwrite in case previous calls failed
# Using "w" here to overwrite in case previous calls failed
h5open
(
xdmf
.
h5_filename
,
"w"
)
do
file
h5open
(
xdmf
.
h5_filename
,
"w"
)
do
file
write_mesh
(
file
,
buffer
,
mesh
)
write_mesh
(
file
,
xdmf
.
buffer
,
mesh
)
end
end
write_mesh_end
(
buffer
)
write_mesh_end
(
xdmf
.
buffer
)
# finish up xml files
# finish up xml files
# Using "w" here to overwrite in case previous calls failed
# Using "w" here to overwrite in case previous calls failed
open
(
xdmf
.
mesh_filename
,
"a"
)
do
file
open
(
xdmf
.
mesh_filename
,
"a"
)
do
file
seekstart
(
buffer
)
write
(
file
,
take!
(
xdmf
.
buffer
))
write
(
file
,
buffer
)
end
end
# updating count after all writes succeeded
# updating count after all writes succeeded
...
@@ -59,13 +61,14 @@ end
...
@@ -59,13 +61,14 @@ end
function
write_mesh_header
(
buffer
)
function
write_mesh_header
(
buffer
)
w
rit
e
(
buffer
,
"""
p
ri
n
t
(
buffer
,
"""
<?xml version="
1.0
" ?>
<?xml version="
1.0
" ?>
<!DOCTYPE Xdmf SYSTEM "
Xdmf
.
dtd
" []>
<!DOCTYPE Xdmf SYSTEM "
Xdmf
.
dtd
" []>
<Xdmf xmlns:xi="
http
://
www
.
w3
.
org
/
2003
/
XInclude
" Version="
3.0
">
<Xdmf xmlns:xi="
http
://
www
.
w3
.
org
/
2003
/
XInclude
" Version="
3.0
">
<Domain>
<Domain>
<Grid GridType="
Collection
" CollectionType="
Spatial
">
<Grid GridType="
Collection
" CollectionType="
Spatial
">
"""
)
"""
)
return
end
end
...
@@ -131,61 +134,67 @@ size3d(sz::NTuple{2}) = (sz[1],sz[2],1)
...
@@ -131,61 +134,67 @@ size3d(sz::NTuple{2}) = (sz[1],sz[2],1)
# size3d(sz::NTuple{3}) = (sz[1],sz[2],sz[3])
# size3d(sz::NTuple{3}) = (sz[1],sz[2],sz[3])
function
xdmf
write_scalar
(
xdmf
::
XDMFOutput
,
mesh
::
Mesh
,
t
,
var
s
,
name
s
)
function
write_scalar
(
h5
,
buffer
,
mesh
::
Mesh
,
dir
,
tcount
,
t
,
var
,
name
)
buffers
=
IOBuffer
[]
write
(
buffer
,
"""
for
_
in
vars
b
=
IOBuffer
()
write
(
b
,
"""
<Grid GridType="
Collection
" CollectionType="
Temporal
"><Grid GridType="
Collection
" CollectionType="
Spatial
">
<Grid GridType="
Collection
" CollectionType="
Temporal
"><Grid GridType="
Collection
" CollectionType="
Spatial
">
<Time Value="
$
t
"/>
<Time Value="
$
t
"/>
"""
)
"""
)
push!
(
buffers
,
b
)
end
tcount
=
xdmf
.
tcount
+
1
variter
=
eachcell
(
mesh
,
var
)
celliters
=
[
eachcell
(
mesh
,
var
)
for
var
in
vars
]
h5open
(
xdmf
.
h5_filename
,
"r+"
)
do
h5
# a group to gather grids from one time step
tgroup
=
create_group
(
h5
,
string
(
tcount
))
# an other group the new time step
h5group_t
=
create_group
(
h5
,
string
(
tcount
))
for
(
i
,
v
)
in
enumerate
(
variter
)
# sub groups for each mesh
nx
,
ny
,
nz
=
size3d
(
size
(
v
))
groups
=
[
create_group
(
h5group_t
,
string
(
i
))
for
i
=
1
:
ncells
(
mesh
)
]
g
=
create_group
(
tgroup
,
string
(
i
))
g
[
name
]
=
v
for
(
cellit
,
name
,
buffer
)
in
zip
(
celliters
,
names
,
buffers
)
write
(
buffer
,
"""
for
(
i
,(
cellvar
,
g
))
in
enumerate
(
zip
(
cellit
,
groups
))
nx
,
ny
,
nz
=
size3d
(
size
(
cellvar
))
g
[
name
]
=
view
(
cellvar
,
:
)
write
(
buffer
,
"""
<Grid Name="
$
i
"><xi:include href="
mesh
.
xdmf
" xpointer="
xpointer
(
/
Xdmf
/
Domain
/
Grid
[
1
]
/
Grid
[
$
i
]
/
child
::*
)
"/><Attribute Name="
$
name
"><DataItem Dimensions="
$
nz
$
ny
$
nx
" Format="
HDF
"> data.h5:/
$
tcount/
$
i/
$
name </DataItem></Attribute></Grid>
<Grid Name="
$
i
"><xi:include href="
mesh
.
xdmf
" xpointer="
xpointer
(
/
Xdmf
/
Domain
/
Grid
[
1
]
/
Grid
[
$
i
]
/
child
::*
)
"/><Attribute Name="
$
name
"><DataItem Dimensions="
$
nz
$
ny
$
nx
" Format="
HDF
"> data.h5:/
$
tcount/
$
i/
$
name </DataItem></Attribute></Grid>
"""
)
"""
)
end
write
(
buffer
,
"</Grid></Grid>"
)
end
end
end
# finish up xml files
write
(
buffer
,
"</Grid></Grid>"
)
for
(
buffer
,
name
)
in
zip
(
buffers
,
names
)
filename
=
joinpath
(
xdmf
.
dir
,
"
$
name.xdmf"
)
filename
=
joinpath
(
dir
,
"
$
name.xdmf"
)
existed
=
isfile
(
filename
)
existed
=
isfile
(
filename
)
open
(
filename
,
"a"
)
do
file
open
(
filename
,
"a"
)
do
file
if
!
existed
if
!
existed
print
(
file
,
"""
print
(
file
,
"""
<?xml version="
1.0
" ?>
<?xml version="
1.0
" ?>
<!DOCTYPE Xdmf SYSTEM "
Xdmf
.
dtd
" []>
<!DOCTYPE Xdmf SYSTEM "
Xdmf
.
dtd
" []>
<Xdmf xmlns:xi="
http
://
www
.
w3
.
org
/
2003
/
XInclude
" Version="
3.0
">
<Xdmf xmlns:xi="
http
://
www
.
w3
.
org
/
2003
/
XInclude
" Version="
3.0
">
<Domain>
<Domain>
"""
)
"""
)
else
else
skip
(
file
,
-
length
(
XDMF_END
))
skip
(
file
,
-
length
(
XDMF_END
))
end
end
seekstart
(
buffer
)
write
(
file
,
take!
(
buffer
))
write
(
file
,
buffer
)
write
(
file
,
XDMF_END
)
write
(
file
,
XDMF_END
)
end
return
end
function
xdmfwrite_scalar
(
xdmf
::
XDMFOutput
,
mesh
::
Mesh
,
t
,
vars
,
names
)
tcount
=
xdmf
.
tcount
+
1
h5
=
h5open
(
xdmf
.
h5_filename
,
"r+"
)
try
for
(
var
,
name
)
in
zip
(
vars
,
names
)
write_scalar
(
h5
,
xdmf
.
buffer
,
mesh
,
xdmf
.
dir
,
tcount
,
t
,
var
,
name
)
end
end
catch
e
rethrow
(
e
)
finally
close
(
h5
)
end
end
# updating count after all writes succeeded
# updating count after all writes succeeded
xdmf
.
tcount
=
tcount
xdmf
.
tcount
=
tcount
return
end
end
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment