ensightOutputFunctions.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
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 "ensightOutputFunctions.H"
27 
28 #include "passiveParticle.H"
29 #include "IOField.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 
33 #include "OFstream.H"
34 #include "IOmanip.H"
35 
36 
37 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
38 
40 (
41  OFstream& caseFile,
42  const string& ensightType,
43  const word& fieldName,
44  const fileName& dataMask,
45  const fileName& local,
46  const label cloudNo,
47  const label timeSet
48 )
49 {
50  caseFile.setf(ios_base::left);
51 
52  fileName dirName(dataMask);
53  if (local.size())
54  {
55  dirName = dirName/local;
56  }
57 
58  if (cloudNo >= 0)
59  {
60  label ts = 1;
61  if (timeSet > ts)
62  {
63  ts = timeSet;
64  }
65 
66  // prefix variables with 'c' (cloud)
67  caseFile
68  << ensightType.c_str()
69  << " per measured node: " << ts << " "
70  << setw(15)
71  << ("c" + Foam::name(cloudNo) + fieldName).c_str()
72  << " "
73  << (dirName/fieldName).c_str()
74  << nl;
75  }
76  else
77  {
78  caseFile
79  << ensightType.c_str()
80  << " per element: "
81  << setw(15) << fieldName
82  << " "
83  << (dirName/fieldName).c_str()
84  << nl;
85  }
86 }
87 
88 
90 (
91  const polyMesh& mesh,
92  const fileName& dataDir,
93  const fileName& subDir,
94  const word& cloudName,
95  IOstream::streamFormat format
96 )
97 {
98  Cloud<passiveParticle> parcels(mesh, cloudName, false);
99 
100  fileName cloudDir = subDir/cloud::prefix/cloudName;
101  fileName postFileName = cloudDir/"positions";
102 
103  // the ITER/lagrangian subdirectory must exist
104  mkDir(dataDir/cloudDir);
105  ensightFile os(dataDir/postFileName, format);
106 
107  // tag binary format (just like geometry files)
108  os.writeBinaryHeader();
109  os.write(postFileName);
110  os.newline();
111  os.write("particle coordinates");
112  os.newline();
113  os.write(parcels.size(), 8); // unusual width
114  os.newline();
115 
116  // binary write is Ensight6 - first ids, then positions
117  if (format == IOstream::BINARY)
118  {
119  forAll(parcels, i)
120  {
121  os.write(i+1);
122  }
123 
124  forAllConstIter(Cloud<passiveParticle>, parcels, elmnt)
125  {
126  const vector& p = elmnt().position();
127 
128  os.write(p.x());
129  os.write(p.y());
130  os.write(p.z());
131  }
132  }
133  else
134  {
135  label nParcels = 0;
136 
137  forAllConstIter(Cloud<passiveParticle>, parcels, elmnt)
138  {
139  const vector& p = elmnt().position();
140 
141  os.write(++nParcels, 8); // unusual width
142  os.write(p.x());
143  os.write(p.y());
144  os.write(p.z());
145  os.newline();
146  }
147  }
148 }
149 
150 
151 
152 template<class Type>
154 (
155  const IOobject& fieldObject,
156  const fileName& dataDir,
157  const fileName& subDir,
158  const word& cloudName,
159  IOstream::streamFormat format
160 )
161 {
162  Info<< " " << fieldObject.name() << flush;
163 
164  fileName cloudDir = subDir/cloud::prefix/cloudName;
165  fileName postFileName = cloudDir/fieldObject.name();
166 
167  string title =
168  postFileName + " with " + pTraits<Type>::typeName + " values";
169 
170  ensightFile os(dataDir/postFileName, format);
171  os.write(title);
172  os.newline();
173 
174  IOField<Type> field(fieldObject);
175 
176  // 6 values per line
177  label count = 0;
178 
179  forAll(field, i)
180  {
181  Type val = field[i];
182 
183  if (mag(val) < 1.0e-90)
184  {
185  val = Zero;
186  }
187 
188  for (direction cmpt=0; cmpt < pTraits<Type>::nComponents; cmpt++)
189  {
190  os.write( component(val, cmpt) );
191  }
192 
193  count += pTraits<Type>::nComponents;
194 
195  if (count % 6 == 0)
196  {
197  os.newline();
198  }
199  }
200 
201  // add final newline if required
202  if (count % 6)
203  {
204  os.newline();
205  }
206 }
207 
208 
209 template<class Type>
211 (
212  const ensightParts& partsList,
213  const IOobject& fieldObject,
214  const fvMesh& mesh,
215  const fileName& dataDir,
216  const fileName& subDir,
217  IOstream::streamFormat format
218 )
219 {
220  Info<< " " << fieldObject.name() << flush;
221 
222  fileName postFileName = subDir/fieldObject.name();
223 
224  ensightFile os(dataDir/postFileName, format);
225  os.write(postFileName);
226  os.newline();
227 
228  // ie, volField<Type>
229  partsList.writeField
230  (
231  os,
232  GeometricField<Type, fvPatchField, volMesh>
233  (
234  fieldObject,
235  mesh
236  )
237  );
238 }
239 
240 
241 // ************************************************************************* //
Foam::surfaceFields.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
uint8_t direction
Definition: direction.H:45
void ensightParticlePositions(const polyMesh &mesh, const fileName &dataDir, const fileName &subDir, const word &cloudName, IOstream::streamFormat format)
void ensightVolField(const ensightParts &partsList, const IOobject &fieldObject, const fvMesh &mesh, const fileName &dataDir, const fileName &subDir, IOstream::streamFormat format)
Write generalised field components.
void ensightCaseEntry(OFstream &caseFile, const string &ensightType, const word &fieldName, const fileName &dataMask, const fileName &local=fileName::null, const label cloudNo=-1, const label timeSet=1)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:294
Miscellaneous collection of functions and template related to Ensight data.
static const zero Zero
Definition: zero.H:97
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Istream and Ostream manipulators taking arguments.
static const char nl
Definition: Ostream.H:260
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Ostream & flush(Ostream &os)
Flush stream.
Definition: Ostream.H:243
messageStream Info
void ensightLagrangianField(const IOobject &fieldObject, const fileName &dataDir, const fileName &subDir, const word &cloudName, IOstream::streamFormat format)
Write lagrangian parcels.
dimensioned< scalar > mag(const dimensioned< Type > &)
const word cloudName(propsDict.lookup("cloudName"))
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
rDeltaTY field()
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
const dimensionedScalar e
Elementary charge.
Definition: doubleScalar.H:105