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/core/ignore_unused.hpp>
41#include <boost/format.hpp>
42
44
45#include "OutputVtkBase.h"
46
47namespace Nektar
48{
49namespace FieldUtils
50{
51
52// Disable the base VTK factory if using the VTK library, this is so we can
53// register the same extension for both.
54#if !NEKTAR_USING_VTK
58 "Writes a VTU file.");
59#endif
60
62{
64}
65
67{
68}
69
70void OutputVtkBase::v_OutputFromPts(po::variables_map &vm)
71{
72 int i, j;
73 LibUtilities::PtsFieldSharedPtr fPts = m_f->m_fieldPts;
74
75 // Extract the output filename and extension
76 string filename = PrepareOutput(vm);
77
78 // Write solution.
79 ofstream outfile(filename.c_str());
80 WriteVtkHeader(outfile);
81 int nfields = 0;
82 int dim = fPts->GetDim();
83
84 int nvert = 1;
85 int vtktype = 1;
86 switch (fPts->GetPtsType())
87 {
90 {
92 "VTK output needs setting up for ePtsFile and ePtsLine");
93 break;
94 }
96 {
98 "VTK output needs setting up for PtsPlane");
99 break;
100 }
102 {
104 "VTK output needs setting up for PtsBox");
105 break;
106 }
108 {
109 nvert = 2;
110 vtktype = 3;
111 break;
112 }
114 {
115 nvert = 3;
116 vtktype = 5;
117 break;
118 }
120 {
121 nvert = 4;
122 vtktype = 10;
123 break;
124 }
125 default:
126 NEKERROR(ErrorUtil::efatal, "ptsType not supported yet.");
127 }
128
129 vector<Array<OneD, int>> ptsConn;
130 fPts->GetConnectivity(ptsConn);
131
132 nfields = fPts->GetNFields();
133
134 int nPts = fPts->GetNpoints();
135 int numBlocks = 0;
136 for (i = 0; i < ptsConn.size(); ++i)
137 {
138 numBlocks += ptsConn[i].size() / nvert;
139 }
140
141 // write out pieces of data.
142 outfile << " <Piece NumberOfPoints=\"" << nPts << "\" NumberOfCells=\""
143 << numBlocks << "\">" << endl;
144 outfile << " <Points>" << endl;
145 outfile << " <DataArray type=\"Float64\" "
146 << "NumberOfComponents=\"" << 3 << "\" format=\"ascii\">" << endl;
147 for (i = 0; i < nPts; ++i)
148 {
149 for (j = 0; j < dim; ++j)
150 {
151 outfile << " " << setprecision(8) << scientific
152 << fPts->GetPointVal(j, i) << " ";
153 }
154 for (j = dim; j < 3; ++j)
155 {
156 // pack to 3D since paraview does not seem to handle 2D
157 outfile << " 0.000000";
158 }
159 outfile << endl;
160 }
161 outfile << " </DataArray>" << endl;
162 outfile << " </Points>" << endl;
163 outfile << " <Cells>" << endl;
164 outfile << " <DataArray type=\"Int32\" "
165 << "Name=\"connectivity\" format=\"ascii\">" << endl;
166
167 // dump connectivity data if it exists
168 outfile << " ";
169 int cnt = 1;
170 for (i = 0; i < ptsConn.size(); ++i)
171 {
172 for (j = 0; j < ptsConn[i].size(); ++j)
173 {
174 outfile << ptsConn[i][j] << " ";
175 if ((!(cnt % nvert)) && cnt)
176 {
177 outfile << std::endl;
178 outfile << " ";
179 }
180 cnt++;
181 }
182 }
183 outfile << " </DataArray>" << endl;
184 outfile << " <DataArray type=\"Int32\" "
185 << "Name=\"offsets\" format=\"ascii\">" << endl;
186
187 outfile << " ";
188 for (i = 0; i < numBlocks; ++i)
189 {
190 outfile << i * nvert + nvert << " ";
191 }
192 outfile << endl;
193 outfile << " </DataArray>" << endl;
194 outfile << " <DataArray type=\"UInt8\" "
195 << "Name=\"types\" format=\"ascii\">" << endl;
196 outfile << " ";
197 for (i = 0; i < numBlocks; ++i)
198 {
199 outfile << vtktype << " ";
200 }
201 outfile << endl;
202 outfile << " </DataArray>" << endl;
203 outfile << " </Cells>" << endl;
204 outfile << " <PointData>" << endl;
205
206 // printing the fields
207 for (j = 0; j < nfields; ++j)
208 {
209 outfile << " <DataArray type=\"Float64\" Name=\""
210 << m_f->m_variables[j] << "\">" << endl;
211 outfile << " ";
212 for (i = 0; i < fPts->GetNpoints(); ++i)
213 {
214 outfile << fPts->GetPointVal(dim + j, i) << " ";
215 }
216 outfile << endl;
217 outfile << " </DataArray>" << endl;
218 }
219
220 outfile << " </PointData>" << endl;
221 outfile << " </Piece>" << endl;
222
223 WriteVtkFooter(outfile);
224 cout << "Written file: " << filename << endl;
225
226 // output parallel outline info if necessary
227 if ((m_f->m_comm->GetSpaceComm()->GetRank() == 0) &&
228 (m_f->m_comm->GetSpaceComm()->GetSize() != 1))
229 {
230 WritePVtu(vm);
231 cout << "Written file: " << filename << endl;
232 }
233}
234
235void OutputVtkBase::v_OutputFromExp(po::variables_map &vm)
236{
237 int i, j;
238 // Extract the output filename and extension
239 string filename = PrepareOutput(vm);
240
241 // Write solution.
242 ofstream outfile(filename.c_str());
243 WriteVtkHeader(outfile);
244 int nfields = m_f->m_variables.size();
245
246 int nstrips;
247 m_f->m_session->LoadParameter("Strip_Z", nstrips, 1);
248
249 // Homogeneous strip variant
250 for (int s = 0; s < nstrips; ++s)
251 {
252 // For each field write out field data for each expansion.
253 for (i = 0; i < m_f->m_exp[0]->GetNumElmts(); ++i)
254 {
255 m_f->m_exp[0]->WriteVtkPieceHeader(outfile, i, s);
256
257 // For this expansion write out each field.
258 for (j = 0; j < nfields; ++j)
259 {
260 m_f->m_exp[s * nfields + j]->WriteVtkPieceData(
261 outfile, i, m_f->m_variables[j]);
262 }
263 m_f->m_exp[0]->WriteVtkPieceFooter(outfile, i);
264 }
265 }
266
267 if (m_f->m_exp[0]->GetNumElmts() == 0)
268 {
269 WriteEmptyVtkPiece(outfile);
270 }
271
272 WriteVtkFooter(outfile);
273 cout << "Written file: " << filename << endl;
274
275 // output parallel outline info if necessary
276 if ((m_f->m_comm->GetSpaceComm()->GetRank() == 0) &&
277 (m_f->m_comm->GetSpaceComm()->GetSize() != 1))
278 {
279 WritePVtu(vm);
280 }
281}
282
283void OutputVtkBase::v_OutputFromData(po::variables_map &vm)
284{
285 boost::ignore_unused(vm);
287 "OutputVtk can't write using only FieldData. You may need "
288 "to add a mesh XML file to your input files.");
289}
290
291fs::path OutputVtkBase::v_GetPath(std::string &filename, po::variables_map &vm)
292{
293 boost::ignore_unused(vm);
294
295 int nprocs = m_f->m_comm->GetSpaceComm()->GetSize();
296 fs::path specPath;
297 if (nprocs == 1)
298 {
299 specPath = fs::path(filename);
300 }
301 else
302 {
303 // replace .vtu by _vtu
304 int dot = filename.find_last_of('.');
305 string path = filename.substr(0, dot) + "_vtu";
306 specPath = fs::path(path);
307 }
308 return fs::path(specPath);
309}
310
311fs::path OutputVtkBase::v_GetFullOutName(std::string &filename,
312 po::variables_map &vm)
313{
314 int nprocs = m_f->m_comm->GetSpaceComm()->GetSize();
315
316 fs::path fulloutname;
317 if (nprocs == 1)
318 {
319 fulloutname = filename;
320 }
321 else
322 {
323 // Guess at filename that might belong to this process.
324 boost::format pad("P%1$07d.%2$s");
325 pad % m_f->m_comm->GetSpaceComm()->GetRank() % "vtu";
326
327 // Generate full path name
328 fs::path specPath = GetPath(filename, vm);
329 fs::path poutfile(pad.str());
330 fulloutname = specPath / poutfile;
331 }
332 return fulloutname;
333}
334
335void OutputVtkBase::WriteVtkHeader(std::ostream &outfile)
336{
337 outfile << "<?xml version=\"1.0\"?>" << endl;
338 outfile << "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" "
339 << "byte_order=\"LittleEndian\">" << endl;
340 outfile << " <UnstructuredGrid>" << endl;
341}
342
343void OutputVtkBase::WriteVtkFooter(std::ostream &outfile)
344{
345 outfile << " </UnstructuredGrid>" << endl;
346 outfile << "</VTKFile>" << endl;
347}
348
349void OutputVtkBase::WriteEmptyVtkPiece(std::ofstream &outfile)
350{
351 // write out empty piece of data.
352 outfile << " <Piece NumberOfPoints=\"" << 0 << "\" NumberOfCells=\"" << 0
353 << "\">" << endl;
354 outfile << " <Points>" << endl;
355 outfile << " <DataArray type=\"Float64\" "
356 << "NumberOfComponents=\"" << 3 << "\" format=\"ascii\">" << endl;
357 outfile << " </DataArray>" << endl;
358 outfile << " </Points>" << endl;
359 outfile << " <Cells>" << endl;
360 outfile << " <DataArray type=\"Int32\" "
361 << "Name=\"connectivity\" format=\"ascii\">" << endl;
362 outfile << " </DataArray>" << endl;
363 outfile << " <DataArray type=\"Int32\" "
364 << "Name=\"offsets\" format=\"ascii\">" << endl;
365
366 outfile << " ";
367 outfile << endl;
368 outfile << " </DataArray>" << endl;
369 outfile << " <DataArray type=\"UInt8\" "
370 << "Name=\"types\" format=\"ascii\">" << endl;
371 outfile << " ";
372 outfile << endl;
373 outfile << " </DataArray>" << endl;
374 outfile << " </Cells>" << endl;
375 outfile << " <PointData>" << endl;
376
377 outfile << " </PointData>" << endl;
378 outfile << " </Piece>" << endl;
379}
380
381void OutputVtkBase::WritePVtu(po::variables_map &vm)
382{
383 string filename = m_config["outfile"].as<string>();
384 int dot = filename.find_last_of('.');
385 string body = filename.substr(0, dot);
386 filename = body + ".pvtu";
387
388 ofstream outfile(filename.c_str());
389
390 int nprocs = m_f->m_comm->GetSpaceComm()->GetSize();
391 string path = LibUtilities::PortablePath(GetPath(filename, vm));
392
393 outfile << "<?xml version=\"1.0\"?>" << endl;
394 outfile << "<VTKFile type=\"PUnstructuredGrid\" version=\"0.1\" "
395 << "byte_order=\"LittleEndian\">" << endl;
396 outfile << "<PUnstructuredGrid GhostLevel=\"0\">" << endl;
397 outfile << "<PPoints> " << endl;
398 outfile << "<PDataArray type=\"Float64\" NumberOfComponents=\"" << 3
399 << "\"/> " << endl;
400 outfile << "</PPoints>" << endl;
401 outfile << "<PCells>" << endl;
402 outfile << "<PDataArray type=\"Int32\" Name=\"connectivity\" "
403 "NumberOfComponents=\"1\"/>"
404 << endl;
405 outfile << "<PDataArray type=\"Int32\" Name=\"offsets\" "
406 "NumberOfComponents=\"1\"/>"
407 << endl;
408 outfile << "<PDataArray type=\"UInt8\" Name=\"types\" "
409 "NumberOfComponents=\"1\"/>"
410 << endl;
411 outfile << "</PCells>" << endl;
412 outfile << "<PPointData Scalars=\"Material\">" << endl;
413 for (int i = 0; i < m_f->m_variables.size(); ++i)
414 {
415 outfile << "<PDataArray type=\"Float64\" Name=\"" << m_f->m_variables[i]
416 << "\"/>" << endl;
417 }
418 outfile << "</PPointData>" << endl;
419
420 for (int i = 0; i < nprocs; ++i)
421 {
422 boost::format pad("P%1$07d.vtu");
423 pad % i;
424 outfile << "<Piece Source=\"" << path << "/" << pad.str() << "\"/>"
425 << endl;
426 }
427 outfile << "</PUnstructuredGrid>" << endl;
428 outfile << "</VTKFile>" << endl;
429
430 cout << "Written file: " << filename << endl;
431}
432
433std::string OutputVtkBase::PrepareOutput(po::variables_map &vm)
434{
435 // Extract the output filename and extension
436 string filename = m_config["outfile"].as<string>();
437
438 fs::path specPath = GetPath(filename, vm);
439 fs::path fulloutname = GetFullOutName(filename, vm);
440 filename = LibUtilities::PortablePath(fulloutname);
441
442 if (m_f->m_comm->GetSpaceComm()->GetSize() != 1)
443 {
444 if (m_f->m_comm->GetSpaceComm()->TreatAsRankZero())
445 {
446 try
447 {
448 fs::create_directory(specPath);
449 }
450 catch (fs::filesystem_error &e)
451 {
452 ASSERTL0(false, "Filesystem error: " + string(e.what()));
453 }
454 cout << "Writing files to directory: " << specPath << endl;
455 }
456 m_f->m_comm->GetSpaceComm()->Block();
457 }
458 else
459 {
460 cout << "Writing: " << specPath << endl;
461 }
462 return filename;
463}
464
465} // namespace FieldUtils
466} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
FieldSharedPtr m_f
Field object.
Definition: Module.h:234
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:263
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)
virtual 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:51
virtual void v_OutputFromData(po::variables_map &vm) override
Write from data to output file.
virtual fs::path v_GetPath(std::string &filename, po::variables_map &vm) override
void WritePVtu(po::variables_map &vm)
void WriteVtkFooter(std::ostream &outfile)
virtual void v_OutputFromExp(po::variables_map &vm) override
Write from m_exp to output file.
virtual 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:198
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:991
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:317
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:49
std::string PortablePath(const boost::filesystem::path &path)
create portable path on different platforms for boost::filesystem path
Definition: FileSystem.cpp:45
std::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:190
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2