Nektar++
OutputVtkBase.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: OutputVtkBase.cpp
4//
5// For more information, please see: http://www.nektar.info/
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: VTK file format output base class.
32//
33////////////////////////////////////////////////////////////////////////////////
34
35#include <iomanip>
36#include <set>
37#include <string>
38using namespace std;
39
40#include <boost/format.hpp>
41
43
44#include "OutputVtkBase.h"
45
46namespace Nektar::FieldUtils
47{
48
49// Disable the base VTK factory if using the VTK library, this is so we can
50// register the same extension for both.
51#if !NEKTAR_USING_VTK
55 "Writes a VTU file.");
56#endif
57
59{
61}
62
64{
65}
66
67void OutputVtkBase::v_OutputFromPts(po::variables_map &vm)
68{
69 int i, j;
70 LibUtilities::PtsFieldSharedPtr fPts = m_f->m_fieldPts;
71
72 // Extract the output filename and extension
73 string filename = PrepareOutput(vm);
74
75 // Write solution.
76 ofstream outfile(filename.c_str());
77 WriteVtkHeader(outfile);
78 int nfields = 0;
79 int dim = fPts->GetDim();
80
81 int nvert = 1;
82 int vtktype = 1;
83 switch (fPts->GetPtsType())
84 {
87 {
89 "VTK output needs setting up for ePtsFile and ePtsLine");
90 break;
91 }
93 {
95 "VTK output needs setting up for PtsPlane");
96 break;
97 }
99 {
101 "VTK output needs setting up for PtsBox");
102 break;
103 }
105 {
106 nvert = 2;
107 vtktype = 3;
108 break;
109 }
111 {
112 nvert = 3;
113 vtktype = 5;
114 break;
115 }
117 {
118 nvert = 4;
119 vtktype = 10;
120 break;
121 }
122 default:
123 NEKERROR(ErrorUtil::efatal, "ptsType not supported yet.");
124 }
125
126 vector<Array<OneD, int>> ptsConn;
127 fPts->GetConnectivity(ptsConn);
128
129 nfields = fPts->GetNFields();
130
131 int nPts = fPts->GetNpoints();
132 int numBlocks = 0;
133 for (i = 0; i < ptsConn.size(); ++i)
134 {
135 numBlocks += ptsConn[i].size() / nvert;
136 }
137
138 // write out pieces of data.
139 outfile << " <Piece NumberOfPoints=\"" << nPts << "\" NumberOfCells=\""
140 << numBlocks << "\">" << endl;
141 outfile << " <Points>" << endl;
142 outfile << " <DataArray type=\"Float64\" "
143 << "NumberOfComponents=\"" << 3 << R"(" format="ascii">)" << endl;
144 for (i = 0; i < nPts; ++i)
145 {
146 for (j = 0; j < dim; ++j)
147 {
148 outfile << " " << setprecision(8) << scientific
149 << fPts->GetPointVal(j, i) << " ";
150 }
151 for (j = dim; j < 3; ++j)
152 {
153 // pack to 3D since paraview does not seem to handle 2D
154 outfile << " 0.000000";
155 }
156 outfile << endl;
157 }
158 outfile << " </DataArray>" << endl;
159 outfile << " </Points>" << endl;
160 outfile << " <Cells>" << endl;
161 outfile << " <DataArray type=\"Int32\" "
162 << R"(Name="connectivity" format="ascii">)" << endl;
163
164 // dump connectivity data if it exists
165 outfile << " ";
166 int cnt = 1;
167 for (i = 0; i < ptsConn.size(); ++i)
168 {
169 for (j = 0; j < ptsConn[i].size(); ++j)
170 {
171 outfile << ptsConn[i][j] << " ";
172 if ((!(cnt % nvert)) && cnt)
173 {
174 outfile << std::endl;
175 outfile << " ";
176 }
177 cnt++;
178 }
179 }
180 outfile << " </DataArray>" << endl;
181 outfile << " <DataArray type=\"Int32\" "
182 << R"(Name="offsets" format="ascii">)" << endl;
183
184 outfile << " ";
185 for (i = 0; i < numBlocks; ++i)
186 {
187 outfile << i * nvert + nvert << " ";
188 }
189 outfile << endl;
190 outfile << " </DataArray>" << endl;
191 outfile << " <DataArray type=\"UInt8\" "
192 << R"(Name="types" format="ascii">)" << endl;
193 outfile << " ";
194 for (i = 0; i < numBlocks; ++i)
195 {
196 outfile << vtktype << " ";
197 }
198 outfile << endl;
199 outfile << " </DataArray>" << endl;
200 outfile << " </Cells>" << endl;
201 outfile << " <PointData>" << endl;
202
203 // printing the fields
204 for (j = 0; j < nfields; ++j)
205 {
206 outfile << R"( <DataArray type="Float64" Name=")"
207 << m_f->m_variables[j] << "\">" << endl;
208 outfile << " ";
209 for (i = 0; i < fPts->GetNpoints(); ++i)
210 {
211 outfile << fPts->GetPointVal(dim + j, i) << " ";
212 }
213 outfile << endl;
214 outfile << " </DataArray>" << endl;
215 }
216
217 outfile << " </PointData>" << endl;
218 outfile << " </Piece>" << endl;
219
220 WriteVtkFooter(outfile);
221 cout << "Written file: " << filename << endl;
222
223 // output parallel outline info if necessary
224 if ((m_f->m_comm->GetSpaceComm()->GetRank() == 0) &&
225 (m_f->m_comm->GetSpaceComm()->GetSize() != 1))
226 {
227 WritePVtu(vm);
228 cout << "Written file: " << filename << endl;
229 }
230}
231
232void OutputVtkBase::v_OutputFromExp(po::variables_map &vm)
233{
234 int i, j;
235 // Extract the output filename and extension
236 string filename = PrepareOutput(vm);
237
238 // Write solution.
239 ofstream outfile(filename.c_str());
240 WriteVtkHeader(outfile);
241 int nfields = m_f->m_variables.size();
242
243 int nstrips;
244 m_f->m_session->LoadParameter("Strip_Z", nstrips, 1);
245
246 // Homogeneous strip variant
247 for (int s = 0; s < nstrips; ++s)
248 {
249 // For each field write out field data for each expansion.
250 for (i = 0; i < m_f->m_exp[0]->GetNumElmts(); ++i)
251 {
252 m_f->m_exp[0]->WriteVtkPieceHeader(outfile, i, s);
253
254 // For this expansion write out each field.
255 for (j = 0; j < nfields; ++j)
256 {
257 m_f->m_exp[s * nfields + j]->WriteVtkPieceData(
258 outfile, i, m_f->m_variables[j]);
259 }
260 m_f->m_exp[0]->WriteVtkPieceFooter(outfile, i);
261 }
262 }
263
264 if (m_f->m_exp[0]->GetNumElmts() == 0)
265 {
266 WriteEmptyVtkPiece(outfile);
267 }
268
269 WriteVtkFooter(outfile);
270 cout << "Written file: " << filename << endl;
271
272 // output parallel outline info if necessary
273 if ((m_f->m_comm->GetSpaceComm()->GetRank() == 0) &&
274 (m_f->m_comm->GetSpaceComm()->GetSize() != 1))
275 {
276 WritePVtu(vm);
277 }
278}
279
280void OutputVtkBase::v_OutputFromData([[maybe_unused]] po::variables_map &vm)
281{
283 "OutputVtk can't write using only FieldData. You may need "
284 "to add a mesh XML file to your input files.");
285}
286
287fs::path OutputVtkBase::v_GetPath(std::string &filename,
288 [[maybe_unused]] po::variables_map &vm)
289{
290 int nprocs = m_f->m_comm->GetSpaceComm()->GetSize();
291 fs::path specPath;
292 if (nprocs == 1)
293 {
294 specPath = fs::path(filename);
295 }
296 else
297 {
298 // replace .vtu by _vtu
299 int dot = filename.find_last_of('.');
300 string path = filename.substr(0, dot) + "_vtu";
301 specPath = fs::path(path);
302 }
303 return fs::path(specPath);
304}
305
306fs::path OutputVtkBase::v_GetFullOutName(std::string &filename,
307 po::variables_map &vm)
308{
309 int nprocs = m_f->m_comm->GetSpaceComm()->GetSize();
310
311 fs::path fulloutname;
312 if (nprocs == 1)
313 {
314 fulloutname = filename;
315 }
316 else
317 {
318 // Guess at filename that might belong to this process.
319 boost::format pad("P%1$07d.%2$s");
320 pad % m_f->m_comm->GetSpaceComm()->GetRank() % "vtu";
321
322 // Generate full path name
323 fs::path specPath = GetPath(filename, vm);
324 fs::path poutfile(pad.str());
325 fulloutname = specPath / poutfile;
326 }
327 return fulloutname;
328}
329
330void OutputVtkBase::WriteVtkHeader(std::ostream &outfile)
331{
332 outfile << "<?xml version=\"1.0\"?>" << endl;
333 outfile << R"(<VTKFile type="UnstructuredGrid" version="0.1" )"
334 << "byte_order=\"LittleEndian\">" << endl;
335 outfile << " <UnstructuredGrid>" << endl;
336}
337
338void OutputVtkBase::WriteVtkFooter(std::ostream &outfile)
339{
340 outfile << " </UnstructuredGrid>" << endl;
341 outfile << "</VTKFile>" << endl;
342}
343
344void OutputVtkBase::WriteEmptyVtkPiece(std::ofstream &outfile)
345{
346 // write out empty piece of data.
347 outfile << " <Piece NumberOfPoints=\"" << 0 << "\" NumberOfCells=\"" << 0
348 << "\">" << endl;
349 outfile << " <Points>" << endl;
350 outfile << " <DataArray type=\"Float64\" "
351 << "NumberOfComponents=\"" << 3 << R"(" format="ascii">)" << endl;
352 outfile << " </DataArray>" << endl;
353 outfile << " </Points>" << endl;
354 outfile << " <Cells>" << endl;
355 outfile << " <DataArray type=\"Int32\" "
356 << R"(Name="connectivity" format="ascii">)" << endl;
357 outfile << " </DataArray>" << endl;
358 outfile << " <DataArray type=\"Int32\" "
359 << R"(Name="offsets" format="ascii">)" << endl;
360
361 outfile << " ";
362 outfile << endl;
363 outfile << " </DataArray>" << endl;
364 outfile << " <DataArray type=\"UInt8\" "
365 << R"(Name="types" format="ascii">)" << endl;
366 outfile << " ";
367 outfile << endl;
368 outfile << " </DataArray>" << endl;
369 outfile << " </Cells>" << endl;
370 outfile << " <PointData>" << endl;
371
372 outfile << " </PointData>" << endl;
373 outfile << " </Piece>" << endl;
374}
375
376void OutputVtkBase::WritePVtu(po::variables_map &vm)
377{
378 string filename = m_config["outfile"].as<string>();
379 int dot = filename.find_last_of('.');
380 string body = filename.substr(0, dot);
381 filename = body + ".pvtu";
382
383 ofstream outfile(filename.c_str());
384
385 int nprocs = m_f->m_comm->GetSpaceComm()->GetSize();
386 string path = LibUtilities::PortablePath(GetPath(filename, vm));
387
388 outfile << "<?xml version=\"1.0\"?>" << endl;
389 outfile << R"(<VTKFile type="PUnstructuredGrid" version="0.1" )"
390 << "byte_order=\"LittleEndian\">" << endl;
391 outfile << "<PUnstructuredGrid GhostLevel=\"0\">" << endl;
392 outfile << "<PPoints> " << endl;
393 outfile << R"(<PDataArray type="Float64" NumberOfComponents=")" << 3
394 << "\"/> " << endl;
395 outfile << "</PPoints>" << endl;
396 outfile << "<PCells>" << endl;
397 outfile << "<PDataArray type=\"Int32\" Name=\"connectivity\" "
398 "NumberOfComponents=\"1\"/>"
399 << endl;
400 outfile << "<PDataArray type=\"Int32\" Name=\"offsets\" "
401 "NumberOfComponents=\"1\"/>"
402 << endl;
403 outfile << "<PDataArray type=\"UInt8\" Name=\"types\" "
404 "NumberOfComponents=\"1\"/>"
405 << endl;
406 outfile << "</PCells>" << endl;
407 outfile << "<PPointData Scalars=\"Material\">" << endl;
408 for (int i = 0; i < m_f->m_variables.size(); ++i)
409 {
410 outfile << R"(<PDataArray type="Float64" Name=")" << m_f->m_variables[i]
411 << "\"/>" << endl;
412 }
413 outfile << "</PPointData>" << endl;
414
415 for (int i = 0; i < nprocs; ++i)
416 {
417 boost::format pad("P%1$07d.vtu");
418 pad % i;
419 outfile << "<Piece Source=\"" << path << "/" << pad.str() << "\"/>"
420 << endl;
421 }
422 outfile << "</PUnstructuredGrid>" << endl;
423 outfile << "</VTKFile>" << endl;
424
425 cout << "Written file: " << filename << endl;
426}
427
428std::string OutputVtkBase::PrepareOutput(po::variables_map &vm)
429{
430 // Extract the output filename and extension
431 string filename = m_config["outfile"].as<string>();
432
433 fs::path specPath = GetPath(filename, vm);
434 fs::path fulloutname = GetFullOutName(filename, vm);
435 filename = LibUtilities::PortablePath(fulloutname);
436
437 if (m_f->m_comm->GetSpaceComm()->GetSize() != 1)
438 {
439 if (m_f->m_comm->GetSpaceComm()->TreatAsRankZero())
440 {
441 try
442 {
443 fs::create_directory(specPath);
444 }
445 catch (fs::filesystem_error &e)
446 {
447 ASSERTL0(false, "Filesystem error: " + string(e.what()));
448 }
449 cout << "Writing files to directory: " << specPath << endl;
450 }
451 m_f->m_comm->GetSpaceComm()->Block();
452 }
453 else
454 {
455 cout << "Writing: " << specPath << endl;
456 }
457 return filename;
458}
459
460} // namespace Nektar::FieldUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
FieldSharedPtr m_f
Field object.
Definition: Module.h:239
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:272
Converter from fld to vtk.
fs::path GetFullOutName(std::string &filename, po::variables_map &vm)
fs::path GetPath(std::string &filename, po::variables_map &vm)
std::string PrepareOutput(po::variables_map &vm)
void WriteEmptyVtkPiece(std::ofstream &outfile)
fs::path v_GetFullOutName(std::string &filename, po::variables_map &vm) override
void WriteVtkHeader(std::ostream &outfile)
static std::shared_ptr< Module > create(FieldSharedPtr f)
Creates an instance of this class.
Definition: OutputVtkBase.h:49
void v_OutputFromData(po::variables_map &vm) override
Write from data to output file.
fs::path v_GetPath(std::string &filename, po::variables_map &vm) override
void WritePVtu(po::variables_map &vm)
void WriteVtkFooter(std::ostream &outfile)
void v_OutputFromExp(po::variables_map &vm) override
Write from m_exp to output file.
void v_OutputFromPts(po::variables_map &vm) override
Write from pts to output file.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:990
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:180
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:47
static std::string PortablePath(const fs::path &path)
create portable path on different platforms for std::filesystem path.
Definition: Filesystem.hpp:56
std::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:184