FacePostProcessing.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-2016 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 "FacePostProcessing.H"
27 #include "Pstream.H"
28 #include "ListListOps.H"
29 #include "surfaceWriter.H"
30 #include "globalIndex.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 template<class CloudType>
36 (
37  const word& zoneName,
38  const label zoneI,
39  const label nFaces,
40  const scalar totArea
41 )
42 {
43  // Create the output file if not already created
44  if (log_)
45  {
46  if (debug)
47  {
48  Info<< "Creating output file." << endl;
49  }
50 
51  if (Pstream::master())
52  {
53  // Create directory if does not exist
54  mkDir(this->writeTimeDir());
55 
56  // Open new file at start up
57  outputFilePtr_.set
58  (
59  zoneI,
60  new OFstream
61  (
62  this->writeTimeDir()/(type() + '_' + zoneName + ".dat")
63  )
64  );
65 
66  outputFilePtr_[zoneI]
67  << "# Source : " << type() << nl
68  << "# Face zone : " << zoneName << nl
69  << "# Faces : " << nFaces << nl
70  << "# Area : " << totArea << nl
71  << "# Time" << tab << "mass" << tab << "massFlowRate" << endl;
72  }
73  }
74 }
75 
76 
77 template<class CloudType>
79 {
80  const fvMesh& mesh = this->owner().mesh();
81  const Time& time = mesh.time();
82  const faceZoneMesh& fzm = mesh.faceZones();
83  scalar timeNew = time.value();
84  scalar timeElapsed = timeNew - timeOld_;
85 
86  totalTime_ += timeElapsed;
87 
88  const scalar alpha = (totalTime_ - timeElapsed)/totalTime_;
89  const scalar beta = timeElapsed/totalTime_;
90 
91  forAll(faceZoneIDs_, zoneI)
92  {
93  massFlowRate_[zoneI] =
94  alpha*massFlowRate_[zoneI] + beta*mass_[zoneI]/timeElapsed;
95  massTotal_[zoneI] += mass_[zoneI];
96  }
97 
98  const label proci = Pstream::myProcNo();
99 
100  Info<< type() << " output:" << nl;
101 
102  List<scalarField> zoneMassTotal(mass_.size());
103  List<scalarField> zoneMassFlowRate(massFlowRate_.size());
104  forAll(faceZoneIDs_, zoneI)
105  {
106  const word& zoneName = fzm[faceZoneIDs_[zoneI]].name();
107 
108  scalarListList allProcMass(Pstream::nProcs());
109  allProcMass[proci] = massTotal_[zoneI];
110  Pstream::gatherList(allProcMass);
111  zoneMassTotal[zoneI] =
112  ListListOps::combine<scalarList>
113  (
114  allProcMass, accessOp<scalarList>()
115  );
116  const scalar sumMassTotal = sum(zoneMassTotal[zoneI]);
117 
118  scalarListList allProcMassFlowRate(Pstream::nProcs());
119  allProcMassFlowRate[proci] = massFlowRate_[zoneI];
120  Pstream::gatherList(allProcMassFlowRate);
121  zoneMassFlowRate[zoneI] =
122  ListListOps::combine<scalarList>
123  (
124  allProcMassFlowRate, accessOp<scalarList>()
125  );
126  const scalar sumMassFlowRate = sum(zoneMassFlowRate[zoneI]);
127 
128  Info<< " " << zoneName
129  << ": total mass = " << sumMassTotal
130  << "; average mass flow rate = " << sumMassFlowRate
131  << nl;
132 
133  if (outputFilePtr_.set(zoneI))
134  {
135  OFstream& os = outputFilePtr_[zoneI];
136  os << time.timeName() << token::TAB << sumMassTotal << token::TAB
137  << sumMassFlowRate<< endl;
138  }
139  }
140 
141  Info<< endl;
142 
143 
144  if (surfaceFormat_ != "none")
145  {
146  forAll(faceZoneIDs_, zoneI)
147  {
148  const faceZone& fZone = fzm[faceZoneIDs_[zoneI]];
149 
150  labelList pointToGlobal;
151  labelList uniqueMeshPointLabels;
152  autoPtr<globalIndex> globalPointsPtr =
153  mesh.globalData().mergePoints
154  (
155  fZone().meshPoints(),
156  fZone().meshPointMap(),
157  pointToGlobal,
158  uniqueMeshPointLabels
159  );
160 
161  pointField uniquePoints(mesh.points(), uniqueMeshPointLabels);
162  List<pointField> allProcPoints(Pstream::nProcs());
163  allProcPoints[proci] = uniquePoints;
164  Pstream::gatherList(allProcPoints);
165 
166  faceList faces(fZone().localFaces());
167  forAll(faces, i)
168  {
169  inplaceRenumber(pointToGlobal, faces[i]);
170  }
171  List<faceList> allProcFaces(Pstream::nProcs());
172  allProcFaces[proci] = faces;
173  Pstream::gatherList(allProcFaces);
174 
175  if (Pstream::master())
176  {
178  (
179  ListListOps::combine<pointField>
180  (
181  allProcPoints, accessOp<pointField>()
182  )
183  );
184 
185  faceList allFaces
186  (
187  ListListOps::combine<faceList>
188  (
189  allProcFaces, accessOp<faceList>()
190  )
191  );
192 
194  (
196  (
197  surfaceFormat_,
198  this->coeffDict().subOrEmptyDict("formatOptions").
199  subOrEmptyDict(surfaceFormat_)
200  )
201  );
202 
203  writer->write
204  (
205  this->writeTimeDir(),
206  fZone.name(),
207  allPoints,
208  allFaces,
209  "massTotal",
210  zoneMassTotal[zoneI],
211  false
212  );
213 
214  writer->write
215  (
216  this->writeTimeDir(),
217  fZone.name(),
218  allPoints,
219  allFaces,
220  "massFlowRate",
221  zoneMassFlowRate[zoneI],
222  false
223  );
224  }
225  }
226  }
227 
228 
229  if (resetOnWrite_)
230  {
231  forAll(faceZoneIDs_, zoneI)
232  {
233  massFlowRate_[zoneI] = 0.0;
234  }
235  timeOld_ = timeNew;
236  totalTime_ = 0.0;
237  }
238 
239  forAll(mass_, zoneI)
240  {
241  mass_[zoneI] = 0.0;
242  }
243 
244  // writeProperties();
245 }
246 
247 
248 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
249 
250 template<class CloudType>
252 (
253  const dictionary& dict,
254  CloudType& owner,
255  const word& modelName
256 )
257 :
258  CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
259  faceZoneIDs_(),
260  surfaceFormat_(this->coeffDict().lookup("surfaceFormat")),
261  resetOnWrite_(this->coeffDict().lookup("resetOnWrite")),
262  totalTime_(0.0),
263  mass_(),
264  massTotal_(),
265  massFlowRate_(),
266  log_(this->coeffDict().lookup("log")),
267  outputFilePtr_(),
268  timeOld_(owner.mesh().time().value())
269 {
270  wordList faceZoneNames(this->coeffDict().lookup("faceZones"));
271  mass_.setSize(faceZoneNames.size());
272  massTotal_.setSize(faceZoneNames.size());
273  massFlowRate_.setSize(faceZoneNames.size());
274 
275  outputFilePtr_.setSize(faceZoneNames.size());
276 
277  DynamicList<label> zoneIDs;
278  const faceZoneMesh& fzm = owner.mesh().faceZones();
279  const surfaceScalarField& magSf = owner.mesh().magSf();
280  const polyBoundaryMesh& pbm = owner.mesh().boundaryMesh();
281  forAll(faceZoneNames, i)
282  {
283  const word& zoneName = faceZoneNames[i];
284  label zoneI = fzm.findZoneID(zoneName);
285  if (zoneI != -1)
286  {
287  zoneIDs.append(zoneI);
288  const faceZone& fz = fzm[zoneI];
289  mass_[i].setSize(fz.size(), 0.0);
290  massTotal_[i].setSize(fz.size(), 0.0);
291  massFlowRate_[i].setSize(fz.size(), 0.0);
292 
293  label nFaces = returnReduce(fz.size(), sumOp<label>());
294  Info<< " " << zoneName << " faces: " << nFaces << nl;
295 
296  scalar totArea = 0.0;
297  forAll(fz, j)
298  {
299  label facei = fz[j];
300  if (facei < owner.mesh().nInternalFaces())
301  {
302  totArea += magSf[fz[j]];
303  }
304  else
305  {
306  label bFacei = facei - owner.mesh().nInternalFaces();
307  label patchi = pbm.patchID()[bFacei];
308  const polyPatch& pp = pbm[patchi];
309 
310  if
311  (
312  !magSf.boundaryField()[patchi].coupled()
313  || refCast<const coupledPolyPatch>(pp).owner()
314  )
315  {
316  label localFacei = pp.whichFace(facei);
317  totArea += magSf.boundaryField()[patchi][localFacei];
318  }
319  }
320  }
321  totArea = returnReduce(totArea, sumOp<scalar>());
322 
323  makeLogFile(zoneName, i, nFaces, totArea);
324  }
325  }
326 
327  faceZoneIDs_.transfer(zoneIDs);
328 
329  // readProperties(); AND initialise mass... fields
330 }
331 
332 
333 template<class CloudType>
335 (
337 )
338 :
340  faceZoneIDs_(pff.faceZoneIDs_),
341  surfaceFormat_(pff.surfaceFormat_),
342  resetOnWrite_(pff.resetOnWrite_),
343  totalTime_(pff.totalTime_),
344  mass_(pff.mass_),
345  massTotal_(pff.massTotal_),
346  massFlowRate_(pff.massFlowRate_),
347  log_(pff.log_),
348  outputFilePtr_(),
349  timeOld_(0.0)
350 {}
351 
352 
353 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
354 
355 template<class CloudType>
357 {}
358 
359 
360 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
361 
362 template<class CloudType>
364 (
365  const parcelType& p,
366  const label facei,
367  bool&
368 )
369 {
370  if
371  (
372  this->owner().solution().output()
373  || this->owner().solution().transient()
374  )
375  {
376  const faceZoneMesh& fzm = this->owner().mesh().faceZones();
377 
378  forAll(faceZoneIDs_, i)
379  {
380  const faceZone& fz = fzm[faceZoneIDs_[i]];
381 
382  label faceId = -1;
383  forAll(fz, j)
384  {
385  if (fz[j] == facei)
386  {
387  faceId = j;
388  break;
389  }
390  }
391 
392  if (faceId != -1)
393  {
394  mass_[i][faceId] += p.mass()*p.nParticle();
395  }
396  }
397  }
398 }
399 
400 
401 // ************************************************************************* //
const labelList & patchID() const
Per boundary face label the patch index.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
label faceId(-1)
autoPtr< globalIndex > mergePoints(labelList &pointToGlobal, labelList &uniquePoints) const
Helper for merging (collocated!) mesh point data.
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:377
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
void write()
Write post-processing info.
static const char tab
Definition: Ostream.H:261
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Output to file stream.
Definition: OFstream.H:81
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
FacePostProcessing(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Base class for graphics format writing. Entry points are.
Definition: writer.H:78
const Type & value() const
Return const reference to value.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
virtual ~FacePostProcessing()
Destructor.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
virtual void postFace(const parcelType &p, const label facei, bool &keepParticle)
Post-face hook.
const word & name() const
Return name.
Definition: zone.H:150
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:341
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A class for handling words, derived from string.
Definition: word.H:59
Records particle face quantities on used-specified face zone.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
Foam::polyBoundaryMesh.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1140
static const char nl
Definition: Ostream.H:262
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:295
label patchi
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
messageStream Info
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
Templated cloud function object base class.
const word & name() const
Return name.
Definition: IOobject.H:260
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243