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-2020 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,
109  DynamicList<word>& names,
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,
136  DynamicList<word>& names,
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,
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 
212  label faceID = -1;
213  label facePatchID = -1;
214  if (mesh_.isInternalFace(facei))
215  {
216  faceID = facei;
217  facePatchID = -1;
218  }
219  else
220  {
221  facePatchID = mesh_.boundaryMesh().whichPatch(facei);
222  const polyPatch& pp = mesh_.boundaryMesh()[facePatchID];
223  if (isA<coupledPolyPatch>(pp))
224  {
225  if (refCast<const coupledPolyPatch>(pp).owner())
226  {
227  faceID = pp.whichFace(facei);
228  }
229  else
230  {
231  faceID = -1;
232  }
233  }
234  else if (!isA<emptyPolyPatch>(pp))
235  {
236  faceID = facei - pp.start();
237  }
238  else
239  {
240  faceID = -1;
241  facePatchID = -1;
242  }
243  }
244 
245  if (faceID >= 0)
246  {
247  // Orientation set by faceZone flip map
248  if (fZone.flipMap()[i])
249  {
250  flips.append(true);
251  }
252  else
253  {
254  flips.append(false);
255  }
256 
257  faceIDs.append(faceID);
258  facePatchIDs.append(facePatchID);
259  }
260  }
261 
262  // could reduce some copying here
263  faceID.append(faceIDs);
264  facePatchID.append(facePatchIDs);
265  faceFlip.append(flips);
266 }
267 
268 
270 (
271  const word& faceZoneName,
272  const vector& dir,
273  DynamicList<word>& names,
275  DynamicList<labelList>& faceID,
276  DynamicList<labelList>& facePatchID,
277  DynamicList<boolList>& faceFlip
278 ) const
279 {
280  const vector refDir = dir/(mag(dir) + ROOTVSMALL);
281 
282  label zonei = mesh_.faceZones().findZoneID(faceZoneName);
283  if (zonei == -1)
284  {
286  << "Unable to find faceZone " << faceZoneName
287  << ". Valid zones: "
288  << mesh_.faceZones().sortedNames() << nl
289  << exit(FatalError);
290  }
291  const faceZone& fZone = mesh_.faceZones()[zonei];
292 
293  names.append(faceZoneName);
294  directions.append(refDir);
295 
296  DynamicList<label> faceIDs(fZone.size());
297  DynamicList<label> facePatchIDs(fZone.size());
298  DynamicList<bool> flips(fZone.size());
299 
300  const surfaceVectorField& Sf = mesh_.Sf();
301  const surfaceScalarField& magSf = mesh_.magSf();
302 
303  vector n(Zero);
304 
305  forAll(fZone, i)
306  {
307  label facei = fZone[i];
308 
309  label faceID = -1;
310  label facePatchID = -1;
311  if (mesh_.isInternalFace(facei))
312  {
313  faceID = facei;
314  facePatchID = -1;
315  }
316  else
317  {
318  facePatchID = mesh_.boundaryMesh().whichPatch(facei);
319  const polyPatch& pp = mesh_.boundaryMesh()[facePatchID];
320  if (isA<coupledPolyPatch>(pp))
321  {
322  if (refCast<const coupledPolyPatch>(pp).owner())
323  {
324  faceID = pp.whichFace(facei);
325  }
326  else
327  {
328  faceID = -1;
329  }
330  }
331  else if (!isA<emptyPolyPatch>(pp))
332  {
333  faceID = facei - pp.start();
334  }
335  else
336  {
337  faceID = -1;
338  facePatchID = -1;
339  }
340  }
341 
342  if (faceID >= 0)
343  {
344  // orientation set by comparison with reference direction
345  if (facePatchID != -1)
346  {
347  n = Sf.boundaryField()[facePatchID][faceID]
348  /(magSf.boundaryField()[facePatchID][faceID] + ROOTVSMALL);
349  }
350  else
351  {
352  n = Sf[faceID]/(magSf[faceID] + ROOTVSMALL);
353  }
354 
355  if ((n & refDir) > tolerance_)
356  {
357  flips.append(false);
358  }
359  else
360  {
361  flips.append(true);
362  }
363 
364  faceIDs.append(faceID);
365  facePatchIDs.append(facePatchID);
366  }
367  }
368 
369  // could reduce copying here
370  faceID.append(faceIDs);
371  facePatchID.append(facePatchIDs);
372  faceFlip.append(flips);
373 }
374 
375 
377 (
378  const word& cellZoneName,
379  const vector& dir,
380  DynamicList<word>& names,
382  DynamicList<labelList>& faceID,
383  DynamicList<labelList>& facePatchID,
384  DynamicList<boolList>& faceFlip
385 ) const
386 {
387  const vector refDir = dir/(mag(dir) + ROOTVSMALL);
388 
389  const label cellZonei = mesh_.cellZones().findZoneID(cellZoneName);
390  if (cellZonei == -1)
391  {
393  << "Unable to find cellZone " << cellZoneName
394  << ". Valid zones: "
395  << mesh_.cellZones().sortedNames() << nl
396  << exit(FatalError);
397  }
398 
399  const label nInternalFaces = mesh_.nInternalFaces();
400  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
401 
402  labelList cellAddr(mesh_.nCells(), -1);
403  const labelList& cellIDs = mesh_.cellZones()[cellZonei];
404  labelUIndList(cellAddr, cellIDs) = identity(cellIDs.size());
405  labelList nbrFaceCellAddr(mesh_.nFaces() - nInternalFaces, -1);
406 
407  forAll(pbm, patchi)
408  {
409  const polyPatch& pp = pbm[patchi];
410 
411  if (pp.coupled())
412  {
413  forAll(pp, i)
414  {
415  label facei = pp.start() + i;
416  label nbrFacei = facei - nInternalFaces;
417  label own = mesh_.faceOwner()[facei];
418  nbrFaceCellAddr[nbrFacei] = cellAddr[own];
419  }
420  }
421  }
422 
423  // Correct boundary values for parallel running
424  syncTools::swapBoundaryFaceList(mesh_, nbrFaceCellAddr);
425 
426  // Collect faces
427  DynamicList<label> faceIDs(floor(0.1*mesh_.nFaces()));
428  DynamicList<label> facePatchIDs(faceIDs.size());
429  DynamicList<label> faceLocalPatchIDs(faceIDs.size());
430  DynamicList<bool> flips(faceIDs.size());
431 
432  // Internal faces
433  for (label facei = 0; facei < nInternalFaces; facei++)
434  {
435  const label own = cellAddr[mesh_.faceOwner()[facei]];
436  const label nbr = cellAddr[mesh_.faceNeighbour()[facei]];
437 
438  if (((own != -1) && (nbr == -1)) || ((own == -1) && (nbr != -1)))
439  {
440  vector n = mesh_.faces()[facei].unitNormal(mesh_.points());
441 
442  if ((n & refDir) > tolerance_)
443  {
444  faceIDs.append(facei);
445  faceLocalPatchIDs.append(facei);
446  facePatchIDs.append(-1);
447  flips.append(false);
448  }
449  else if ((n & -refDir) > tolerance_)
450  {
451  faceIDs.append(facei);
452  faceLocalPatchIDs.append(facei);
453  facePatchIDs.append(-1);
454  flips.append(true);
455  }
456  }
457  }
458 
459  // Loop over boundary faces
460  forAll(pbm, patchi)
461  {
462  const polyPatch& pp = pbm[patchi];
463 
464  forAll(pp, localFacei)
465  {
466  const label facei = pp.start() + localFacei;
467  const label own = cellAddr[mesh_.faceOwner()[facei]];
468  const label nbr = nbrFaceCellAddr[facei - nInternalFaces];
469 
470  if ((own != -1) && (nbr == -1))
471  {
472  vector n = mesh_.faces()[facei].unitNormal(mesh_.points());
473 
474  if ((n & refDir) > tolerance_)
475  {
476  faceIDs.append(facei);
477  faceLocalPatchIDs.append(localFacei);
478  facePatchIDs.append(patchi);
479  flips.append(false);
480  }
481  else if ((n & -refDir) > tolerance_)
482  {
483  faceIDs.append(facei);
484  faceLocalPatchIDs.append(localFacei);
485  facePatchIDs.append(patchi);
486  flips.append(true);
487  }
488  }
489  }
490  }
491 
492  // Convert into primitivePatch for convenience
494  (
495  IndirectList<face>(mesh_.faces(), faceIDs),
496  mesh_.points()
497  );
498 
499  if (debug)
500  {
501  OBJstream os(mesh_.time().path()/"patch.obj");
502  faceList faces(patch);
503  os.write(faces, mesh_.points(), false);
504  }
505 
506 
507  // Data on all edges and faces
508  List<edgeTopoDistanceData<label>> allEdgeInfo(patch.nEdges());
509  List<edgeTopoDistanceData<label>> allFaceInfo(patch.size());
510 
511  bool search = true;
512 
513  DebugInfo
514  << "initialiseCellZoneAndDirection: "
515  << "Starting walk to split patch into faceZones"
516  << endl;
517 
518  globalIndex globalFaces(patch.size());
519 
520  label oldFaceID = 0;
521  label regioni = 0;
522 
523  // Dummy tracking data
524  bool dummyData{false};
525 
526  while (search)
527  {
528  DynamicList<label> changedEdges;
530 
531  label seedFacei = labelMax;
532  for (; oldFaceID < patch.size(); oldFaceID++)
533  {
534  if (!allFaceInfo[oldFaceID].valid<bool>(dummyData))
535  {
536  seedFacei = globalFaces.toGlobal(oldFaceID);
537  break;
538  }
539  }
540  reduce(seedFacei, minOp<label>());
541 
542  if (seedFacei == labelMax)
543  {
544  break;
545  }
546 
547  if (globalFaces.isLocal(seedFacei))
548  {
549  const label localFacei = globalFaces.toLocal(seedFacei);
550  const labelList& fEdges = patch.faceEdges()[localFacei];
551 
552  for (const label edgei : fEdges)
553  {
554  if (allEdgeInfo[edgei].valid<bool>(dummyData))
555  {
557  << "Problem in edge face wave: attempted to assign a "
558  << "value to an edge that has already been visited. "
559  << "Edge info: " << allEdgeInfo[edgei]
560  << endl;
561  }
562 
563  changedEdges.append(edgei);
564  changedInfo.append
565  (
567  (
568  0, // distance
569  regioni
570  )
571  );
572  }
573  }
574 
575 
577  <
580  > calc
581  (
582  mesh_,
583  patch,
584  changedEdges,
585  changedInfo,
586  allEdgeInfo,
587  allFaceInfo,
588  returnReduce(patch.nEdges(), sumOp<label>())
589  );
590 
591  if (debug)
592  {
593  label nCells = 0;
594  forAll(allFaceInfo, facei)
595  {
596  if (allFaceInfo[facei].data() == regioni)
597  {
598  nCells++;
599  }
600  }
601 
602  Info<< "*** region:" << regioni
603  << " found:" << returnReduce(nCells, sumOp<label>())
604  << " faces" << endl;
605  }
606 
607  regioni++;
608  }
609 
610  // Collect the data per region
611  const label nRegion = regioni;
612 
613  #ifdef FULLDEBUG
614  // May wish to have enabled always
615  if (nRegion == 0)
616  {
618  << "Region split failed" << nl
619  << exit(FatalError);
620  }
621  #endif
622 
623  List<DynamicList<label>> regionFaceIDs(nRegion);
624  List<DynamicList<label>> regionFacePatchIDs(nRegion);
625  List<DynamicList<bool>> regionFaceFlips(nRegion);
626 
627  forAll(allFaceInfo, facei)
628  {
629  regioni = allFaceInfo[facei].data();
630 
631  regionFaceIDs[regioni].append(faceLocalPatchIDs[facei]);
632  regionFacePatchIDs[regioni].append(facePatchIDs[facei]);
633  regionFaceFlips[regioni].append(flips[facei]);
634  }
635 
636  // Transfer to persistent storage
637  forAll(regionFaceIDs, regioni)
638  {
639  const word zoneName = cellZoneName + ":faceZone" + Foam::name(regioni);
640  names.append(zoneName);
641  directions.append(refDir);
642  faceID.append(regionFaceIDs[regioni]);
643  facePatchID.append(regionFacePatchIDs[regioni]);
644  faceFlip.append(regionFaceFlips[regioni]);
645 
646  // Write OBj of faces to file
647  if (debug)
648  {
649  OBJstream os(mesh_.time().path()/zoneName + ".obj");
650  faceList faces(mesh_.faces(), regionFaceIDs[regioni]);
651  os.write(faces, mesh_.points(), false);
652  }
653  }
654 
655  if (log)
656  {
657  Info<< type() << " " << name() << " output:" << nl
658  << " Created " << faceID.size()
659  << " separate face zones from cell zone " << cellZoneName << nl;
660 
661  forAll(names, i)
662  {
663  label nFaces = returnReduce(faceID[i].size(), sumOp<label>());
664  Info<< " " << names[i] << ": "
665  << nFaces << " faces" << nl;
666  }
667 
668  Info<< endl;
669  }
670 }
671 
672 
674 (
675  const label idx
676 ) const
677 {
678  scalar sumMagSf = 0;
679 
680  if (isSurfaceMode())
681  {
682  const polySurface& s =
683  storedObjects().lookupObject<polySurface>(zoneNames_[idx]);
684 
685  sumMagSf = sum(s.magSf());
686  }
687  else
688  {
689  const surfaceScalarField& magSf = mesh_.magSf();
690 
691  const labelList& faceIDs = faceID_[idx];
692  const labelList& facePatchIDs = facePatchID_[idx];
693 
694  forAll(faceIDs, i)
695  {
696  label facei = faceIDs[i];
697 
698  if (facePatchIDs[i] == -1)
699  {
700  sumMagSf += magSf[facei];
701  }
702  else
703  {
704  label patchi = facePatchIDs[i];
705  sumMagSf += magSf.boundaryField()[patchi][facei];
706  }
707  }
708  }
709 
710  return returnReduce(sumMagSf, sumOp<scalar>());
711 }
712 
713 
715 {
716  for (const word& surfName : zoneNames_)
717  {
718  const polySurface& s =
719  storedObjects().lookupObject<polySurface>(surfName);
720 
721  const auto& phi = s.lookupObject<polySurfaceVectorField>(phiName_);
722 
723  Log << type() << ' ' << name() << ' '
724  << checkFlowType(phi.dimensions(), phi.name()) << " write:" << nl;
725  }
726 
727 
728  forAll(zoneNames_, surfi)
729  {
730  const polySurface& s =
731  storedObjects().lookupObject<polySurface>(zoneNames_[surfi]);
732 
733  const auto& phi = s.lookupObject<polySurfaceVectorField>(phiName_);
734 
735  checkFlowType(phi.dimensions(), phi.name());
736 
737  const boolList& flips = faceFlip_[surfi];
738 
739  scalar phiPos(0);
740  scalar phiNeg(0);
741 
742  tmp<scalarField> tphis = phi & s.Sf();
743  const scalarField& phis = tphis();
744 
745  forAll(s, i)
746  {
747  scalar phif = phis[i];
748  if (flips[i])
749  {
750  phif *= -1;
751  }
752 
753  if (phif > 0)
754  {
755  phiPos += phif;
756  }
757  else
758  {
759  phiNeg += phif;
760  }
761  }
762 
763  reduce(phiPos, sumOp<scalar>());
764  reduce(phiNeg, sumOp<scalar>());
765 
766  phiPos *= scaleFactor_;
767  phiNeg *= scaleFactor_;
768 
769  scalar netFlux = phiPos + phiNeg;
770  scalar absoluteFlux = phiPos - phiNeg;
771 
772  Log << " surface " << zoneNames_[surfi] << ':' << nl
773  << " positive : " << phiPos << nl
774  << " negative : " << phiNeg << nl
775  << " net : " << netFlux << nl
776  << " absolute : " << absoluteFlux
777  << nl << endl;
778 
779  if (writeToFile())
780  {
781  filePtrs_[surfi]
782  << time_.value() << token::TAB
783  << phiPos << token::TAB
784  << phiNeg << token::TAB
785  << netFlux << token::TAB
786  << absoluteFlux
787  << endl;
788  }
789  }
790 
791  Log << endl;
792 
793  return true;
794 }
795 
796 
798 {
799  if (!needsUpdate_)
800  {
801  return false;
802  }
803 
804  // Initialise with capacity == number of input names
805  DynamicList<word> faceZoneName(zoneNames_.size());
806  DynamicList<vector> refDir(faceZoneName.capacity());
807  DynamicList<labelList> faceID(faceZoneName.capacity());
808  DynamicList<labelList> facePatchID(faceZoneName.capacity());
809  DynamicList<boolList> faceFlips(faceZoneName.capacity());
810 
811  switch (mode_)
812  {
813  case mdFaceZone:
814  {
815  forAll(zoneNames_, zonei)
816  {
817  initialiseFaceZone
818  (
819  zoneNames_[zonei],
820  faceZoneName,
821  refDir, // fill with dummy value
822  faceID,
823  facePatchID,
824  faceFlips
825  );
826  }
827  break;
828  }
829  case mdFaceZoneAndDirection:
830  {
831  forAll(zoneNames_, zonei)
832  {
833  initialiseFaceZoneAndDirection
834  (
835  zoneNames_[zonei],
836  zoneDirections_[zonei],
837  faceZoneName,
838  refDir,
839  faceID,
840  facePatchID,
841  faceFlips
842  );
843  }
844  break;
845  }
846  case mdCellZoneAndDirection:
847  {
848  forAll(zoneNames_, zonei)
849  {
850  initialiseCellZoneAndDirection
851  (
852  zoneNames_[zonei],
853  zoneDirections_[zonei],
854  faceZoneName,
855  refDir,
856  faceID,
857  facePatchID,
858  faceFlips
859  );
860  }
861  break;
862  }
863  case mdSurface:
864  {
865  forAll(zoneNames_, zonei)
866  {
867  initialiseSurface
868  (
869  zoneNames_[zonei],
870  faceZoneName,
871  refDir,
872  faceFlips
873  );
874  }
875  break;
876  }
877  case mdSurfaceAndDirection:
878  {
879  forAll(zoneNames_, zonei)
880  {
881  initialiseSurfaceAndDirection
882  (
883  zoneNames_[zonei],
884  zoneDirections_[zonei],
885  faceZoneName,
886  refDir,
887  faceFlips
888  );
889  }
890  break;
891  }
892 
893  // Compiler warning if we forgot an enumeration
894  }
895 
896  zoneNames_.transfer(faceZoneName);
897  faceID_.transfer(faceID);
898  facePatchID_.transfer(facePatchID);
899  faceFlip_.transfer(faceFlips);
900 
901  Info<< type() << " " << name() << " output:" << nl;
902 
903  // Calculate and report areas
904  List<scalar> areas(zoneNames_.size());
905  forAll(zoneNames_, zonei)
906  {
907  const word& zoneName = zoneNames_[zonei];
908  areas[zonei] = totalArea(zonei);
909 
910  if (isSurfaceMode())
911  {
912  Info<< " Surface: " << zoneName
913  << ", area: " << areas[zonei] << nl;
914  }
915  else
916  {
917  Info<< " Zone: " << zoneName
918  << ", area: " << areas[zonei] << nl;
919  }
920  }
921  Info<< endl;
922 
923  if (writeToFile())
924  {
925  filePtrs_.resize(zoneNames_.size());
926 
927  forAll(filePtrs_, zonei)
928  {
929  const word& zoneName = zoneNames_[zonei];
930  filePtrs_.set(zonei, createFile(zoneName));
931  writeFileHeader
932  (
933  zoneName,
934  areas[zonei],
935  refDir[zonei],
936  filePtrs_[zonei]
937  );
938  }
939  }
940 
941  Info<< endl;
942 
943  needsUpdate_ = false;
944 
945  return true;
946 }
947 
948 
949 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
950 
952 (
953  const word& name,
954  const Time& runTime,
955  const dictionary& dict
956 )
957 :
959  writeFile(obr_, name),
960  needsUpdate_(true),
961  mode_(mdFaceZone),
962  scaleFactor_(1),
963  phiName_("phi"),
964  zoneNames_(),
965  zoneDirections_(),
966  faceID_(),
967  facePatchID_(),
968  faceFlip_(),
969  filePtrs_(),
970  tolerance_(0.8)
971 {
972  read(dict);
973 }
974 
975 
976 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
977 
979 {
982 
983  needsUpdate_ = true;
984  mode_ = modeTypeNames_.get("mode", dict);
985  phiName_ = dict.getOrDefault<word>("phi", "phi");
986  scaleFactor_ = dict.getOrDefault<scalar>("scaleFactor", 1);
987  tolerance_ = dict.getOrDefault<scalar>("tolerance", 0.8);
988 
989  zoneNames_.clear();
990  zoneDirections_.clear();
991 
992  List<Tuple2<word, vector>> nameAndDirection;
993 
994  switch (mode_)
995  {
996  case mdFaceZone:
997  {
998  dict.readEntry("faceZones", zoneNames_);
999  break;
1000  }
1001  case mdFaceZoneAndDirection:
1002  {
1003  dict.readEntry("faceZoneAndDirection", nameAndDirection);
1004  break;
1005  }
1006  case mdCellZoneAndDirection:
1007  {
1008  dict.readEntry("cellZoneAndDirection", nameAndDirection);
1009  break;
1010  }
1011  case mdSurface:
1012  {
1013  dict.readEntry("surfaces", zoneNames_);
1014  break;
1015  }
1016  case mdSurfaceAndDirection:
1017  {
1018  dict.readEntry("surfaceAndDirection", nameAndDirection);
1019  break;
1020  }
1021  default:
1022  {
1024  << "unhandled enumeration " << modeTypeNames_[mode_]
1025  << abort(FatalIOError);
1026  }
1027  }
1028 
1029 
1030  // Split name/vector into separate lists
1031  if (nameAndDirection.size())
1032  {
1033  zoneNames_.resize(nameAndDirection.size());
1034  zoneDirections_.resize(nameAndDirection.size());
1035 
1036  label zonei = 0;
1037 
1038  for (const Tuple2<word, vector>& nameDirn : nameAndDirection)
1039  {
1040  zoneNames_[zonei] = nameDirn.first();
1041  zoneDirections_[zonei] = nameDirn.second();
1042  ++zonei;
1043  }
1044 
1045  nameAndDirection.clear();
1046  }
1047 
1048 
1049  Info<< type() << ' ' << name() << " ("
1050  << modeTypeNames_[mode_] << ") with selection:\n "
1051  << flatOutput(zoneNames_) << endl;
1052 
1053  return !zoneNames_.empty();
1054 }
1055 
1056 
1059  const word& zoneName,
1060  const scalar area,
1061  const vector& refDir,
1062  Ostream& os
1063 ) const
1064 {
1065  writeHeader(os, "Flux summary");
1066  if (isSurfaceMode())
1067  {
1068  writeHeaderValue(os, "Surface", zoneName);
1069  }
1070  else
1071  {
1072  writeHeaderValue(os, "Face zone", zoneName);
1073  }
1074  writeHeaderValue(os, "Total area", area);
1075 
1076  switch (mode_)
1077  {
1078  case mdFaceZoneAndDirection:
1079  case mdCellZoneAndDirection:
1080  case mdSurfaceAndDirection:
1081  {
1082  writeHeaderValue(os, "Reference direction", refDir);
1083  break;
1084  }
1085  default:
1086  {}
1087  }
1088 
1089  writeHeaderValue(os, "Scale factor", scaleFactor_);
1090 
1091  writeCommented(os, "Time");
1092  os << tab << "positive"
1093  << tab << "negative"
1094  << tab << "net"
1095  << tab << "absolute"
1096  << endl;
1097 }
1098 
1099 
1101 {
1102  return true;
1103 }
1104 
1105 
1107 {
1108  update();
1109 
1110  if (isSurfaceMode())
1111  {
1112  return surfaceModeWrite();
1113  }
1114 
1115  const surfaceScalarField& phi = lookupObject<surfaceScalarField>(phiName_);
1116 
1117  Log << type() << ' ' << name() << ' '
1118  << checkFlowType(phi.dimensions(), phi.name()) << " write:" << nl;
1119 
1120  forAll(zoneNames_, zonei)
1121  {
1122  const labelList& faceID = faceID_[zonei];
1123  const labelList& facePatchID = facePatchID_[zonei];
1124  const boolList& faceFlips = faceFlip_[zonei];
1125 
1126  scalar phiPos(0);
1127  scalar phiNeg(0);
1128  scalar phif(0);
1129 
1130  forAll(faceID, i)
1131  {
1132  label facei = faceID[i];
1133  label patchi = facePatchID[i];
1134 
1135  if (patchi != -1)
1136  {
1137  phif = phi.boundaryField()[patchi][facei];
1138  }
1139  else
1140  {
1141  phif = phi[facei];
1142  }
1143 
1144  if (faceFlips[i])
1145  {
1146  phif *= -1;
1147  }
1148 
1149  if (phif > 0)
1150  {
1151  phiPos += phif;
1152  }
1153  else
1154  {
1155  phiNeg += phif;
1156  }
1157  }
1158 
1159  reduce(phiPos, sumOp<scalar>());
1160  reduce(phiNeg, sumOp<scalar>());
1161 
1162  phiPos *= scaleFactor_;
1163  phiNeg *= scaleFactor_;
1164 
1165  scalar netFlux = phiPos + phiNeg;
1166  scalar absoluteFlux = phiPos - phiNeg;
1167 
1168  Log << " faceZone " << zoneNames_[zonei] << ':' << nl
1169  << " positive : " << phiPos << nl
1170  << " negative : " << phiNeg << nl
1171  << " net : " << netFlux << nl
1172  << " absolute : " << absoluteFlux
1173  << nl << endl;
1174 
1175  if (writeToFile())
1176  {
1177  filePtrs_[zonei]
1178  << time_.value() << token::TAB
1179  << phiPos << token::TAB
1180  << phiNeg << token::TAB
1181  << netFlux << token::TAB
1182  << absoluteFlux
1183  << endl;
1184  }
1185  }
1186 
1187  Log << endl;
1188 
1189  return true;
1190 }
1191 
1192 
1193 // ************************************************************************* //
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:377
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:62
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:270
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:62
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:714
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:329
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:439
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:69
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:182
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
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:49
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
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:1100
Foam::functionObjects::fluxSummary::totalArea
scalar totalArea(const label idx) const
Calculate the total area for the surface or derived faceZone.
Definition: fluxSummary.C:674
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:54
Foam::functionObjects::writeFile::read
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:214
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:1058
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
Foam::Field< scalar >
Foam::polyPatch::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:307
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:65
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::functionObjects::fluxSummary::fluxSummary
fluxSummary(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: fluxSummary.C:952
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:474
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:1106
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:121
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
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::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
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:217
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
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:313
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:381
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:359
Foam::polyPatch::whichFace
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:400
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:384
Foam::nl
constexpr char nl
Definition: Ostream.H:385
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:978
Foam::token::TAB
Tab [isspace].
Definition: token.H:118
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:576
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
Foam::fieldTypes::area
const wordList area
Standard area field types (scalar, vector, tensor, etc)
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:401
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: Tuple2.H:57
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::GeometricField< vector, fvsPatchField, surfaceMesh >
Foam::faceZone::flipMap
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:271
Foam::functionObjects::fluxSummary::update
bool update()
Initialise - after read(), before write()
Definition: fluxSummary.C:797
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
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:85
Foam::labelUIndList
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: UIndirectList.H:58