ensightSurfaceWriterTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2014 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "IOmanip.H"
27 #include "IFstream.H"
28 #include "OFstream.H"
29 #include "OSspecific.H"
30 #include "ensightPartFaces.H"
31 #include "ensightPTraits.H"
32 #include "OStringStream.H"
33 #include "regExp.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 template<class Type>
39 (
40  const fileName& outputDir,
41  const fileName& surfaceName,
42  const pointField& points,
43  const faceList& faces,
44  const word& fieldName,
45  const Field<Type>& values,
46  const bool isNodeValues,
47  const bool verbose
48 ) const
49 {
50  if (!isDir(outputDir/fieldName))
51  {
52  mkDir(outputDir/fieldName);
53  }
54 
55  // const scalar timeValue = Foam::name(this->mesh().time().timeValue());
56  const scalar timeValue = 0.0;
57 
58  OFstream osCase(outputDir/fieldName/surfaceName + ".case");
59  ensightGeoFile osGeom
60  (
61  outputDir/fieldName/surfaceName + ".0000.mesh",
62  writeFormat_
63  );
64  ensightFile osField
65  (
66  outputDir/fieldName/surfaceName + ".0000." + fieldName,
67  writeFormat_
68  );
69 
70  if (verbose)
71  {
72  Info<< "Writing case file to " << osCase.name() << endl;
73  }
74 
75  osCase
76  << "FORMAT" << nl
77  << "type: ensight gold" << nl
78  << nl
79  << "GEOMETRY" << nl
80  << "model: 1 " << osGeom.name().name() << nl
81  << nl
82  << "VARIABLE" << nl
83  << ensightPTraits<Type>::typeName << " per "
84  << word(isNodeValues ? "node:" : "element:") << setw(10) << 1
85  << " " << fieldName
86  << " " << surfaceName.c_str() << ".****." << fieldName << nl
87  << nl
88  << "TIME" << nl
89  << "time set: 1" << nl
90  << "number of steps: 1" << nl
91  << "filename start number: 0" << nl
92  << "filename increment: 1" << nl
93  << "time values:" << nl
94  << timeValue << nl
95  << nl;
96 
97  ensightPartFaces ensPart(0, osGeom.name().name(), points, faces, true);
98  osGeom << ensPart;
99 
100  // Write field
102  ensPart.writeField(osField, values, isNodeValues);
103 
104  return osCase.name();
105 }
106 
107 
108 template<class Type>
110 (
111  const fileName& outputDir,
112  const fileName& surfaceName,
113  const pointField& points,
114  const faceList& faces,
115  const word& fieldName,
116  const Field<Type>& values,
117  const bool isNodeValues,
118  const bool verbose
119 ) const
120 {
121 
122  const fileName baseDir = outputDir.path()/surfaceName;
123  const fileName timeDir = outputDir.name();
124 
125  if (!isDir(baseDir))
126  {
127  mkDir(baseDir);
128  }
129 
130  const fileName meshFile(baseDir/surfaceName + ".0000.mesh");
131  const scalar timeValue = readScalar(IStringStream(timeDir)());
132  label timeIndex = 0;
133 
134 
135  // Do case file
136  {
138  scalarList times;
139  bool stateChanged = false;
140 
141  if (isFile(baseDir/"fieldsDict"))
142  {
143  IFstream is(baseDir/"fieldsDict");
144  if (is.good() && dict.read(is))
145  {
146  dict.lookup("times") >> times;
147  const scalar timeValue = readScalar(IStringStream(timeDir)());
148  label index = findLower(times, timeValue);
149  timeIndex = index+1;
150  }
151  }
152 
153 
154  // Update stored times list
155  times.setSize(timeIndex+1, -1);
156 
157  if (times[timeIndex] != timeValue)
158  {
159  stateChanged = true;
160  }
161  times[timeIndex] = timeValue;
162 
163 
164  // Add my information to dictionary
165  {
166  dict.set("times", times);
167  if (dict.found("fields"))
168  {
169  dictionary& fieldsDict = dict.subDict("fields");
170  if (!fieldsDict.found(fieldName))
171  {
172  dictionary fieldDict;
173  fieldDict.set("type", ensightPTraits<Type>::typeName);
174  fieldsDict.set(fieldName, fieldDict);
175 
176  stateChanged = true;
177  }
178  }
179  else
180  {
181  dictionary fieldDict;
182  fieldDict.set("type", ensightPTraits<Type>::typeName);
183 
184  dictionary fieldsDict;
185  fieldsDict.set(fieldName, fieldDict);
186 
187  dict.set("fields", fieldsDict);
188 
189  stateChanged = true;
190  }
191  }
192 
193 
194  if (stateChanged)
195  {
196  if (verbose)
197  {
198  Info<< "Writing state file to fieldsDict" << endl;
199  }
200  OFstream os(baseDir/"fieldsDict");
201  os << dict;
202 
203 
204  OFstream osCase(baseDir/surfaceName + ".case");
205 
206  if (verbose)
207  {
208  Info<< "Writing case file to " << osCase.name() << endl;
209  }
210 
211  osCase
212  << "FORMAT" << nl
213  << "type: ensight gold" << nl
214  << nl
215  << "GEOMETRY" << nl
216  << "model: 1 " << meshFile.name() << nl
217  << nl
218  << "VARIABLE" << nl;
219  const dictionary& fieldsDict = dict.subDict("fields");
220  forAllConstIter(dictionary, fieldsDict, iter)
221  {
222  const word& fieldName = iter().keyword();
223  const word fieldType(iter().dict().lookup("type"));
224 
225  osCase
226  << fieldType << " per "
227  << word(isNodeValues ? "node:" : "element:")
228  << setw(10) << 1
229  << setw(15) << fieldName
230  << " " << surfaceName.c_str() << ".****." << fieldName
231  << nl;
232  }
233  osCase << nl;
234 
235  osCase
236  << "TIME" << nl
237  << "time set: 1" << nl
238  << "number of steps: " << timeIndex+1 << nl
239  << "filename start number: 0" << nl
240  << "filename increment: 1" << nl
241  << "time values:" << nl;
242  forAll(times, timeI)
243  {
244  osCase << setw(12) << times[timeI] << " ";
245 
246  if (timeI != 0 && (timeI % 6) == 0)
247  {
248  osCase << nl;
249  }
250  }
251  osCase << nl;
252  }
253  }
254 
255 
256  // Write geometry
257  ensightPartFaces ensPart(0, meshFile.name(), points, faces, true);
258  if (!exists(meshFile))
259  {
260  if (verbose)
261  {
262  Info<< "Writing mesh file to " << meshFile.name() << endl;
263  }
264  ensightGeoFile osGeom(meshFile, writeFormat_);
265  osGeom << ensPart;
266  }
267 
268 
269  // Get string representation
270  string timeString;
271  {
272  OStringStream os;
273  os.stdStream().fill('0');
274  os << setw(4) << timeIndex;
275  timeString = os.str();
276  }
277 
278  // Write field
279  ensightFile osField
280  (
281  baseDir/surfaceName + "." + timeString + "." + fieldName,
282  writeFormat_
283  );
284  if (verbose)
285  {
286  Info<< "Writing field file to " << osField.name() << endl;
287  }
289  ensPart.writeField(osField, values, isNodeValues);
290 
291  return baseDir/surfaceName + ".case";
292 }
293 
294 
295 template<class Type>
297 (
298  const fileName& outputDir,
299  const fileName& surfaceName,
300  const pointField& points,
301  const faceList& faces,
302  const word& fieldName,
303  const Field<Type>& values,
304  const bool isNodeValues,
305  const bool verbose
306 ) const
307 {
308  if (collateTimes_)
309  {
310  return writeCollated
311  (
312  outputDir,
313  surfaceName,
314  points,
315  faces,
316  fieldName,
317  values,
318  isNodeValues,
319  verbose
320  );
321  }
322  else
323  {
324  return writeUncollated
325  (
326  outputDir,
327  surfaceName,
328  points,
329  faces,
330  fieldName,
331  values,
332  isNodeValues,
333  verbose
334  );
335  }
336 }
337 
338 
339 // ************************************************************************* //
Foam::ensightSurfaceWriter::writeCollated
fileName writeCollated(const fileName &outputDir, const fileName &surfaceName, const pointField &points, const faceList &faces, const word &fieldName, const Field< Type > &values, const bool isNodeValues, const bool verbose) const
Templated write operation - one file per timestep.
Foam::ensightSurfaceWriter::writeTemplate
fileName writeTemplate(const fileName &outputDir, const fileName &surfaceName, const pointField &points, const faceList &faces, const word &fieldName, const Field< Type > &values, const bool isNodeValues, const bool verbose) const
Templated write operation.
OSspecific.H
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::exists
bool exists(const fileName &, const bool checkGzip=true)
Does the name exist (as DIRECTORY or FILE) in the file system?
Definition: POSIX.C:608
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::IFstream
Input from file stream.
Definition: IFstream.H:81
ensightPartFaces.H
Foam::ensightSurfaceWriter::writeUncollated
fileName writeUncollated(const fileName &outputDir, const fileName &surfaceName, const pointField &points, const faceList &faces, const word &fieldName, const Field< Type > &values, const bool isNodeValues, const bool verbose) const
Templated write operation - all time steps in single file.
Foam::fileName::path
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:293
Foam::OSstream::stdStream
virtual ostream & stdStream()
Access to underlying std::ostream.
Definition: OSstream.H:178
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::OStringStream::str
string str() const
Return the string.
Definition: OStringStream.H:107
Foam::fileName::name
word name() const
Return file name (part beyond last /)
Definition: fileName.C:212
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
OFstream.H
Foam::ensightPTraits
Conversion of OpenFOAM pTraits into the Ensight equivalent.
Definition: ensightPTraits.H:48
Foam::ensightGeoFile
Specialized Ensight output with extra geometry file header.
Definition: ensightGeoFile.H:45
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
ensightPTraits.H
OStringStream.H
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::findLower
label findLower(const ListType &, typename ListType::const_reference, const label start, const BinaryOp &bop)
Find last element < given value in sorted list and return index,.
IFstream.H
Foam::ensightFile::writeKeyword
virtual Ostream & writeKeyword(const string &key)
Write element keyword with trailing newline, optionally with undef.
Definition: ensightFile.C:275
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::ensightFile
Ensight output with specialized write() for strings, integers and floats. Correctly handles binary wr...
Definition: ensightFile.H:47
Foam::IStringStream
Input from memory buffer stream.
Definition: IStringStream.H:49
Foam::isFile
bool isFile(const fileName &, const bool checkGzip=true)
Does the name exist as a FILE in the file system?
Definition: POSIX.C:622
Foam::isDir
bool isDir(const fileName &)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:615
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::OFstream
Output to file stream.
Definition: OFstream.H:81
readScalar
#define readScalar
Definition: doubleScalar.C:38
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::OStringStream
Output to memory buffer stream.
Definition: OStringStream.H:49
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::ensightPartFaces
An implementation of ensightPart to hold volume mesh faces.
Definition: ensightPartFaces.H:48
regExp.H
Foam::OFstream::name
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:118
Foam::IOstream::good
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
Foam::ensightPart::writeField
void writeField(ensightFile &, const Field< Type > &, const bool perNode=false) const
Write generalized field components.
Definition: ensightPartTemplates.C:35
Foam::mkDir
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:419
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress
Foam::dictionary::set
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:856