Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
FieldConvert/OutputModules/OutputVtk.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: OutputVtk.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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: VTK file format output.
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 
36 #include <set>
37 #include <string>
38 using namespace std;
39 
40 #include <iomanip>
41 #include "OutputVtk.h"
42 #include <boost/format.hpp>
44 
45 namespace Nektar
46 {
47 namespace Utilities
48 {
49 
50 ModuleKey OutputVtk::m_className =
52  ModuleKey(eOutputModule, "vtu"), OutputVtk::create,
53  "Writes a VTU file.");
54 
55 OutputVtk::OutputVtk(FieldSharedPtr f) : OutputModule(f)
56 {
57  m_requireEquiSpaced = true;
58 }
59 
61 {
62 }
63 
64 void OutputVtk::Process(po::variables_map &vm)
65 {
66  LibUtilities::PtsFieldSharedPtr fPts = m_f->m_fieldPts;
67 
68  // Do nothing if no expansion defined
69  if (fPts == LibUtilities::NullPtsField && !m_f->m_exp.size())
70  {
71  return;
72  }
73 
74  int i, j;
75  if (m_f->m_verbose)
76  {
77  cout << "OutputVtk: Writing file..." << endl;
78  }
79 
80  // Extract the output filename and extension
81  string filename = m_config["outfile"].as<string>();
82  string path;
83 
84  // amend for parallel output if required
85  if(m_f->m_session->GetComm()->GetSize() != 1)
86  {
87  int dot = filename.find_last_of('.');
88  string ext = filename.substr(dot,filename.length()-dot);
89  string start = filename.substr(0,dot);
90  path = start + "_vtu";
91 
92  boost::format pad("P%1$07d.vtu");
93  pad % m_f->m_session->GetComm()->GetRank();
94  filename = pad.str();
95 
96  fs::path poutfile(filename.c_str());
97  fs::path specPath(path.c_str());
98 
99  if(m_f->m_comm->GetRank() == 0)
100  {
101  try
102  {
103  fs::create_directory(specPath);
104  }
105  catch (fs::filesystem_error& e)
106  {
107  ASSERTL0(false, "Filesystem error: " + string(e.what()));
108  }
109  cout << "Writing files to directory: " << specPath << endl;
110  }
111 
112  fs::path fulloutname = specPath / poutfile;
113  filename = LibUtilities::PortablePath(fulloutname);
114  m_f->m_comm->Block();
115  }
116  else
117  {
118  fs::path specPath(filename.c_str());
119  cout << "Writing: " << specPath << endl;
120  filename = LibUtilities::PortablePath(specPath);
121  }
122 
123  // Write solution.
124  ofstream outfile(filename.c_str());
125  m_f->m_exp[0]->WriteVtkHeader(outfile);
126  int nfields = 0;
127  int dim = 0;
128 
129  vector<string> fieldname;
130  if(fPts == LibUtilities::NullPtsField) // standard output in collapsed coordinates
131  {
132  dim = m_f->m_exp[0]->GetExp(0)->GetCoordim();
133 
134  int nstrips;
135  if (m_f->m_fielddef.size() == 0)
136  {
137  nfields = 0;
138  }
139  else
140  {
141  nfields = m_f->m_fielddef[0]->m_fields.size();
142  }
143  m_f->m_session->LoadParameter("Strip_Z", nstrips, 1);
144 
145  // Homogeneous strip variant
146  for(int s = 0; s < nstrips; ++s)
147  {
148  // For each field write out field data for each expansion.
149  for (i = 0; i < m_f->m_exp[0]->GetNumElmts(); ++i)
150  {
151  m_f->m_exp[0]->WriteVtkPieceHeader(outfile,i,s);
152 
153  // For this expansion write out each field.
154  for (j = 0; j < nfields; ++j)
155  {
156  m_f->m_exp[s*nfields+j]->WriteVtkPieceData(
157  outfile, i, m_f->m_fielddef[0]->m_fields[j]);
158  }
159  m_f->m_exp[0]->WriteVtkPieceFooter(outfile, i);
160  }
161  }
162 
163  // save field names for parallel output
164  for(i = 0; i < nfields; ++i)
165  {
166  fieldname.push_back(m_f->m_fielddef[0]->m_fields[i]);
167  }
168  }
169  else // write out data stored in fPts (for example if equispaced output is called).
170  {
171  int i = 0;
172  int j = 0;
173 
174  dim = fPts->GetDim();
175 
176  if(fPts->GetNpoints() == 0)
177  {
178  return;
179  }
180 
181  int nvert, vtktype;
182  switch(fPts->GetPtsType())
183  {
186  {
187  ASSERTL0(false,"VTK output needs settig up for PtsFile or Pts Line");
188  break;
189  }
191  {
192  ASSERTL0(false,"VTK output needs settig up for PtsPlane");
193  break;
194  }
196  {
197  nvert = 3;
198  vtktype = 5;
199  break;
200  }
202  {
203  nvert = 4;
204  vtktype = 10;
205  break;
206  }
207  default:
208  ASSERTL0(false, "ptsType not supported yet.");
209  }
210 
211  vector<Array<OneD, int> > ptsConn;
212  fPts->GetConnectivity(ptsConn);
213 
214  nfields = fPts->GetNFields();
215 
216  int nPts = fPts->GetNpoints();
217  int numBlocks = 0;
218  for(i = 0; i < ptsConn.size(); ++i)
219  {
220  numBlocks += ptsConn[i].num_elements()/nvert;
221  }
222 
223  // write out pieces of data.
224  outfile << " <Piece NumberOfPoints=\""
225  << nPts << "\" NumberOfCells=\""
226  << numBlocks << "\">" << endl;
227  outfile << " <Points>" << endl;
228  outfile << " <DataArray type=\"Float64\" "
229  << "NumberOfComponents=\""<<3<<"\" format=\"ascii\">" << endl;
230  for(i = 0; i < nPts; ++i)
231  {
232  for(j = 0; j < dim; ++j)
233  {
234  outfile << " " << setprecision(8) << scientific
235  << fPts->GetPointVal(j, i) << " ";
236  }
237  for(j = dim; j < 3; ++j) // pack to 3D since paraview does not seem to handle 2D
238  {
239  outfile << " 0.000000" ;
240  }
241  outfile << endl;
242  }
243  outfile << " </DataArray>" << endl;
244  outfile << " </Points>" << endl;
245  outfile << " <Cells>" << endl;
246  outfile << " <DataArray type=\"Int32\" "
247  << "Name=\"connectivity\" format=\"ascii\">" << endl;
248 
249  // dump connectivity data if it exists
250  outfile << " ";
251  int cnt = 1;
252  for(i = 0; i < ptsConn.size();++i)
253  {
254  for(j = 0; j < ptsConn[i].num_elements(); ++j)
255  {
256  outfile << ptsConn[i][j] << " ";
257  if( (!(cnt % nvert)) && cnt )
258  {
259  outfile << std::endl;
260  outfile << " ";
261  }
262  cnt ++;
263  }
264  }
265  outfile << " </DataArray>" << endl;
266  outfile << " <DataArray type=\"Int32\" "
267  << "Name=\"offsets\" format=\"ascii\">" << endl;
268 
269  outfile << " ";
270  for (i = 0; i < numBlocks; ++i)
271  {
272  outfile << i*nvert+nvert << " ";
273  }
274  outfile << endl;
275  outfile << " </DataArray>" << endl;
276  outfile << " <DataArray type=\"UInt8\" "
277  << "Name=\"types\" format=\"ascii\">" << endl;
278  outfile << " ";
279  for (i = 0; i < numBlocks; ++i)
280  {
281  outfile << vtktype <<" ";
282  }
283  outfile << endl;
284  outfile << " </DataArray>" << endl;
285  outfile << " </Cells>" << endl;
286  outfile << " <PointData>" << endl;
287 
288  // printing the fields
289  for(j = 0; j < nfields; ++j)
290  {
291  fieldname.push_back(fPts->GetFieldName(j));
292  outfile << " <DataArray type=\"Float64\" Name=\""
293  << fPts->GetFieldName(j) << "\">" << endl;
294  outfile << " ";
295  for(i = 0; i < fPts->GetNpoints(); ++i)
296  {
297  outfile << fPts->GetPointVal(dim+j, i) << " ";
298  }
299  outfile << endl;
300  outfile << " </DataArray>" << endl;
301  }
302 
303  outfile << " </PointData>" << endl;
304  outfile << " </Piece>" << endl;
305  }
306 
307  m_f->m_exp[0]->WriteVtkFooter(outfile);
308  cout << "Written file: " << filename << endl;
309 
310  // output parallel outline info if necessary
311  if(m_f->m_comm->GetRank() == 0)
312  {
313  ASSERTL1(fieldname.size() == nfields, "fieldname not the same size as nfields");
314  int nprocs = m_f->m_comm->GetSize();
315  if(nprocs != 1)
316  {
317  filename = m_config["outfile"].as<string>();
318  int dot = filename.find_last_of('.');
319  string body = filename.substr(0,dot);
320  filename = body + ".pvtu";
321 
322  ofstream outfile(filename.c_str());
323 
324  outfile << "<?xml version=\"1.0\"?>" << endl;
325  outfile << "<VTKFile type=\"PUnstructuredGrid\" version=\"0.1\" "
326  << "byte_order=\"LittleEndian\">" << endl;
327  outfile << "<PUnstructuredGrid GhostLevel=\"0\">" << endl;
328  outfile << "<PPoints> " << endl;
329  outfile << "<PDataArray type=\"Float64\" NumberOfComponents=\""
330  << 3 << "\"/> " << endl;
331  outfile << "</PPoints>" << endl;
332  outfile << "<PCells>" << endl;
333  outfile << "<PDataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\"/>" << endl;
334  outfile << "<PDataArray type=\"Int32\" Name=\"offsets\" NumberOfComponents=\"1\"/>" << endl;
335  outfile << "<PDataArray type=\"UInt8\" Name=\"types\" NumberOfComponents=\"1\"/>" << endl;
336  outfile << "</PCells>" << endl;
337  outfile << "<PPointData Scalars=\"Material\">" << endl;
338  for(int i = 0; i < nfields; ++i)
339  {
340  outfile << "<PDataArray type=\"Float64\" Name=\"" <<
341  fieldname[i] << "\"/>" << endl;
342  }
343  outfile << "</PPointData>" << endl;
344 
345  for(int i = 0; i < nprocs; ++i)
346  {
347  boost::format pad("P%1$07d.vtu");
348  pad % i;
349  outfile << "<Piece Source=\"" << path << "/" << pad.str() << "\"/>" <<endl;
350  }
351  outfile << "</PUnstructuredGrid>" << endl;
352  outfile << "</VTKFile>" << endl;
353  cout << "Written file: " << filename << endl;
354  }
355  }
356 }
357 
358 }
359 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
pair< ModuleType, string > ModuleKey
Abstract base class for output modules.
map< string, ConfigOption > m_config
List of configuration values.
STL namespace.
FieldSharedPtr m_f
Field object.
boost::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:262
virtual void Process()
Write mesh to output file.
std::string PortablePath(const boost::filesystem::path &path)
create portable path on different platforms for boost::filesystem path
Definition: FileSystem.cpp:41
boost::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:695
static PtsFieldSharedPtr NullPtsField
Definition: PtsField.h:263
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
ModuleFactory & GetModuleFactory()
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215