forked from mdelorme/fv3d
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathIOManager.h
195 lines (167 loc) · 7.39 KB
/
IOManager.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
#pragma once
#include <highfive/H5Easy.hpp>
#include <ostream>
#include <iomanip>
// https://visit-sphinx-github-user-manual.readthedocs.io/en/3.4rc/data_into_visit/XdmfFormat.html
#include "SimInfo.h"
using namespace H5Easy;
namespace fv3d {
// xdmf strings
namespace {
char str_xdmf_header[] = R"xml(<?xml version="1.0" ?>
<!DOCTYPE Xdmf SYSTEM "Xdmf.dtd" []>
<Xdmf Version="2.0">
<Domain CollectionType="Temporal">
<Grid Name="MainTimeSeries" GridType="Collection" CollectionType="Temporal">
<Topology Name="Main Topology" TopologyType="3DSMesh" NumberOfElements="%d %d %d"/>
<Geometry Name="Main Geometry" GeometryType="X_Y_Z">
<DataItem Dimensions="%d %d %d" NumberType="Float" Precision="8" Format="HDF">%s:/x</DataItem>
<DataItem Dimensions="%d %d %d" NumberType="Float" Precision="8" Format="HDF">%s:/y</DataItem>
<DataItem Dimensions="%d %d %d" NumberType="Float" Precision="8" Format="HDF">%s:/z</DataItem>
</Geometry>
)xml";
#define format_xdmf_header(params) \
params.Nz + 1, params.Ny + 1, params.Nx + 1, \
params.Nz + 1, params.Ny + 1, params.Nx + 1, (params.filename_out + ".h5").c_str(), \
params.Nz + 1, params.Ny + 1, params.Nx + 1, (params.filename_out + ".h5").c_str(), \
params.Nz + 1, params.Ny + 1, params.Nx + 1, (params.filename_out + ".h5").c_str()
char str_xdmf_footer[] =
R"xml(
</Grid>
</Domain>
</Xdmf>)xml";
char str_xdmf_ite_header[] =
R"xml(
<Grid Name="Cells" GridType="Uniform">
<Time TimeType="Single" Value="%lf" />
<Topology Reference="//Topology[@Name='Main Topology']" />
<Geometry Reference="//Geometry[@Name='Main Geometry']" />)xml";
char str_xdmf_scalar_field[] =
R"xml(
<Attribute Name="%s" AttributeType="Scalar" Center="Cell">
<DataItem Dimensions="%d %d %d" NumberType="Float" Precision="8" Format="HDF">%s:/%s/%s</DataItem>
</Attribute>)xml";
#define format_xdmf_scalar_field(params, group, field) \
field, \
params.Nz, params.Ny, params.Nx, \
(params.filename_out + ".h5").c_str(), group.c_str(), field
char str_xdmf_vector_field[] =
R"xml(
<Attribute Name="%s" AttributeType="Vector" Center="Cell">
<DataItem Dimensions="%d %d %d 3" ItemType="Function" Function="JOIN($0, $1, $2)">
<DataItem Dimensions="%d %d %d" NumberType="Float" Precision="8" Format="HDF">%s:/%s/%s</DataItem>
<DataItem Dimensions="%d %d %d" NumberType="Float" Precision="8" Format="HDF">%s:/%s/%s</DataItem>
<DataItem Dimensions="%d %d %d" NumberType="Float" Precision="8" Format="HDF">%s:/%s/%s</DataItem>
</DataItem>
</Attribute>)xml";
#define format_xdmf_vector_field(params, group, name, field_x, field_y, field_z) \
name, \
params.Nz, params.Ny, params.Nx, \
params.Nz, params.Ny, params.Nx, \
(params.filename_out + ".h5").c_str(), group.c_str(), field_x, \
params.Nz, params.Ny, params.Nx, \
(params.filename_out + ".h5").c_str(), group.c_str(), field_y, \
params.Nz, params.Ny, params.Nx, \
(params.filename_out + ".h5").c_str(), group.c_str(), field_z
char str_xdmf_ite_footer[] =
R"xml(
</Grid>
)xml";
} // anonymous namespace
class IOManager {
public:
Params params;
IOManager(Params ¶ms)
: params(params) {};
~IOManager() = default;
void saveSolution(const Array &Q, int iteration, real_t t, real_t dt) {
std::ostringstream oss;
std::setw(4);
std::setfill('0');
oss << "ite_" << iteration;
std::string path = oss.str();
auto flag_h5 = (iteration == 0 ? File::Truncate : File::ReadWrite);
auto flag_xdmf = (iteration == 0 ? "w+" : "r+");
File file(params.filename_out + ".h5", flag_h5);
FILE* xdmf_fd = fopen((params.filename_out + ".xdmf").c_str(), flag_xdmf);
if (iteration == 0) {
file.createAttribute("Ntx", params.Ntx);
file.createAttribute("Nty", params.Nty);
file.createAttribute("Ntz", params.Ntz);
file.createAttribute("Nx", params.Nx);
file.createAttribute("Ny", params.Ny);
file.createAttribute("Nz", params.Nz);
file.createAttribute("ibeg", params.ibeg);
file.createAttribute("iend", params.iend);
file.createAttribute("jbeg", params.jbeg);
file.createAttribute("jend", params.jend);
file.createAttribute("kbeg", params.kbeg);
file.createAttribute("kend", params.kend);
file.createAttribute("problem", params.problem);
std::vector<real_t> x, y, z;
// -- vertex pos
for (int k=params.kbeg; k <= params.kend; ++k) {
for (int j=params.jbeg; j <= params.jend; ++j) {
for (int i=params.ibeg; i <= params.iend; ++i) {
x.push_back((i-params.ibeg) * params.dx);
y.push_back((j-params.jbeg) * params.dy);
z.push_back((k-params.kbeg) * params.dz);
}
}
}
file.createDataSet("x", x);
file.createDataSet("y", y);
file.createDataSet("z", z);
fprintf(xdmf_fd, str_xdmf_header, format_xdmf_header(params));
fprintf(xdmf_fd, "%s", str_xdmf_footer);
}
using Table = std::vector<std::vector<std::vector<real_t>>>;
auto Qhost = Kokkos::create_mirror(Q);
Kokkos::deep_copy(Qhost, Q);
Table trho, tu, tv, tw, tprs;
for (int k=params.kbeg; k<params.kend; ++k) {
std::vector<std::vector<real_t>> rcrho, rcu, rcv, rcw, rcprs;
for (int j=params.jbeg; j<params.jend; ++j) {
std::vector<real_t> rrho, ru, rv, rw, rprs;
for (int i=params.ibeg; i<params.iend; ++i) {
real_t rho = Qhost(k, j, i, IR);
real_t u = Qhost(k, j, i, IU);
real_t v = Qhost(k, j, i, IV);
real_t w = Qhost(k, j, i, IW);
real_t p = Qhost(k, j, i, IP);
rrho.push_back(rho);
ru.push_back(u);
rv.push_back(v);
rw.push_back(w);
rprs.push_back(p);
}
rcrho.push_back(rrho);
rcu.push_back(ru);
rcv.push_back(rv);
rcw.push_back(rw);
rcprs.push_back(rprs);
}
trho.push_back(rcrho);
tu.push_back(rcu);
tv.push_back(rcv);
tw.push_back(rcw);
tprs.push_back(rcprs);
}
auto group = file.createGroup(path);
group.createDataSet("rho", trho);
group.createDataSet("u", tu);
group.createDataSet("v", tv);
group.createDataSet("w", tw);
group.createDataSet("prs", tprs);
group.createAttribute("time", t);
fseek(xdmf_fd, -sizeof(str_xdmf_footer), SEEK_END);
fprintf(xdmf_fd, str_xdmf_ite_header, t);
fprintf(xdmf_fd, str_xdmf_scalar_field, format_xdmf_scalar_field(params, path, "rho"));
fprintf(xdmf_fd, str_xdmf_vector_field, format_xdmf_vector_field(params, path, "velocity", "u", "v", "w"));
fprintf(xdmf_fd, str_xdmf_scalar_field, format_xdmf_scalar_field(params, path, "prs"));
fprintf(xdmf_fd, "%s", str_xdmf_ite_footer);
fprintf(xdmf_fd, "%s", str_xdmf_footer);
fclose(xdmf_fd);
}
};
}