fluxSummary.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) 2015 OpenFOAM Foundation
9  Copyright (C) 2015-2022 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 "fluxSummary.H"
30 #include "surfaceFields.H"
31 #include "polySurfaceFields.H"
32 #include "dictionary.H"
33 #include "Time.H"
34 #include "syncTools.H"
35 #include "meshTools.H"
36 #include "PatchEdgeFaceWave.H"
37 #include "edgeTopoDistanceData.H"
38 #include "globalIndex.H"
39 #include "OBJstream.H"
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 namespace functionObjects
47 {
48  defineTypeNameAndDebug(fluxSummary, 0);
49  addToRunTimeSelectionTable(functionObject, fluxSummary, dictionary);
50 }
51 }
52 
53 const Foam::Enum
54 <
56 >
58 ({
59  { modeType::mdFaceZone , "faceZone" },
60  { modeType::mdFaceZoneAndDirection, "faceZoneAndDirection" },
61  { modeType::mdCellZoneAndDirection, "cellZoneAndDirection" },
62  { modeType::mdSurface, "functionObjectSurface" },
63  { modeType::mdSurface, "surface" },
64  { modeType::mdSurfaceAndDirection, "surfaceAndDirection" },
65 });
66 
67 
68 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
69 
71 {
72  return (mdSurface == mode_ || mdSurfaceAndDirection == mode_);
73 }
74 
75 
77 (
78  const dimensionSet& fieldDims,
79  const word& fieldName
80 ) const
81 {
82  // Surfaces are multiplied by their area, so account for that
83  // in the dimension checking
84  const dimensionSet dims =
85  (fieldDims * (isSurfaceMode() ? dimTime*dimArea : dimTime));
86 
87  if (dims == dimVolume)
88  {
89  return "volumetric";
90  }
91  else if (dims == dimMass)
92  {
93  return "mass";
94  }
95 
97  << "Unsupported flux field " << fieldName << " with dimensions "
98  << fieldDims
99  << ". Expected either mass flow or volumetric flow rate."
101 
102  return word::null;
103 }
104 
105 
107 (
108  const word& surfName,
111  DynamicList<boolList>& faceFlip
112 ) const
113 {
114  const polySurface* surfptr =
115  storedObjects().cfindObject<polySurface>(surfName);
116 
117  if (!surfptr)
118  {
120  << "Unable to find surface " << surfName
121  << ". Valid surfaces: "
122  << storedObjects().sortedNames<polySurface>() << nl
123  << exit(FatalError);
124  }
125 
126  names.append(surfName);
127  directions.append(Zero); // Dummy value
128  faceFlip.append(boolList()); // No flip-map
129 }
130 
131 
133 (
134  const word& surfName,
135  const vector& dir,
138  DynamicList<boolList>& faceFlip
139 ) const
140 {
141  const polySurface* surfptr =
142  storedObjects().cfindObject<polySurface>(surfName);
143 
144  if (!surfptr)
145  {
147  << "Unable to find surface " << surfName
148  << ". Valid surfaces: "
149  << storedObjects().sortedNames<polySurface>() << nl
150  << exit(FatalError);
151  }
152 
153  const auto& s = *surfptr;
154  const vector refDir = dir/(mag(dir) + ROOTVSMALL);
155 
156  names.append(surfName);
157  directions.append(refDir);
158  faceFlip.append(boolList()); // No flip-map
159 
160  boolList& flips = faceFlip[faceFlip.size()-1];
161  flips.setSize(s.size(), false);
162 
163  forAll(s, i)
164  {
165  // Orientation set by comparison with reference direction
166  const vector& n = s.faceNormals()[i];
167 
168  if ((n & refDir) > tolerance_)
169  {
170  flips[i] = false;
171  }
172  else
173  {
174  flips[i] = true;
175  }
176  }
177 }
178 
179 
181 (
182  const word& faceZoneName,
183  DynamicList<word>& names,
184  DynamicList<vector>& directions,
185  DynamicList<labelList>& faceID,
186  DynamicList<labelList>& facePatchID,
187  DynamicList<boolList>& faceFlip
188 ) const
189 {
190  label zonei = mesh_.faceZones().findZoneID(faceZoneName);
191  if (zonei == -1)
192  {
194  << "Unable to find faceZone " << faceZoneName
195  << ". Valid zones: "
196  << mesh_.faceZones().sortedNames() << nl
197  << exit(FatalError);
198  }
199  const faceZone& fZone = mesh_.faceZones()[zonei];
200 
201  names.append(faceZoneName);
202  directions.append(Zero); // dummy value
203 
204  labelList faceIds(fZone.size());
205  labelList facePatchIds(fZone.size());
206  boolList faceFlips(fZone.size());
207 
208  // Total number of faces selected
209  label numFaces = 0;
210 
211  forAll(fZone, i)
212  {
213  const label meshFacei = fZone[i];
214  const bool isFlip = fZone.flipMap()[i];
215 
216  // Internal faces
217  label faceId = meshFacei;
218  label facePatchId = -1;
219 
220  // Boundary faces
221  if (!mesh_.isInternalFace(meshFacei))
222  {
223  facePatchId = mesh_.boundaryMesh().whichPatch(meshFacei);
224  const polyPatch& pp = mesh_.boundaryMesh()[facePatchId];
225 
226  if (isA<emptyPolyPatch>(pp))
227  {
228  continue; // Ignore empty patch
229  }
230 
231  const auto* cpp = isA<coupledPolyPatch>(pp);
232 
233  if (cpp && !cpp->owner())
234  {
235  continue; // Ignore neighbour side
236  }
237 
238  faceId = pp.whichFace(meshFacei);
239  }
240 
241  if (faceId >= 0)
242  {
243  // Orientation set by faceZone flip map
244  faceIds[numFaces] = faceId;
245  facePatchIds[numFaces] = facePatchId;
246  faceFlips[numFaces] = isFlip;
247 
248  ++numFaces;
249  }
250  }
251 
252  // Shrink to size used
253  faceIds.resize(numFaces);
254  facePatchIds.resize(numFaces);
255  faceFlips.resize(numFaces);
256 
257  faceID.append(std::move(faceIds));
258  facePatchID.append(std::move(facePatchIds));
259  faceFlip.append(std::move(faceFlips));
260 }
261 
262 
264 (
265  const word& faceZoneName,
266  const vector& dir,
269  DynamicList<labelList>& faceID,
270  DynamicList<labelList>& facePatchID,
271  DynamicList<boolList>& faceFlip
272 ) const
273 {
274  const vector refDir = dir/(mag(dir) + ROOTVSMALL);
275 
276  label zonei = mesh_.faceZones().findZoneID(faceZoneName);
277  if (zonei == -1)
278  {
280  << "Unable to find faceZone " << faceZoneName
281  << ". Valid zones: "
282  << mesh_.faceZones().sortedNames() << nl
283  << exit(FatalError);
284  }
285  const faceZone& fZone = mesh_.faceZones()[zonei];
286 
287  names.append(faceZoneName);
288  directions.append(refDir);
289 
290  labelList faceIds(fZone.size());
291  labelList facePatchIds(fZone.size());
292  boolList faceFlips(fZone.size());
293 
294  const surfaceVectorField& Sf = mesh_.Sf();
295  const surfaceScalarField& magSf = mesh_.magSf();
296 
297  vector n(Zero);
298 
299  // Total number of faces selected
300  label numFaces = 0;
301 
302  forAll(fZone, i)
303  {
304  const label meshFacei = fZone[i];
305 
306  // Internal faces
307  label faceId = meshFacei;
308  label facePatchId = -1;
309 
310  // Boundary faces
311  if (!mesh_.isInternalFace(meshFacei))
312  {
313  facePatchId = mesh_.boundaryMesh().whichPatch(meshFacei);
314  const polyPatch& pp = mesh_.boundaryMesh()[facePatchId];
315 
316  if (isA<emptyPolyPatch>(pp))
317  {
318  continue; // Ignore empty patch
319  }
320 
321  const auto* cpp = isA<coupledPolyPatch>(pp);
322 
323  if (cpp && !cpp->owner())
324  {
325  continue; // Ignore neighbour side
326  }
327 
328  faceId = pp.whichFace(meshFacei);
329  }
330 
331  if (faceId >= 0)
332  {
333  // Orientation set by comparison with reference direction
334  if (facePatchId != -1)
335  {
336  n = Sf.boundaryField()[facePatchId][faceId]
337  /(magSf.boundaryField()[facePatchId][faceId] + ROOTVSMALL);
338  }
339  else
340  {
341  n = Sf[faceId]/(magSf[faceId] + ROOTVSMALL);
342  }
343 
344  if ((n & refDir) > tolerance_)
345  {
346  faceFlips[numFaces] = false;
347  }
348  else
349  {
350  faceFlips[numFaces] = true;
351  }
352 
353  faceIds[numFaces] = faceId;
354  facePatchIds[numFaces] = facePatchId;
355 
356  ++numFaces;
357  }
358  }
359 
360  // Shrink to size used
361  faceIds.resize(numFaces);
362  facePatchIds.resize(numFaces);
363  faceFlips.resize(numFaces);
364 
365  faceID.append(std::move(faceIds));
366  facePatchID.append(std::move(facePatchIds));
367  faceFlip.append(std::move(faceFlips));
368 }
369 
370 
372 (
373  const word& cellZoneName,
374  const vector& dir,
377  DynamicList<labelList>& faceID,
378  DynamicList<labelList>& facePatchID,
379  DynamicList<boolList>& faceFlip
380 ) const
381 {
382  const vector refDir = dir/(mag(dir) + ROOTVSMALL);
383 
384  const label cellZonei = mesh_.cellZones().findZoneID(cellZoneName);
385  if (cellZonei == -1)
386  {
388  << "Unable to find cellZone " << cellZoneName
389  << ". Valid zones: "
390  << mesh_.cellZones().sortedNames() << nl
391  << exit(FatalError);
392  }
393 
394  const label nInternalFaces = mesh_.nInternalFaces();
395  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
396 
397  labelList cellAddr(mesh_.nCells(), -1);
398  const labelList& cellIDs = mesh_.cellZones()[cellZonei];
399  labelUIndList(cellAddr, cellIDs) = identity(cellIDs.size());
400  labelList nbrFaceCellAddr(mesh_.nBoundaryFaces(), -1);
401 
402  forAll(pbm, patchi)
403  {
404  const polyPatch& pp = pbm[patchi];
405 
406  if (pp.coupled())
407  {
408  forAll(pp, i)
409  {
410  label facei = pp.start() + i;
411  label nbrFacei = facei - nInternalFaces;
412  label own = mesh_.faceOwner()[facei];
413  nbrFaceCellAddr[nbrFacei] = cellAddr[own];
414  }
415  }
416  }
417 
418  // Correct boundary values for parallel running
419  syncTools::swapBoundaryFaceList(mesh_, nbrFaceCellAddr);
420 
421  // Collect faces
422  DynamicList<label> faceIDs(floor(0.1*mesh_.nFaces()));
423  DynamicList<label> facePatchIDs(faceIDs.size());
424  DynamicList<label> faceLocalPatchIDs(faceIDs.size());
425  DynamicList<bool> flips(faceIDs.size());
426 
427  // Internal faces
428  for (label facei = 0; facei < nInternalFaces; facei++)
429  {
430  const label own = cellAddr[mesh_.faceOwner()[facei]];
431  const label nbr = cellAddr[mesh_.faceNeighbour()[facei]];
432 
433  if (((own != -1) && (nbr == -1)) || ((own == -1) && (nbr != -1)))
434  {
435  vector n = mesh_.faces()[facei].unitNormal(mesh_.points());
436 
437  if ((n & refDir) > tolerance_)
438  {
439  faceIDs.append(facei);
440  faceLocalPatchIDs.append(facei);
441  facePatchIDs.append(-1);
442  flips.append(false);
443  }
444  else if ((n & -refDir) > tolerance_)
445  {
446  faceIDs.append(facei);
447  faceLocalPatchIDs.append(facei);
448  facePatchIDs.append(-1);
449  flips.append(true);
450  }
451  }
452  }
453 
454  // Loop over boundary faces
455  forAll(pbm, patchi)
456  {
457  const polyPatch& pp = pbm[patchi];
458 
459  forAll(pp, localFacei)
460  {
461  const label facei = pp.start() + localFacei;
462  const label own = cellAddr[mesh_.faceOwner()[facei]];
463  const label nbr = nbrFaceCellAddr[facei - nInternalFaces];
464 
465  if ((own != -1) && (nbr == -1))
466  {
467  vector n = mesh_.faces()[facei].unitNormal(mesh_.points());
468 
469  if ((n & refDir) > tolerance_)
470  {
471  faceIDs.append(facei);
472  faceLocalPatchIDs.append(localFacei);
473  facePatchIDs.append(patchi);
474  flips.append(false);
475  }
476  else if ((n & -refDir) > tolerance_)
477  {
478  faceIDs.append(facei);
479  faceLocalPatchIDs.append(localFacei);
480  facePatchIDs.append(patchi);
481  flips.append(true);
482  }
483  }
484  }
485  }
486 
487  // Convert into primitivePatch for convenience
489  (
490  IndirectList<face>(mesh_.faces(), faceIDs),
491  mesh_.points()
492  );
493 
494  if (debug)
495  {
496  OBJstream os(mesh_.time().path()/"patch.obj");
497  faceList faces(patch);
498  os.write(faces, mesh_.points(), false);
499  }
500 
501 
502  // Data on all edges and faces
503  List<edgeTopoDistanceData<label>> allEdgeInfo(patch.nEdges());
504  List<edgeTopoDistanceData<label>> allFaceInfo(patch.size());
505 
506  bool search = true;
507 
508  DebugInfo
509  << "initialiseCellZoneAndDirection: "
510  << "Starting walk to split patch into faceZones"
511  << endl;
512 
513  const globalIndex globalFaces(patch.size());
514 
515  label oldFaceID = 0;
516  label regioni = 0;
517 
518  // Dummy tracking data
519  bool dummyData{false};
520 
521  while (search)
522  {
523  DynamicList<label> changedEdges;
524  DynamicList<edgeTopoDistanceData<label>> changedInfo;
525 
526  label seedFacei = labelMax;
527  for (; oldFaceID < patch.size(); oldFaceID++)
528  {
529  if (!allFaceInfo[oldFaceID].valid<bool>(dummyData))
530  {
531  seedFacei = globalFaces.toGlobal(oldFaceID);
532  break;
533  }
534  }
535  reduce(seedFacei, minOp<label>());
536 
537  if (seedFacei == labelMax)
538  {
539  break;
540  }
541 
542  if (globalFaces.isLocal(seedFacei))
543  {
544  const label localFacei = globalFaces.toLocal(seedFacei);
545  const labelList& fEdges = patch.faceEdges()[localFacei];
546 
547  for (const label edgei : fEdges)
548  {
549  if (allEdgeInfo[edgei].valid<bool>(dummyData))
550  {
552  << "Problem in edge face wave: attempted to assign a "
553  << "value to an edge that has already been visited. "
554  << "Edge info: " << allEdgeInfo[edgei]
555  << endl;
556  }
557 
558  changedEdges.append(edgei);
559  changedInfo.append
560  (
561  edgeTopoDistanceData<label>
562  (
563  0, // distance
564  regioni
565  )
566  );
567  }
568  }
569 
570 
571  PatchEdgeFaceWave
572  <
574  edgeTopoDistanceData<label>
575  > calc
576  (
577  mesh_,
578  patch,
579  changedEdges,
580  changedInfo,
581  allEdgeInfo,
582  allFaceInfo,
583  returnReduce(patch.nEdges(), sumOp<label>())
584  );
585 
586  if (debug)
587  {
588  label nCells = 0;
589  forAll(allFaceInfo, facei)
590  {
591  if (allFaceInfo[facei].data() == regioni)
592  {
593  nCells++;
594  }
595  }
596 
597  Info<< "*** region:" << regioni
598  << " found:" << returnReduce(nCells, sumOp<label>())
599  << " faces" << endl;
600  }
601 
602  regioni++;
603  }
604 
605  // Collect the data per region
606  const label nRegion = regioni;
607 
608  #ifdef FULLDEBUG
609  // May wish to have enabled always
610  if (nRegion == 0)
611  {
613  << "Region split failed" << nl
614  << exit(FatalError);
615  }
616  #endif
617 
618  List<DynamicList<label>> regionFaceIDs(nRegion);
619  List<DynamicList<label>> regionFacePatchIDs(nRegion);
620  List<DynamicList<bool>> regionFaceFlips(nRegion);
621 
622  forAll(allFaceInfo, facei)
623  {
624  regioni = allFaceInfo[facei].data();
625 
626  regionFaceIDs[regioni].append(faceLocalPatchIDs[facei]);
627  regionFacePatchIDs[regioni].append(facePatchIDs[facei]);
628  regionFaceFlips[regioni].append(flips[facei]);
629  }
630 
631  // Transfer to persistent storage
632  forAll(regionFaceIDs, regioni)
633  {
634  const word zoneName = cellZoneName + ":faceZone" + Foam::name(regioni);
635  names.append(zoneName);
636  directions.append(refDir);
637  faceID.append(regionFaceIDs[regioni]);
638  facePatchID.append(regionFacePatchIDs[regioni]);
639  faceFlip.append(regionFaceFlips[regioni]);
640 
641  // Write OBj of faces to file
642  if (debug)
643  {
644  OBJstream os(mesh_.time().path()/zoneName + ".obj");
645  faceList faces(mesh_.faces(), regionFaceIDs[regioni]);
646  os.write(faces, mesh_.points(), false);
647  }
648  }
649 
650  if (log)
651  {
652  Info<< type() << " " << name() << " output:" << nl
653  << " Created " << faceID.size()
654  << " separate face zones from cell zone " << cellZoneName << nl;
655 
656  forAll(names, i)
657  {
658  label nFaces = returnReduce(faceID[i].size(), sumOp<label>());
659  Info<< " " << names[i] << ": "
660  << nFaces << " faces" << nl;
661  }
663  Info<< endl;
664  }
665 }
666 
667 
669 (
670  const label idx
671 ) const
672 {
673  scalar sumMagSf = 0;
674 
675  if (isSurfaceMode())
676  {
677  const polySurface& s =
678  storedObjects().lookupObject<polySurface>(zoneNames_[idx]);
679 
680  sumMagSf = sum(s.magSf());
681  }
682  else
683  {
684  const surfaceScalarField& magSf = mesh_.magSf();
685 
686  const labelList& faceIDs = faceID_[idx];
687  const labelList& facePatchIDs = facePatchID_[idx];
688 
689  forAll(faceIDs, i)
690  {
691  label facei = faceIDs[i];
692 
693  if (facePatchIDs[i] == -1)
694  {
695  sumMagSf += magSf[facei];
696  }
697  else
698  {
699  label patchi = facePatchIDs[i];
700  sumMagSf += magSf.boundaryField()[patchi][facei];
701  }
702  }
703  }
704 
705  return returnReduce(sumMagSf, sumOp<scalar>());
706 }
707 
708 
710 {
711  for (const word& surfName : zoneNames_)
712  {
713  const polySurface& s =
714  storedObjects().lookupObject<polySurface>(surfName);
715 
716  const auto& phi = s.lookupObject<polySurfaceVectorField>(phiName_);
717 
718  Log << type() << ' ' << name() << ' '
719  << checkFlowType(phi.dimensions(), phi.name()) << " write:" << nl;
720  }
721 
722 
723  forAll(zoneNames_, surfi)
724  {
725  const polySurface& s =
726  storedObjects().lookupObject<polySurface>(zoneNames_[surfi]);
727 
728  const auto& phi = s.lookupObject<polySurfaceVectorField>(phiName_);
729 
730  checkFlowType(phi.dimensions(), phi.name());
731 
732  const boolList& flips = faceFlip_[surfi];
733 
734  scalar phiPos(0);
735  scalar phiNeg(0);
736 
737  tmp<scalarField> tphis = phi & s.Sf();
738  const scalarField& phis = tphis();
739 
740  forAll(s, i)
741  {
742  scalar phif = phis[i];
743  if (flips[i])
744  {
745  phif *= -1;
746  }
747 
748  if (phif > 0)
749  {
750  phiPos += phif;
751  }
752  else
753  {
754  phiNeg += phif;
755  }
756  }
757 
758  reduce(phiPos, sumOp<scalar>());
759  reduce(phiNeg, sumOp<scalar>());
760 
761  phiPos *= scaleFactor_;
762  phiNeg *= scaleFactor_;
763 
764  scalar netFlux = phiPos + phiNeg;
765  scalar absoluteFlux = phiPos - phiNeg;
766 
767  Log << " surface " << zoneNames_[surfi] << ':' << nl
768  << " positive : " << phiPos << nl
769  << " negative : " << phiNeg << nl
770  << " net : " << netFlux << nl
771  << " absolute : " << absoluteFlux
772  << nl << endl;
773 
774  if (writeToFile())
775  {
776  filePtrs_[surfi]
777  << time_.value() << token::TAB
778  << phiPos << token::TAB
779  << phiNeg << token::TAB
780  << netFlux << token::TAB
781  << absoluteFlux
782  << endl;
783  }
784  }
786  Log << endl;
787 
788  return true;
789 }
790 
791 
793 {
794  if (!needsUpdate_)
795  {
796  return false;
797  }
798 
799  // Initialise with capacity == number of input names
800  DynamicList<word> faceZoneName(zoneNames_.size());
801  DynamicList<vector> refDir(faceZoneName.capacity());
802  DynamicList<labelList> faceID(faceZoneName.capacity());
803  DynamicList<labelList> facePatchID(faceZoneName.capacity());
804  DynamicList<boolList> faceFlips(faceZoneName.capacity());
805 
806  switch (mode_)
807  {
808  case mdFaceZone:
809  {
810  forAll(zoneNames_, zonei)
811  {
812  initialiseFaceZone
813  (
814  zoneNames_[zonei],
815  faceZoneName,
816  refDir, // fill with dummy value
817  faceID,
818  facePatchID,
819  faceFlips
820  );
821  }
822  break;
823  }
824  case mdFaceZoneAndDirection:
825  {
826  forAll(zoneNames_, zonei)
827  {
828  initialiseFaceZoneAndDirection
829  (
830  zoneNames_[zonei],
831  zoneDirections_[zonei],
832  faceZoneName,
833  refDir,
834  faceID,
835  facePatchID,
836  faceFlips
837  );
838  }
839  break;
840  }
841  case mdCellZoneAndDirection:
842  {
843  forAll(zoneNames_, zonei)
844  {
845  initialiseCellZoneAndDirection
846  (
847  zoneNames_[zonei],
848  zoneDirections_[zonei],
849  faceZoneName,
850  refDir,
851  faceID,
852  facePatchID,
853  faceFlips
854  );
855  }
856  break;
857  }
858  case mdSurface:
859  {
860  forAll(zoneNames_, zonei)
861  {
862  initialiseSurface
863  (
864  zoneNames_[zonei],
865  faceZoneName,
866  refDir,
867  faceFlips
868  );
869  }
870  break;
871  }
872  case mdSurfaceAndDirection:
873  {
874  forAll(zoneNames_, zonei)
875  {
876  initialiseSurfaceAndDirection
877  (
878  zoneNames_[zonei],
879  zoneDirections_[zonei],
880  faceZoneName,
881  refDir,
882  faceFlips
883  );
884  }
885  break;
886  }
887 
888  // Compiler warning if we forgot an enumeration
889  }
890 
891  zoneNames_.transfer(faceZoneName);
892  faceID_.transfer(faceID);
893  facePatchID_.transfer(facePatchID);
894  faceFlip_.transfer(faceFlips);
895 
896  Info<< type() << " " << name() << " output:" << nl;
897 
898  // Calculate and report areas
899  List<scalar> areas(zoneNames_.size());
900  forAll(zoneNames_, zonei)
901  {
902  const word& zoneName = zoneNames_[zonei];
903  areas[zonei] = totalArea(zonei);
904 
905  if (isSurfaceMode())
906  {
907  Info<< " Surface: " << zoneName
908  << ", area: " << areas[zonei] << nl;
909  }
910  else
911  {
912  Info<< " Zone: " << zoneName
913  << ", area: " << areas[zonei] << nl;
914  }
915  }
916  Info<< endl;
917 
918  if (writeToFile())
919  {
920  filePtrs_.resize(zoneNames_.size());
921 
922  forAll(filePtrs_, zonei)
923  {
924  const word& zoneName = zoneNames_[zonei];
925  filePtrs_.set(zonei, newFileAtStartTime(zoneName));
926  writeFileHeader
927  (
928  zoneName,
929  areas[zonei],
930  refDir[zonei],
931  filePtrs_[zonei]
932  );
933  }
934  }
935 
936  Info<< endl;
937 
938  needsUpdate_ = false;
939 
940  return true;
941 }
942 
943 
944 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
945 
947 (
948  const word& name,
949  const Time& runTime,
950  const dictionary& dict
951 )
952 :
954  writeFile(obr_, name),
955  needsUpdate_(true),
956  mode_(mdFaceZone),
957  scaleFactor_(1),
958  phiName_("phi"),
959  zoneNames_(),
960  zoneDirections_(),
961  faceID_(),
962  facePatchID_(),
963  faceFlip_(),
964  filePtrs_(),
965  tolerance_(0.8)
966 {
967  read(dict);
968 }
969 
970 
971 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
972 
974 {
977 
978  needsUpdate_ = true;
979  mode_ = modeTypeNames_.get("mode", dict);
980  phiName_ = dict.getOrDefault<word>("phi", "phi");
981  scaleFactor_ = dict.getOrDefault<scalar>("scaleFactor", 1);
982  tolerance_ = dict.getOrDefault<scalar>("tolerance", 0.8);
983 
984  zoneNames_.clear();
985  zoneDirections_.clear();
986 
987  List<Tuple2<word, vector>> nameAndDirection;
988 
989  switch (mode_)
990  {
991  case mdFaceZone:
992  {
993  dict.readEntry("faceZones", zoneNames_);
994  break;
995  }
996  case mdFaceZoneAndDirection:
997  {
998  dict.readEntry("faceZoneAndDirection", nameAndDirection);
999  break;
1000  }
1001  case mdCellZoneAndDirection:
1002  {
1003  dict.readEntry("cellZoneAndDirection", nameAndDirection);
1004  break;
1005  }
1006  case mdSurface:
1007  {
1008  dict.readEntry("surfaces", zoneNames_);
1009  break;
1010  }
1011  case mdSurfaceAndDirection:
1012  {
1013  dict.readEntry("surfaceAndDirection", nameAndDirection);
1014  break;
1015  }
1016  default:
1017  {
1019  << "unhandled enumeration " << modeTypeNames_[mode_]
1020  << abort(FatalIOError);
1021  }
1022  }
1023 
1024 
1025  // Split name/vector into separate lists
1026  if (nameAndDirection.size())
1027  {
1028  zoneNames_.resize(nameAndDirection.size());
1029  zoneDirections_.resize(nameAndDirection.size());
1030 
1031  label zonei = 0;
1032 
1033  for (const Tuple2<word, vector>& nameDirn : nameAndDirection)
1034  {
1035  zoneNames_[zonei] = nameDirn.first();
1036  zoneDirections_[zonei] = nameDirn.second();
1037  ++zonei;
1038  }
1039 
1040  nameAndDirection.clear();
1041  }
1042 
1043 
1044  Info<< type() << ' ' << name() << " ("
1045  << modeTypeNames_[mode_] << ") with selection:\n "
1046  << flatOutput(zoneNames_) << endl;
1047 
1048  return !zoneNames_.empty();
1049 }
1050 
1051 
1053 (
1054  const word& zoneName,
1055  const scalar area,
1056  const vector& refDir,
1057  Ostream& os
1058 ) const
1059 {
1060  writeHeader(os, "Flux summary");
1061  if (isSurfaceMode())
1062  {
1063  writeHeaderValue(os, "Surface", zoneName);
1064  }
1065  else
1066  {
1067  writeHeaderValue(os, "Face zone", zoneName);
1068  }
1069  writeHeaderValue(os, "Total area", area);
1070 
1071  switch (mode_)
1072  {
1073  case mdFaceZoneAndDirection:
1074  case mdCellZoneAndDirection:
1075  case mdSurfaceAndDirection:
1076  {
1077  writeHeaderValue(os, "Reference direction", refDir);
1078  break;
1079  }
1080  default:
1081  {}
1082  }
1083 
1084  writeHeaderValue(os, "Scale factor", scaleFactor_);
1085 
1086  writeCommented(os, "Time");
1087  os << tab << "positive"
1088  << tab << "negative"
1089  << tab << "net"
1090  << tab << "absolute"
1091  << endl;
1092 }
1093 
1096 {
1097  return true;
1098 }
1099 
1100 
1102 {
1103  update();
1104 
1105  if (isSurfaceMode())
1106  {
1107  return surfaceModeWrite();
1108  }
1109 
1110  const surfaceScalarField& phi = lookupObject<surfaceScalarField>(phiName_);
1111 
1112  Log << type() << ' ' << name() << ' '
1113  << checkFlowType(phi.dimensions(), phi.name()) << " write:" << nl;
1114 
1115  forAll(zoneNames_, zonei)
1116  {
1117  const labelList& faceID = faceID_[zonei];
1118  const labelList& facePatchID = facePatchID_[zonei];
1119  const boolList& faceFlips = faceFlip_[zonei];
1120 
1121  scalar phiPos(0);
1122  scalar phiNeg(0);
1123  scalar phif(0);
1124 
1125  forAll(faceID, i)
1126  {
1127  label facei = faceID[i];
1128  label patchi = facePatchID[i];
1129 
1130  if (patchi != -1)
1131  {
1132  phif = phi.boundaryField()[patchi][facei];
1133  }
1134  else
1135  {
1136  phif = phi[facei];
1137  }
1138 
1139  if (faceFlips[i])
1140  {
1141  phif *= -1;
1142  }
1143 
1144  if (phif > 0)
1145  {
1146  phiPos += phif;
1147  }
1148  else
1149  {
1150  phiNeg += phif;
1151  }
1152  }
1153 
1154  reduce(phiPos, sumOp<scalar>());
1155  reduce(phiNeg, sumOp<scalar>());
1156 
1157  phiPos *= scaleFactor_;
1158  phiNeg *= scaleFactor_;
1159 
1160  scalar netFlux = phiPos + phiNeg;
1161  scalar absoluteFlux = phiPos - phiNeg;
1162 
1163  Log << " faceZone " << zoneNames_[zonei] << ':' << nl
1164  << " positive : " << phiPos << nl
1165  << " negative : " << phiNeg << nl
1166  << " net : " << netFlux << nl
1167  << " absolute : " << absoluteFlux
1168  << nl << endl;
1169 
1170  if (writeToFile())
1171  {
1172  filePtrs_[zonei]
1173  << time_.value() << token::TAB
1174  << phiPos << token::TAB
1175  << phiNeg << token::TAB
1176  << netFlux << token::TAB
1177  << absoluteFlux
1178  << endl;
1179  }
1180  }
1181 
1182  Log << endl;
1183 
1184  return true;
1185 }
1186 
1187 
1188 // ************************************************************************* //
A surface mesh consisting of general polygon faces and capable of holding fields. ...
Definition: polySurface.H:62
Foam::surfaceFields.
const polyBoundaryMesh & pbm
modeType mode_
Mode for face determination/to generate faces to test.
Definition: fluxSummary.H:227
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
defineTypeNameAndDebug(ObukhovLength, 0)
virtual bool write()
Write the fluxSummary.
Definition: fluxSummary.C:1094
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
label faceId(-1)
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
List of names generated by calling name() for each list item and filtered for matches.
static void writeHeader(Ostream &os, const word &fieldName)
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalar totalArea(const label idx) const
Calculate the total area for the surface or derived faceZone.
Definition: fluxSummary.C:662
virtual bool read(const dictionary &dict)
Read the field fluxSummary data.
Definition: fluxSummary.C:966
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Tab [isspace].
Definition: token.H:129
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:517
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Fields (face and point) for polySurface.
Set of directions for each cell in the mesh. Either uniform and size=1 or one set of directions per c...
Definition: directions.H:63
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
fluxSummary(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: fluxSummary.C:940
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
bool update()
Initialise - after read(), before write()
Definition: fluxSummary.C:785
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
void initialiseCellZoneAndDirection(const word &cellZoneName, const vector &refDir, DynamicList< word > &names, DynamicList< vector > &dir, DynamicList< labelList > &faceID, DynamicList< labelList > &facePatchID, DynamicList< boolList > &faceFlip) const
Initialise face set from cell zone and direction.
Definition: fluxSummary.C:365
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field...
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.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
void setSize(const label n)
Alias for resize()
Definition: List.H:316
bool surfaceModeWrite()
Specialized write for surfaces.
Definition: fluxSummary.C:702
void initialiseFaceZoneAndDirection(const word &faceZoneName, const vector &refDir, DynamicList< word > &names, DynamicList< vector > &dir, DynamicList< labelList > &faceID, DynamicList< labelList > &facePatchID, DynamicList< boolList > &faceFlip) const
Initialise face set from face zone and direction.
Definition: fluxSummary.C:257
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
bool isSurfaceMode() const
Check if surface mode instead of zone mode.
Definition: fluxSummary.C:63
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
const wordList area
Standard area field types (scalar, vector, tensor, etc)
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static const word null
An empty word.
Definition: word.H:84
Foam::DimensionedField< vector, polySurfaceGeoMesh > polySurfaceVectorField
Vector< scalar > vector
Definition: vector.H:57
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:584
errorManip< error > abort(error &err)
Definition: errorManip.H:139
#define DebugInfo
Report an information message using Foam::Info.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
virtual bool execute()
Execute, currently does nothing.
Definition: fluxSummary.C:1088
edgeScalarField phis(IOobject("phis", runTime.timeName(), aMesh.thisDb(), IOobject::NO_READ, IOobject::NO_WRITE), linearEdgeInterpolate(Us) &aMesh.Le())
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
static const Enum< modeType > modeTypeNames_
Face mode names.
Definition: fluxSummary.H:212
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
mesh update()
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:241
#define WarningInFunction
Report a warning using Foam::Warning.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
#define Log
Definition: PDRblock.C:28
const std::string patch
OpenFOAM patch number as a std::string.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual void writeFileHeader(const word &zoneName, const scalar area, const vector &refDir, Ostream &os) const
Output file header information.
Definition: fluxSummary.C:1046
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
messageStream Info
Information stream (stdout output on master, null elsewhere)
constexpr label labelMax
Definition: label.H:55
virtual bool read(const dictionary &dict)
Read optional controls.
label n
fileName search(const word &file, const fileName &directory)
Recursively search the given directory for the file.
Definition: fileName.C:642
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
void initialiseSurface(const word &surfName, DynamicList< word > &names, DynamicList< vector > &dir, DynamicList< boolList > &faceFlip) const
Initialise for given surface name.
Definition: fluxSummary.C:100
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
List< label > labelList
A List of labels.
Definition: List.H:62
void initialiseFaceZone(const word &faceZoneName, DynamicList< word > &names, DynamicList< vector > &dir, DynamicList< labelList > &faceID, DynamicList< labelList > &facePatchID, DynamicList< boolList > &faceFlip) const
Initialise face set from face zone.
Definition: fluxSummary.C:174
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Base class for writing single files from the function objects.
Definition: writeFile.H:112
word checkFlowType(const dimensionSet &fieldDims, const word &fieldName) const
Check flowType (mass or volume)
Definition: fluxSummary.C:70
List< bool > boolList
A List of bools.
Definition: List.H:60
labelList cellIDs
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:485
void initialiseSurfaceAndDirection(const word &surfName, const vector &refDir, DynamicList< word > &names, DynamicList< vector > &dir, DynamicList< boolList > &faceFlip) const
Initialise for given surface name and direction.
Definition: fluxSummary.C:126
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127