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  if(m_f->m_comm->TreatAsRankZero())
78  {
79  cout << "OutputVtk: Writing file..." << endl;
80  }
81  }
82 
83  // Extract the output filename and extension
84  string filename = m_config["outfile"].as<string>();
85  string path;
86 
87  // amend for parallel output if required
88  if(m_f->m_session->GetComm()->GetSize() != 1)
89  {
90  int dot = filename.find_last_of('.');
91  string ext = filename.substr(dot,filename.length()-dot);
92  string start = filename.substr(0,dot);
93  path = start + "_vtu";
94 
95  boost::format pad("P%1$07d.vtu");
96  pad % m_f->m_session->GetComm()->GetRank();
97  filename = pad.str();
98 
99  fs::path poutfile(filename.c_str());
100  fs::path specPath(path.c_str());
101 
102  if(m_f->m_comm->TreatAsRankZero())
103  {
104  try
105  {
106  fs::create_directory(specPath);
107  }
108  catch (fs::filesystem_error& e)
109  {
110  ASSERTL0(false, "Filesystem error: " + string(e.what()));
111  }
112  cout << "Writing files to directory: " << specPath << endl;
113  }
114 
115  fs::path fulloutname = specPath / poutfile;
116  filename = LibUtilities::PortablePath(fulloutname);
117  m_f->m_comm->Block();
118  }
119  else
120  {
121  fs::path specPath(filename.c_str());
122  cout << "Writing: " << specPath << endl;
123  filename = LibUtilities::PortablePath(specPath);
124  }
125 
126  // Write solution.
127  ofstream outfile(filename.c_str());
128  m_f->m_exp[0]->WriteVtkHeader(outfile);
129  int nfields = 0;
130  int dim = 0;
131 
132  vector<string> fieldname;
133  if(fPts == LibUtilities::NullPtsField) // standard output in collapsed coordinates
134  {
135  dim = m_f->m_exp[0]->GetExp(0)->GetCoordim();
136 
137  int nstrips;
138  if (m_f->m_fielddef.size() == 0)
139  {
140  nfields = 0;
141  }
142  else
143  {
144  nfields = m_f->m_fielddef[0]->m_fields.size();
145  }
146  m_f->m_session->LoadParameter("Strip_Z", nstrips, 1);
147 
148  // Homogeneous strip variant
149  for(int s = 0; s < nstrips; ++s)
150  {
151  // For each field write out field data for each expansion.
152  for (i = 0; i < m_f->m_exp[0]->GetNumElmts(); ++i)
153  {
154  m_f->m_exp[0]->WriteVtkPieceHeader(outfile,i,s);
155 
156  // For this expansion write out each field.
157  for (j = 0; j < nfields; ++j)
158  {
159  m_f->m_exp[s*nfields+j]->WriteVtkPieceData(
160  outfile, i, m_f->m_fielddef[0]->m_fields[j]);
161  }
162  m_f->m_exp[0]->WriteVtkPieceFooter(outfile, i);
163  }
164  }
165 
166  // save field names for parallel output
167  for(i = 0; i < nfields; ++i)
168  {
169  fieldname.push_back(m_f->m_fielddef[0]->m_fields[i]);
170  }
171  }
172  else // write out data stored in fPts (for example if equispaced output is called).
173  {
174  int i = 0;
175  int j = 0;
176 
177  dim = fPts->GetDim();
178 
179  if(fPts->GetNpoints() == 0)
180  {
181  return;
182  }
183 
184  int nvert, vtktype;
185  switch(fPts->GetPtsType())
186  {
189  {
190  ASSERTL0(false,"VTK output needs settig up for PtsFile or Pts Line");
191  break;
192  }
194  {
195  ASSERTL0(false,"VTK output needs settig up for PtsPlane");
196  break;
197  }
199  {
200  nvert = 3;
201  vtktype = 5;
202  break;
203  }
205  {
206  nvert = 4;
207  vtktype = 10;
208  break;
209  }
210  default:
211  ASSERTL0(false, "ptsType not supported yet.");
212  }
213 
214  vector<Array<OneD, int> > ptsConn;
215  fPts->GetConnectivity(ptsConn);
216 
217  nfields = fPts->GetNFields();
218 
219  int nPts = fPts->GetNpoints();
220  int numBlocks = 0;
221  for(i = 0; i < ptsConn.size(); ++i)
222  {
223  numBlocks += ptsConn[i].num_elements()/nvert;
224  }
225 
226  // write out pieces of data.
227  outfile << " <Piece NumberOfPoints=\""
228  << nPts << "\" NumberOfCells=\""
229  << numBlocks << "\">" << endl;
230  outfile << " <Points>" << endl;
231  outfile << " <DataArray type=\"Float64\" "
232  << "NumberOfComponents=\""<<3<<"\" format=\"ascii\">" << endl;
233  for(i = 0; i < nPts; ++i)
234  {
235  for(j = 0; j < dim; ++j)
236  {
237  outfile << " " << setprecision(8) << scientific
238  << fPts->GetPointVal(j, i) << " ";
239  }
240  for(j = dim; j < 3; ++j) // pack to 3D since paraview does not seem to handle 2D
241  {
242  outfile << " 0.000000" ;
243  }
244  outfile << endl;
245  }
246  outfile << " </DataArray>" << endl;
247  outfile << " </Points>" << endl;
248  outfile << " <Cells>" << endl;
249  outfile << " <DataArray type=\"Int32\" "
250  << "Name=\"connectivity\" format=\"ascii\">" << endl;
251 
252  // dump connectivity data if it exists
253  outfile << " ";
254  int cnt = 1;
255  for(i = 0; i < ptsConn.size();++i)
256  {
257  for(j = 0; j < ptsConn[i].num_elements(); ++j)
258  {
259  outfile << ptsConn[i][j] << " ";
260  if( (!(cnt % nvert)) && cnt )
261  {
262  outfile << std::endl;
263  outfile << " ";
264  }
265  cnt ++;
266  }
267  }
268  outfile << " </DataArray>" << endl;
269  outfile << " <DataArray type=\"Int32\" "
270  << "Name=\"offsets\" format=\"ascii\">" << endl;
271 
272  outfile << " ";
273  for (i = 0; i < numBlocks; ++i)
274  {
275  outfile << i*nvert+nvert << " ";
276  }
277  outfile << endl;
278  outfile << " </DataArray>" << endl;
279  outfile << " <DataArray type=\"UInt8\" "
280  << "Name=\"types\" format=\"ascii\">" << endl;
281  outfile << " ";
282  for (i = 0; i < numBlocks; ++i)
283  {
284  outfile << vtktype <<" ";
285  }
286  outfile << endl;
287  outfile << " </DataArray>" << endl;
288  outfile << " </Cells>" << endl;
289  outfile << " <PointData>" << endl;
290 
291  // printing the fields
292  for(j = 0; j < nfields; ++j)
293  {
294  fieldname.push_back(fPts->GetFieldName(j));
295  outfile << " <DataArray type=\"Float64\" Name=\""
296  << fPts->GetFieldName(j) << "\">" << endl;
297  outfile << " ";
298  for(i = 0; i < fPts->GetNpoints(); ++i)
299  {
300  outfile << fPts->GetPointVal(dim+j, i) << " ";
301  }
302  outfile << endl;
303  outfile << " </DataArray>" << endl;
304  }
305 
306  outfile << " </PointData>" << endl;
307  outfile << " </Piece>" << endl;
308  }
309 
310  m_f->m_exp[0]->WriteVtkFooter(outfile);
311  cout << "Written file: " << filename << endl;
312 
313  // output parallel outline info if necessary
314  if(m_f->m_comm->GetRank() == 0)
315  {
316  ASSERTL1(fieldname.size() == nfields, "fieldname not the same size as nfields");
317  int nprocs = m_f->m_comm->GetSize();
318  if(nprocs != 1)
319  {
320  filename = m_config["outfile"].as<string>();
321  int dot = filename.find_last_of('.');
322  string body = filename.substr(0,dot);
323  filename = body + ".pvtu";
324 
325  ofstream outfile(filename.c_str());
326 
327  outfile << "<?xml version=\"1.0\"?>" << endl;
328  outfile << "<VTKFile type=\"PUnstructuredGrid\" version=\"0.1\" "
329  << "byte_order=\"LittleEndian\">" << endl;
330  outfile << "<PUnstructuredGrid GhostLevel=\"0\">" << endl;
331  outfile << "<PPoints> " << endl;
332  outfile << "<PDataArray type=\"Float64\" NumberOfComponents=\""
333  << 3 << "\"/> " << endl;
334  outfile << "</PPoints>" << endl;
335  outfile << "<PCells>" << endl;
336  outfile << "<PDataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\"/>" << endl;
337  outfile << "<PDataArray type=\"Int32\" Name=\"offsets\" NumberOfComponents=\"1\"/>" << endl;
338  outfile << "<PDataArray type=\"UInt8\" Name=\"types\" NumberOfComponents=\"1\"/>" << endl;
339  outfile << "</PCells>" << endl;
340  outfile << "<PPointData Scalars=\"Material\">" << endl;
341  for(int i = 0; i < nfields; ++i)
342  {
343  outfile << "<PDataArray type=\"Float64\" Name=\"" <<
344  fieldname[i] << "\"/>" << endl;
345  }
346  outfile << "</PPointData>" << endl;
347 
348  for(int i = 0; i < nprocs; ++i)
349  {
350  boost::format pad("P%1$07d.vtu");
351  pad % i;
352  outfile << "<Piece Source=\"" << path << "/" << pad.str() << "\"/>" <<endl;
353  }
354  outfile << "</PUnstructuredGrid>" << endl;
355  outfile << "</VTKFile>" << endl;
356  cout << "Written file: " << filename << endl;
357  }
358  }
359 }
360 
361 }
362 }
#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:698
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