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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "FacePostProcessing.H"
30 #include "Pstream.H"
31 #include "ListListOps.H"
32 #include "globalIndex.H"
33 #include "surfaceWriter.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 template<class CloudType>
39 (
40  const word& zoneName,
41  const label zoneI,
42  const label nFaces,
43  const scalar totArea
44 )
45 {
46  // Create the output file if not already created
47  if (logToFile_)
48  {
49  DebugInfo << "Creating output file." << endl;
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  Log_<< 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  Log_<< " " << 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  Log_<< 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  surfaceFormat_,
199  (
200  this->coeffDict(),
201  surfaceFormat_
202  )
203  );
204 
205  if (debug)
206  {
207  writer->verbose(true);
208  }
209 
210  writer->open
211  (
212  allPoints,
213  allFaces,
214  (this->writeTimeDir() / fZone.name()),
215  false // serial - already merged
216  );
217 
218  writer->nFields(2); // Legacy VTK
219  writer->write("massTotal", zoneMassTotal[zoneI]);
220  writer->write("massFlowRate", zoneMassFlowRate[zoneI]);
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().getWord("surfaceFormat")),
258  resetOnWrite_(this->coeffDict().getBool("resetOnWrite")),
259  logToFile_(this->coeffDict().getBool("log")),
260  totalTime_(0.0),
261  mass_(),
262  massTotal_(),
263  massFlowRate_(),
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 
275  const faceZoneMesh& fzm = owner.mesh().faceZones();
276  const surfaceScalarField& magSf = owner.mesh().magSf();
278 
279  forAll(faceZoneNames, i)
280  {
281  const word& zoneName = faceZoneNames[i];
282  label zoneI = fzm.findZoneID(zoneName);
283 
284  if (zoneI != -1)
285  {
286  zoneIDs.append(zoneI);
287  const faceZone& fz = fzm[zoneI];
288  mass_[i].setSize(fz.size(), 0.0);
289  massTotal_[i].setSize(fz.size(), 0.0);
290  massFlowRate_[i].setSize(fz.size(), 0.0);
291 
292  label nFaces = returnReduce(fz.size(), sumOp<label>());
293  Info<< " " << zoneName << " faces: " << nFaces << nl;
294 
295  scalar totArea = 0.0;
296  for (const label facei : fz)
297  {
298  if (facei < owner.mesh().nInternalFaces())
299  {
300  totArea += magSf[facei];
301  }
302  else
303  {
304  const label patchi = pbm.patchID(facei);
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);
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  logToFile_(pff.logToFile_),
341  totalTime_(pff.totalTime_),
342  mass_(pff.mass_),
343  massTotal_(pff.massTotal_),
344  massFlowRate_(pff.massFlowRate_),
345  outputFilePtr_(),
346  timeOld_(0.0)
347 {}
348 
349 
350 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
351 
352 template<class CloudType>
354 (
355  const parcelType& p,
356  const typename parcelType::trackingData& td
357 )
358 {
359  if
360  (
361  this->owner().solution().output()
362  || this->owner().solution().transient()
363  )
364  {
365  const faceZoneMesh& fzm = this->owner().mesh().faceZones();
366 
367  forAll(faceZoneIDs_, i)
368  {
369  const faceZone& fz = fzm[faceZoneIDs_[i]];
370 
371  label faceId = -1;
372  forAll(fz, j)
373  {
374  if (fz[j] == p.face())
375  {
376  faceId = j;
377  break;
378  }
379  }
380 
381  if (faceId != -1)
382  {
383  mass_[i][faceId] += p.mass()*p.nParticle();
384  }
385  }
386  }
387 
388  return true;
389 }
390 
391 
392 // ************************************************************************* //
const polyBoundaryMesh & pbm
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
const labelIOList & zoneIDs
Definition: correctPhi.H:59
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
type
Types of root.
Definition: Roots.H:52
label faceId(-1)
void write()
Write post-processing info.
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/"finiteArea-edgesCentres"))
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
FacePostProcessing(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
Lookup type of boundary radiation properties.
Definition: lookup.H:57
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const CloudType & owner() const
Return const access to the owner cloud.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void inplaceRenumber(const labelUList &oldToNew, IntListType &lists)
Inplace renumber the values (not the indices) of a list of lists.
Definition: ensightOutput.H:92
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:614
#define Log_
Report write to Foam::Info if the class log switch is true.
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:629
List< scalarList > scalarListList
List of scalarList.
Definition: scalarList.H:35
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
Records particle face quantities on used-specified face zone.
label nInternalFaces() const noexcept
Number of internal faces.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:122
#define DebugInfo
Report an information message using Foam::Info.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:670
virtual bool postFace(const parcelType &p, const typename parcelType::trackingData &td)
Post-face hook.
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:37
messageStream Info
Information stream (stdout output on master, null elsewhere)
static Ostream & output(Ostream &os, const IntRange< T > &range)
Definition: IntRanges.C:44
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:92
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
Templated cloud function object base class.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const dictionary formatOptions(propsDict.subOrEmptyDict("formatOptions", keyType::LITERAL))