PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
io.cpp
Go to the documentation of this file.
1
8// STL
9#include <fstream>
10#include <sstream>
11#include <vector>
12
13// PETSc
14#include <petscviewerhdf5.h>
15
16// PetIBM
17#include <petibm/io.h>
18
19namespace petibm
20{
21namespace io
22{
23PetscErrorCode readLagrangianPoints(const std::string &file, PetscInt &nPts,
24 type::RealVec2D &coords)
25{
26 PetscFunctionBeginUser;
27
28 std::ifstream inFile(file);
29 std::string line;
30 PetscInt c = 0;
31 PetscInt dim = 0;
32
33 // check if we successfully open the file
34 if (!inFile.good())
35 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_FILE_READ,
36 "Opening or reading body file %s failed!", file.c_str());
37
38 // read the total number of points
39 if (std::getline(inFile, line))
40 {
41 std::stringstream sline(line);
42
43 if (!(sline >> nPts))
44 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_FILE_READ,
45 "Can't read the total number of points in file %s !\n",
46 file.c_str());
47
48 if (sline.peek() != EOF)
49 SETERRQ1(
50 PETSC_COMM_WORLD, PETSC_ERR_FILE_READ,
51 "The first line in file %s contains more than one integer. "
52 "Please check the format.\n",
53 file.c_str());
54 }
55 else
56 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_FILE_READ,
57 "Error while reading the first line in file %s !\n",
58 file.c_str());
59
60 // get the dimension from the first coordinate set; initialize coords
61 {
62 std::getline(inFile, line);
63 std::stringstream sline(line);
64 PetscReal temp;
65 while (sline >> temp) dim += 1;
66
67 if (dim == 0)
68 SETERRQ1(
69 PETSC_COMM_WORLD, PETSC_ERR_FILE_READ,
70 "Could not calculate the dimension from the first coordinate "
71 "set in the file %s!\n",
72 file.c_str());
73
74 // initialize the size of coordinate array
75 coords = type::RealVec2D(nPts, type::RealVec1D(dim, 0.0));
76
77 // read again to get first coordinate set
78 sline.str(line);
79 sline.clear();
80 for (PetscInt d = 0; d < dim; ++d) sline >> coords[0][d];
81
82 // increase c
83 c += 1;
84 }
85
86 // loop through all points (from the second coordinate set)
87 while (std::getline(inFile, line))
88 {
89 std::stringstream sline(line);
90
91 for (PetscInt d = 0; d < dim; ++d)
92 if (!(sline >> coords[c][d]))
93 SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_FILE_READ,
94 "The number of doubles at line %d in file %s does not "
95 "match the dimension.\n",
96 c + 2, file.c_str());
97
98 if (sline.peek() != EOF)
99 SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_FILE_READ,
100 "The number of doubles at line %d in file %s does not "
101 "match the dimension.\n",
102 c + 2, file.c_str());
103
104 c += 1;
105 }
106
107 // close the file
108 inFile.close();
109
110 // check the total number of points read
111 if (c != nPts)
112 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_FILE_READ,
113 "The total number of coordinates read in does not match the "
114 "number specified at the first line in file %s !\n",
115 file.c_str());
116
117 PetscFunctionReturn(0);
118} // readLagrangianPoints
119
120PetscErrorCode print(const std::string &info)
121{
122 PetscFunctionBeginUser;
123
124 PetscErrorCode ierr;
125
126 ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, info.c_str());
127 CHKERRQ(ierr);
128
129 ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT);
130 CHKERRQ(ierr);
131
132 ierr = PetscPrintf(PETSC_COMM_WORLD, "\n"); CHKERRQ(ierr);
133
134 PetscFunctionReturn(0);
135} // print
136
137PetscErrorCode writeHDF5Vecs(const MPI_Comm comm, const std::string &filePath,
138 const std::string &loc,
139 const std::vector<std::string> &names,
140 const std::vector<Vec> &vecs,
141 const PetscFileMode mode)
142{
143 PetscErrorCode ierr;
144 PetscViewer viewer;
145
146 PetscFunctionBeginUser;
147
148 // create viewer
149 ierr = PetscViewerCreate(comm, &viewer); CHKERRQ(ierr);
150 ierr = PetscViewerSetType(viewer, PETSCVIEWERHDF5); CHKERRQ(ierr);
151 ierr = PetscViewerFileSetMode(viewer, mode); CHKERRQ(ierr);
152 ierr = PetscViewerFileSetName(viewer, filePath.c_str()); CHKERRQ(ierr);
153
154 ierr = PetscViewerHDF5PushGroup(viewer, loc.c_str()); CHKERRQ(ierr);
155
156 for (unsigned int i = 0; i < vecs.size(); ++i)
157 {
158 ierr = PetscObjectSetName((PetscObject)vecs[i], names[i].c_str());
159 CHKERRQ(ierr);
160
161 ierr = VecView(vecs[i], viewer); CHKERRQ(ierr);
162 }
163
164 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
165
166 PetscFunctionReturn(0);
167} // writeHDF5Vecs
168
169PetscErrorCode writeHDF5Vecs(const MPI_Comm comm, const std::string &filePath,
170 const std::string &loc,
171 const std::vector<std::string> &names,
172 const std::vector<PetscInt> &n,
173 const std::vector<PetscReal *> &vecs,
174 const PetscFileMode mode)
175{
176 PetscErrorCode ierr;
177 PetscViewer viewer;
178
179 PetscFunctionBeginUser;
180
181 ierr = PetscViewerCreate(comm, &viewer); CHKERRQ(ierr);
182 ierr = PetscViewerSetType(viewer, PETSCVIEWERHDF5); CHKERRQ(ierr);
183 ierr = PetscViewerFileSetMode(viewer, mode); CHKERRQ(ierr);
184 ierr = PetscViewerFileSetName(viewer, filePath.c_str()); CHKERRQ(ierr);
185
186 ierr = PetscViewerHDF5PushGroup(viewer, loc.c_str()); CHKERRQ(ierr);
187
188 for (unsigned int i = 0; i < vecs.size(); i++)
189 {
190 Vec temp;
191 ierr =
192 VecCreateMPIWithArray(comm, 1, n[i], PETSC_DECIDE, nullptr, &temp);
193 CHKERRQ(ierr);
194 ierr = PetscObjectSetName((PetscObject)temp, names[i].c_str());
195 CHKERRQ(ierr);
196 ierr = VecPlaceArray(temp, vecs[i]); CHKERRQ(ierr);
197 ierr = VecView(temp, viewer); CHKERRQ(ierr);
198 ierr = VecResetArray(temp); CHKERRQ(ierr);
199 ierr = VecDestroy(&temp); CHKERRQ(ierr);
200 }
201
202 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
203
204 PetscFunctionReturn(0);
205} // writeHDF5Vecs
206
207PetscErrorCode writeHDF5Vecs(const MPI_Comm comm, const std::string &filePath,
208 const std::string &loc,
209 const std::vector<std::string> &names,
210 const type::RealVec2D &vecs,
211 const PetscFileMode mode)
212{
213 PetscErrorCode ierr;
214 PetscViewer viewer;
215
216 PetscFunctionBeginUser;
217
218 ierr = PetscViewerCreate(comm, &viewer); CHKERRQ(ierr);
219 ierr = PetscViewerSetType(viewer, PETSCVIEWERHDF5); CHKERRQ(ierr);
220 ierr = PetscViewerFileSetMode(viewer, mode); CHKERRQ(ierr);
221 ierr = PetscViewerFileSetName(viewer, filePath.c_str()); CHKERRQ(ierr);
222
223 ierr = PetscViewerHDF5PushGroup(viewer, loc.c_str()); CHKERRQ(ierr);
224
225 for (unsigned int i = 0; i < vecs.size(); i++)
226 {
227 Vec temp;
228 ierr = VecCreateMPIWithArray(comm, 1, vecs[i].size(), PETSC_DECIDE,
229 nullptr, &temp); CHKERRQ(ierr);
230 ierr = PetscObjectSetName((PetscObject)temp, names[i].c_str());
231 CHKERRQ(ierr);
232 ierr = VecPlaceArray(temp, vecs[i].data()); CHKERRQ(ierr);
233 ierr = VecView(temp, viewer); CHKERRQ(ierr);
234 ierr = VecResetArray(temp); CHKERRQ(ierr);
235 ierr = VecDestroy(&temp); CHKERRQ(ierr);
236 }
237
238 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
239
240 PetscFunctionReturn(0);
241} // writeHDF5Vecs
242
243PetscErrorCode readHDF5Vecs(const MPI_Comm comm, const std::string &filePath,
244 const std::string &loc,
245 const std::vector<std::string> &names,
246 std::vector<Vec> &vecs)
247{
248 PetscErrorCode ierr;
249 PetscViewer viewer;
250
251 PetscFunctionBeginUser;
252
253 // create viewer
254 ierr = PetscViewerCreate(comm, &viewer); CHKERRQ(ierr);
255 ierr = PetscViewerSetType(viewer, PETSCVIEWERHDF5); CHKERRQ(ierr);
256 ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ); CHKERRQ(ierr);
257 ierr = PetscViewerFileSetName(viewer, filePath.c_str()); CHKERRQ(ierr);
258
259 ierr = PetscViewerHDF5PushGroup(viewer, loc.c_str()); CHKERRQ(ierr);
260
261 for (unsigned int i = 0; i < vecs.size(); ++i)
262 {
263 ierr = PetscObjectSetName((PetscObject)vecs[i], names[i].c_str());
264 CHKERRQ(ierr);
265
266 ierr = VecLoad(vecs[i], viewer); CHKERRQ(ierr);
267 }
268
269 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
270
271 PetscFunctionReturn(0);
272} // readHDF5Vecs
273
274PetscErrorCode writePetscLog(const MPI_Comm comm, const std::string &filePath)
275{
276 PetscErrorCode ierr;
277 PetscViewer viewerLog;
278
279 PetscFunctionBeginUser;
280
281 ierr = PetscViewerCreate(comm, &viewerLog); CHKERRQ(ierr);
282 ierr = PetscViewerSetType(viewerLog, PETSCVIEWERASCII); CHKERRQ(ierr);
283 ierr = PetscViewerFileSetMode(viewerLog, FILE_MODE_WRITE); CHKERRQ(ierr);
284 ierr = PetscViewerFileSetName(viewerLog, filePath.c_str()); CHKERRQ(ierr);
285 ierr = PetscLogView(viewerLog); CHKERRQ(ierr);
286 ierr = PetscViewerDestroy(&viewerLog); CHKERRQ(ierr);
287
288 PetscFunctionReturn(0);
289} // writePetscLog
290
291} // end of namespace io
292} // end of namespace petibm
PetscErrorCode writeHDF5Vecs(const MPI_Comm comm, const std::string &filePath, const std::string &loc, const std::vector< std::string > &names, const std::vector< Vec > &vecs, const PetscFileMode mode=FILE_MODE_WRITE)
Write a vector of Vec objects to a HDF5 file.
Definition: io.cpp:137
PetscErrorCode readHDF5Vecs(const MPI_Comm comm, const std::string &filePath, const std::string &loc, const std::vector< std::string > &names, std::vector< Vec > &vecs)
Read a vector of Vec objects from a HDF5 file.
Definition: io.cpp:243
PetscErrorCode readLagrangianPoints(const std::string &file, PetscInt &nPts, type::RealVec2D &coords)
Read the number Lagrangian points and their coordinates from a file.
Definition: io.cpp:23
PetscErrorCode writePetscLog(const MPI_Comm comm, const std::string &filePath)
Write a summary of the PETSc logging into a ASCII file.
Definition: io.cpp:274
PetscErrorCode print(const std::string &info)
Print information of a parallel object to standard output.
Definition: io.cpp:120
std::vector< RealVec1D > RealVec2D
2D std::vector holding PetscReal.
Definition: type.h:142
std::vector< PetscReal > RealVec1D
1D std::vector holding PetscReal.
Definition: type.h:140
Prototypes of I/O functions.
A toolbox for building flow solvers.
Definition: bodypack.h:52