FacePostProcessing.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-2020 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 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
78 
79 template<class CloudType>
81 {
82  const fvMesh& mesh = this->owner().mesh();
83  const Time& time = mesh.time();
84  const faceZoneMesh& fzm = mesh.faceZones();
85  scalar timeNew = time.value();
86  scalar timeElapsed = timeNew - timeOld_;
87 
88  totalTime_ += timeElapsed;
89 
90  const scalar alpha = (totalTime_ - timeElapsed)/totalTime_;
91  const scalar beta = timeElapsed/totalTime_;
92 
93  forAll(faceZoneIDs_, zoneI)
94  {
95  massFlowRate_[zoneI] =
96  alpha*massFlowRate_[zoneI] + beta*mass_[zoneI]/timeElapsed;
97  massTotal_[zoneI] += mass_[zoneI];
98  }
99 
100  const label proci = Pstream::myProcNo();
101 
102  Info<< type() << " output:" << nl;
103 
104  List<scalarField> zoneMassTotal(mass_.size());
105  List<scalarField> zoneMassFlowRate(massFlowRate_.size());
106  forAll(faceZoneIDs_, zoneI)
107  {
108  const word& zoneName = fzm[faceZoneIDs_[zoneI]].name();
109 
110  scalarListList allProcMass(Pstream::nProcs());
111  allProcMass[proci] = massTotal_[zoneI];
112  Pstream::gatherList(allProcMass);
113  zoneMassTotal[zoneI] =
114  ListListOps::combine<scalarList>
115  (
116  allProcMass, accessOp<scalarList>()
117  );
118  const scalar sumMassTotal = sum(zoneMassTotal[zoneI]);
119 
120  scalarListList allProcMassFlowRate(Pstream::nProcs());
121  allProcMassFlowRate[proci] = massFlowRate_[zoneI];
122  Pstream::gatherList(allProcMassFlowRate);
123  zoneMassFlowRate[zoneI] =
124  ListListOps::combine<scalarList>
125  (
126  allProcMassFlowRate, accessOp<scalarList>()
127  );
128  const scalar sumMassFlowRate = sum(zoneMassFlowRate[zoneI]);
129 
130  Info<< " " << zoneName
131  << ": total mass = " << sumMassTotal
132  << "; average mass flow rate = " << sumMassFlowRate
133  << nl;
134 
135  if (outputFilePtr_.set(zoneI))
136  {
137  OFstream& os = outputFilePtr_[zoneI];
138  os << time.timeName() << token::TAB << sumMassTotal << token::TAB
139  << sumMassFlowRate<< endl;
140  }
141  }
142 
143  Info<< endl;
144 
145 
146  if (surfaceFormat_ != "none")
147  {
148  forAll(faceZoneIDs_, zoneI)
149  {
150  const faceZone& fZone = fzm[faceZoneIDs_[zoneI]];
151 
152  labelList pointToGlobal;
153  labelList uniqueMeshPointLabels;
154  autoPtr<globalIndex> globalPointsPtr =
155  mesh.globalData().mergePoints
156  (
157  fZone().meshPoints(),
158  fZone().meshPointMap(),
159  pointToGlobal,
160  uniqueMeshPointLabels
161  );
162 
163  pointField uniquePoints(mesh.points(), uniqueMeshPointLabels);
164  List<pointField> allProcPoints(Pstream::nProcs());
165  allProcPoints[proci] = uniquePoints;
166  Pstream::gatherList(allProcPoints);
167 
168  faceList faces(fZone().localFaces());
169  forAll(faces, i)
170  {
171  inplaceRenumber(pointToGlobal, faces[i]);
172  }
173  List<faceList> allProcFaces(Pstream::nProcs());
174  allProcFaces[proci] = faces;
175  Pstream::gatherList(allProcFaces);
176 
177  if (Pstream::master())
178  {
180  (
181  ListListOps::combine<pointField>
182  (
183  allProcPoints, accessOp<pointField>()
184  )
185  );
186 
187  faceList allFaces
188  (
189  ListListOps::combine<faceList>
190  (
191  allProcFaces, accessOp<faceList>()
192  )
193  );
194 
196  (
197  surfaceWriter::New(surfaceFormat_, this->coeffDict())
198  );
199 
200  writer->write
201  (
202  this->writeTimeDir(),
203  fZone.name(),
204  allPoints,
205  allFaces,
206  "massTotal",
207  zoneMassTotal[zoneI],
208  false
209  );
210 
211  writer->write
212  (
213  this->writeTimeDir(),
214  fZone.name(),
215  allPoints,
216  allFaces,
217  "massFlowRate",
218  zoneMassFlowRate[zoneI],
219  false
220  );
221  }
222  }
223  }
224 
225 
226  if (resetOnWrite_)
227  {
228  forAll(faceZoneIDs_, zoneI)
229  {
230  massFlowRate_[zoneI] = 0.0;
231  }
232  timeOld_ = timeNew;
233  totalTime_ = 0.0;
234  }
235 
236  forAll(mass_, zoneI)
237  {
238  mass_[zoneI] = 0.0;
239  }
240 
241  // writeProperties();
242 }
243 
244 
245 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
246 
247 template<class CloudType>
249 (
250  const dictionary& dict,
251  CloudType& owner,
252  const word& modelName
253 )
254 :
255  CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
256  faceZoneIDs_(),
257  surfaceFormat_(this->coeffDict().lookup("surfaceFormat")),
258  resetOnWrite_(this->coeffDict().lookup("resetOnWrite")),
259  totalTime_(0.0),
260  mass_(),
261  massTotal_(),
262  massFlowRate_(),
263  log_(this->coeffDict().lookup("log")),
264  outputFilePtr_(),
265  timeOld_(owner.mesh().time().value())
266 {
267  wordList faceZoneNames(this->coeffDict().lookup("faceZones"));
268  mass_.setSize(faceZoneNames.size());
269  massTotal_.setSize(faceZoneNames.size());
270  massFlowRate_.setSize(faceZoneNames.size());
271 
272  outputFilePtr_.setSize(faceZoneNames.size());
273 
274  DynamicList<label> zoneIDs;
275  const faceZoneMesh& fzm = owner.mesh().faceZones();
276  const surfaceScalarField& magSf = owner.mesh().magSf();
277  const polyBoundaryMesh& pbm = owner.mesh().boundaryMesh();
278  forAll(faceZoneNames, i)
279  {
280  const word& zoneName = faceZoneNames[i];
281  label zoneI = fzm.findZoneID(zoneName);
282  if (zoneI != -1)
283  {
284  zoneIDs.append(zoneI);
285  const faceZone& fz = fzm[zoneI];
286  mass_[i].setSize(fz.size(), 0.0);
287  massTotal_[i].setSize(fz.size(), 0.0);
288  massFlowRate_[i].setSize(fz.size(), 0.0);
289 
290  label nFaces = returnReduce(fz.size(), sumOp<label>());
291  Info<< " " << zoneName << " faces: " << nFaces << nl;
292 
293  scalar totArea = 0.0;
294  forAll(fz, j)
295  {
296  label facei = fz[j];
297  if (facei < owner.mesh().nInternalFaces())
298  {
299  totArea += magSf[fz[j]];
300  }
301  else
302  {
303  label bFacei = facei - owner.mesh().nInternalFaces();
304  label patchi = pbm.patchID()[bFacei];
305  const polyPatch& pp = pbm[patchi];
306 
307  if
308  (
309  !magSf.boundaryField()[patchi].coupled()
310  || refCast<const coupledPolyPatch>(pp).owner()
311  )
312  {
313  label localFacei = pp.whichFace(facei);
314  totArea += magSf.boundaryField()[patchi][localFacei];
315  }
316  }
317  }
318  totArea = returnReduce(totArea, sumOp<scalar>());
319 
320  makeLogFile(zoneName, i, nFaces, totArea);
321  }
322  }
323 
324  faceZoneIDs_.transfer(zoneIDs);
325 
326  // readProperties(); AND initialise mass... fields
327 }
328 
329 
330 template<class CloudType>
332 (
334 )
335 :
337  faceZoneIDs_(pff.faceZoneIDs_),
338  surfaceFormat_(pff.surfaceFormat_),
339  resetOnWrite_(pff.resetOnWrite_),
340  totalTime_(pff.totalTime_),
341  mass_(pff.mass_),
342  massTotal_(pff.massTotal_),
343  massFlowRate_(pff.massFlowRate_),
344  log_(pff.log_),
345  outputFilePtr_(),
346  timeOld_(0.0)
347 {}
348 
349 
350 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
351 
352 template<class CloudType>
354 {}
355 
356 
357 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
358 
359 template<class CloudType>
360 void Foam::FacePostProcessing<CloudType>::postFace(const parcelType& p, bool&)
361 {
362  if
363  (
364  this->owner().solution().output()
365  || this->owner().solution().transient()
366  )
367  {
368  const faceZoneMesh& fzm = this->owner().mesh().faceZones();
369 
370  forAll(faceZoneIDs_, i)
371  {
372  const faceZone& fz = fzm[faceZoneIDs_[i]];
373 
374  label faceId = -1;
375  forAll(fz, j)
376  {
377  if (fz[j] == p.face())
378  {
379  faceId = j;
380  break;
381  }
382  }
383 
384  if (faceId != -1)
385  {
386  mass_[i][faceId] += p.mass()*p.nParticle();
387  }
388  }
389  }
390 }
391 
392 
393 // ************************************************************************* //
dictionary dict
#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
label faceId(-1)
const word & name() const
Return name.
Definition: IOobject.H:303
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:476
void write()
Write post-processing info.
static const char tab
Definition: Ostream.H:259
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Output to file stream.
Definition: OFstream.H:82
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const labelList & patchID() const
Per boundary face label the patch index.
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
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:622
autoPtr< globalIndex > mergePoints(labelList &pointToGlobal, labelList &uniquePoints) const
Helper for merging (collocated!) mesh point data.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
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:1131
virtual void postFace(const parcelType &p, bool &keepParticle)
Post-face hook.
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:341
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
Records particle face quantities on used-specified face zone.
const Type & value() const
Return const reference to value.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1394
const word & name() const
Return name.
Definition: zone.H:147
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
Foam::polyBoundaryMesh.
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
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
messageStream Info
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
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:52
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
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
volScalarField & p
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
Templated cloud function object base class.
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:383