40#include <boost/format.hpp>
51#ifndef NEKTAR_USING_VTK
55 "Writes a VTU file.");
76 ofstream outfile(filename.c_str());
79 int dim = fPts->GetDim();
83 switch (fPts->GetPtsType())
89 "VTK output needs setting up for ePtsFile and ePtsLine");
95 "VTK output needs setting up for PtsPlane");
101 "VTK output needs setting up for PtsBox");
126 vector<Array<OneD, int>> ptsConn;
127 fPts->GetConnectivity(ptsConn);
129 nfields = fPts->GetNFields();
131 int nPts = fPts->GetNpoints();
133 for (i = 0; i < ptsConn.size(); ++i)
135 numBlocks += ptsConn[i].size() / nvert;
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)
146 for (j = 0; j < dim; ++j)
148 outfile <<
" " << setprecision(8) << scientific
149 << fPts->GetPointVal(j, i) <<
" ";
151 for (j = dim; j < 3; ++j)
154 outfile <<
" 0.000000";
158 outfile <<
" </DataArray>" << endl;
159 outfile <<
" </Points>" << endl;
160 outfile <<
" <Cells>" << endl;
161 outfile <<
" <DataArray type=\"Int32\" "
162 << R
"(Name="connectivity" format="ascii">)" << endl;
167 for (i = 0; i < ptsConn.size(); ++i)
169 for (j = 0; j < ptsConn[i].size(); ++j)
171 outfile << ptsConn[i][j] <<
" ";
172 if ((!(cnt % nvert)) && cnt)
174 outfile << std::endl;
180 outfile <<
" </DataArray>" << endl;
181 outfile <<
" <DataArray type=\"Int32\" "
182 << R
"(Name="offsets" format="ascii">)" << endl;
185 for (i = 0; i < numBlocks; ++i)
187 outfile << i * nvert + nvert <<
" ";
190 outfile <<
" </DataArray>" << endl;
191 outfile <<
" <DataArray type=\"UInt8\" "
192 << R
"(Name="types" format="ascii">)" << endl;
194 for (i = 0; i < numBlocks; ++i)
196 outfile << vtktype <<
" ";
199 outfile <<
" </DataArray>" << endl;
200 outfile <<
" </Cells>" << endl;
201 outfile <<
" <PointData>" << endl;
204 for (j = 0; j < nfields; ++j)
206 outfile << R
"( <DataArray type="Float64" Name=")"
207 << m_f->m_variables[j] << "\">" << endl;
209 for (i = 0; i < fPts->GetNpoints(); ++i)
211 outfile << fPts->GetPointVal(dim + j, i) <<
" ";
214 outfile <<
" </DataArray>" << endl;
217 outfile <<
" </PointData>" << endl;
218 outfile <<
" </Piece>" << endl;
221 cout <<
"Written file: " << filename << endl;
224 if ((
m_f->m_comm->GetSpaceComm()->GetRank() == 0) &&
225 (
m_f->m_comm->GetSpaceComm()->GetSize() != 1))
228 cout <<
"Written file: " << filename << endl;
239 ofstream outfile(filename.c_str());
241 int nfields =
m_f->m_variables.size();
244 m_f->m_session->LoadParameter(
"Strip_Z", nstrips, 1);
247 for (
int s = 0; s < nstrips; ++s)
250 for (i = 0; i <
m_f->m_exp[0]->GetNumElmts(); ++i)
252 m_f->m_exp[0]->WriteVtkPieceHeader(outfile, i, s);
255 for (j = 0; j < nfields; ++j)
257 m_f->m_exp[s * nfields + j]->WriteVtkPieceData(
258 outfile, i,
m_f->m_variables[j]);
260 m_f->m_exp[0]->WriteVtkPieceFooter(outfile, i);
264 if (
m_f->m_exp[0]->GetNumElmts() == 0)
270 cout <<
"Written file: " << filename << endl;
273 if ((
m_f->m_comm->GetSpaceComm()->GetRank() == 0) &&
274 (
m_f->m_comm->GetSpaceComm()->GetSize() != 1))
283 "OutputVtk can't write using only FieldData. You may need "
284 "to add a mesh XML file to your input files.");
288 [[maybe_unused]] po::variables_map &vm)
290 int nprocs =
m_f->m_comm->GetSpaceComm()->GetSize();
294 specPath = fs::path(filename);
299 int dot = filename.find_last_of(
'.');
300 string path = filename.substr(0, dot) +
"_vtu";
301 specPath = fs::path(path);
303 return fs::path(specPath);
307 po::variables_map &vm)
309 int nprocs =
m_f->m_comm->GetSpaceComm()->GetSize();
311 fs::path fulloutname;
314 fulloutname = filename;
320 pad %
m_f->m_comm->GetSpaceComm()->GetRank() %
"vtu";
323 fs::path specPath =
GetPath(filename, vm);
324 fs::path poutfile(pad.str());
325 fulloutname = specPath / poutfile;
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;
340 outfile <<
" </UnstructuredGrid>" << endl;
341 outfile <<
"</VTKFile>" << endl;
347 outfile <<
" <Piece NumberOfPoints=\"" << 0 <<
"\" NumberOfCells=\"" << 0
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;
363 outfile <<
" </DataArray>" << endl;
364 outfile <<
" <DataArray type=\"UInt8\" "
365 << R
"(Name="types" format="ascii">)" << endl;
368 outfile <<
" </DataArray>" << endl;
369 outfile <<
" </Cells>" << endl;
370 outfile <<
" <PointData>" << endl;
372 outfile <<
" </PointData>" << endl;
373 outfile <<
" </Piece>" << endl;
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";
383 ofstream outfile(filename.c_str());
385 int nprocs =
m_f->m_comm->GetSpaceComm()->GetSize();
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
395 outfile <<
"</PPoints>" << endl;
396 outfile <<
"<PCells>" << endl;
397 outfile <<
"<PDataArray type=\"Int32\" Name=\"connectivity\" "
398 "NumberOfComponents=\"1\"/>"
400 outfile <<
"<PDataArray type=\"Int32\" Name=\"offsets\" "
401 "NumberOfComponents=\"1\"/>"
403 outfile <<
"<PDataArray type=\"UInt8\" Name=\"types\" "
404 "NumberOfComponents=\"1\"/>"
406 outfile <<
"</PCells>" << endl;
407 outfile <<
"<PPointData Scalars=\"Material\">" << endl;
408 for (
int i = 0; i <
m_f->m_variables.size(); ++i)
410 outfile << R
"(<PDataArray type="Float64" Name=")" << m_f->m_variables[i]
413 outfile <<
"</PPointData>" << endl;
415 for (
int i = 0; i < nprocs; ++i)
419 outfile <<
"<Piece Source=\"" << path <<
"/" << pad.str() <<
"\"/>"
422 outfile <<
"</PUnstructuredGrid>" << endl;
423 outfile <<
"</VTKFile>" << endl;
425 cout <<
"Written file: " << filename << endl;
431 string filename =
m_config[
"outfile"].as<
string>();
433 ASSERTL0(filename !=
"",
"Legacy VTK output requires a filename.");
435 fs::path specPath =
GetPath(filename, vm);
439 if (
m_f->m_comm->GetSpaceComm()->GetSize() != 1)
441 if (
m_f->m_comm->GetSpaceComm()->TreatAsRankZero())
445 fs::create_directory(specPath);
447 catch (fs::filesystem_error &e)
449 ASSERTL0(
false,
"Filesystem error: " +
string(e.what()));
451 cout <<
"Writing files to directory: " << specPath << endl;
453 m_f->m_comm->GetSpaceComm()->Block();
457 cout <<
"Writing: " << specPath << endl;
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
FieldSharedPtr m_f
Field object.
std::map< std::string, ConfigOption > m_config
List of configuration values.
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 ModuleKey m_className
static std::shared_ptr< Module > create(FieldSharedPtr f)
Creates an instance of this class.
void v_OutputFromData(po::variables_map &vm) override
Write from data to output file.
OutputVtkBase(FieldSharedPtr f)
~OutputVtkBase() override
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.
std::shared_ptr< Field > FieldSharedPtr
std::pair< ModuleType, std::string > ModuleKey
ModuleFactory & GetModuleFactory()
static std::string PortablePath(const fs::path &path)
create portable path on different platforms for std::filesystem path.
std::shared_ptr< PtsField > PtsFieldSharedPtr