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