isoAdvection.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) 2016-2017 DHI
9  Modified code Copyright (C) 2016-2017 OpenCFD Ltd.
10  Modified code Copyright (C) 2019-2020 DLR
11  Modified code Copyright (C) 2018, 2021 Johan Roenby
12 -------------------------------------------------------------------------------
13 License
14  This file is part of OpenFOAM.
15 
16  OpenFOAM is free software: you can redistribute it and/or modify it
17  under the terms of the GNU General Public License as published by
18  the Free Software Foundation, either version 3 of the License, or
19  (at your option) any later version.
20 
21  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24  for more details.
25 
26  You should have received a copy of the GNU General Public License
27  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28 
29 \*---------------------------------------------------------------------------*/
30 
31 #include "isoAdvection.H"
32 #include "volFields.H"
33 #include "interpolationCellPoint.H"
34 #include "volPointInterpolation.H"
35 #include "fvcSurfaceIntegrate.H"
36 #include "fvcGrad.H"
37 #include "upwind.H"
38 #include "cellSet.H"
39 #include "meshTools.H"
40 #include "OBJstream.H"
41 #include "syncTools.H"
42 #include "profiling.H"
43 
45 
46 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50  defineTypeNameAndDebug(isoAdvection, 0);
51 }
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
55 Foam::isoAdvection::isoAdvection
56 (
58  const surfaceScalarField& phi,
59  const volVectorField& U
60 )
61 :
62  mesh_(alpha1.mesh()),
63  dict_(mesh_.solverDict(alpha1.name())),
64  alpha1_(alpha1),
65  alpha1In_(alpha1.ref()),
66  phi_(phi),
67  U_(U),
68  dVf_
69  (
70  IOobject
71  (
72  "dVf_",
73  mesh_.time().timeName(),
74  mesh_,
75  IOobject::NO_READ,
76  IOobject::NO_WRITE
77  ),
78  mesh_,
80  ),
81  alphaPhi_
82  (
83  IOobject
84  (
85  "alphaPhi_",
86  mesh_.time().timeName(),
87  mesh_,
88  IOobject::NO_READ,
89  IOobject::NO_WRITE
90  ),
91  mesh_,
93  ),
94  advectionTime_(0),
95  timeIndex_(-1),
96 
97  // Tolerances and solution controls
98  nAlphaBounds_(dict_.getOrDefault<label>("nAlphaBounds", 3)),
99  isoFaceTol_(dict_.getOrDefault<scalar>("isoFaceTol", 1e-8)),
100  surfCellTol_(dict_.getOrDefault<scalar>("surfCellTol", 1e-8)),
101  writeIsoFacesToFile_(dict_.getOrDefault("writeIsoFaces", false)),
102 
103  // Cell cutting data
104  surfCells_(label(0.2*mesh_.nCells())),
105  surf_(reconstructionSchemes::New(alpha1_, phi_, U_, dict_)),
106  advectFace_(alpha1.mesh(), alpha1),
107  bsFaces_(label(0.2*mesh_.nBoundaryFaces())),
108  bsx0_(bsFaces_.size()),
109  bsn0_(bsFaces_.size()),
110  bsUn0_(bsFaces_.size()),
111 
112  // Porosity
113  porosityEnabled_(dict_.getOrDefault<bool>("porosityEnabled", false)),
114  porosityPtr_(nullptr),
115 
116  // Parallel run data
117  procPatchLabels_(mesh_.boundary().size()),
118  surfaceCellFacesOnProcPatches_(0)
119 {
121 
122  // Prepare lists used in parallel runs
123  if (Pstream::parRun())
124  {
125  // Force calculation of required demand driven data (else parallel
126  // communication may crash)
127  mesh_.cellCentres();
128  mesh_.cellVolumes();
129  mesh_.faceCentres();
130  mesh_.faceAreas();
131  mesh_.magSf();
132  mesh_.boundaryMesh().patchID();
133  mesh_.cellPoints();
134  mesh_.cellCells();
135  mesh_.cells();
136 
137  // Get boundary mesh and resize the list for parallel comms
138  setProcessorPatches();
139  }
140 
141  // Reading porosity properties from constant directory
142  IOdictionary porosityProperties
143  (
144  IOobject
145  (
146  "porosityProperties",
147  mesh_.time().constant(),
148  mesh_,
151  )
152  );
153 
154  porosityEnabled_ =
155  porosityProperties.getOrDefault<bool>("porosityEnabled", false);
156 
157  if (porosityEnabled_)
158  {
159  if (mesh_.foundObject<volScalarField>("porosity"))
160  {
161  porosityPtr_ = mesh_.getObjectPtr<volScalarField>("porosity");
162 
163  if
164  (
165  gMin(porosityPtr_->primitiveField()) <= 0
166  || gMax(porosityPtr_->primitiveField()) > 1 + SMALL
167  )
168  {
170  << "Porosity field has values <= 0 or > 1"
171  << exit(FatalError);
172  }
173  }
174  else
175  {
177  << "Porosity enabled in constant/porosityProperties "
178  << "but no porosity field is found in object registry."
179  << exit(FatalError);
180  }
181  }
182 }
183 
184 
185 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
186 
187 void Foam::isoAdvection::setProcessorPatches()
188 {
189  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
190  surfaceCellFacesOnProcPatches_.clear();
191  surfaceCellFacesOnProcPatches_.resize(patches.size());
192 
193  // Append all processor patch labels to the list
194  procPatchLabels_.clear();
195  forAll(patches, patchi)
196  {
197  if
198  (
199  isA<processorPolyPatch>(patches[patchi])
200  && !patches[patchi].empty()
201  )
202  {
203  procPatchLabels_.append(patchi);
204  }
205  }
206 }
207 
208 
209 void Foam::isoAdvection::extendMarkedCells
210 (
211  bitSet& markedCell
212 ) const
213 {
214  // Mark faces using any marked cell
215  bitSet markedFace(mesh_.nFaces());
216 
217  for (const label celli : markedCell)
218  {
219  markedFace.set(mesh_.cells()[celli]); // set multiple faces
220  }
221 
222  syncTools::syncFaceList(mesh_, markedFace, orEqOp<unsigned int>());
223 
224  // Update cells using any markedFace
225  for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
226  {
227  if (markedFace.test(facei))
228  {
229  markedCell.set(mesh_.faceOwner()[facei]);
230  markedCell.set(mesh_.faceNeighbour()[facei]);
231  }
232  }
233  for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); ++facei)
234  {
235  if (markedFace.test(facei))
236  {
237  markedCell.set(mesh_.faceOwner()[facei]);
238  }
239  }
240 }
241 
242 
243 void Foam::isoAdvection::timeIntegratedFlux()
244 {
245  addProfilingInFunction(geometricVoF);
246  // Get time step
247  const scalar dt = mesh_.time().deltaTValue();
248 
249  // Create object for interpolating velocity to isoface centres
250  interpolationCellPoint<vector> UInterp(U_);
251 
252  // For each downwind face of each surface cell we "isoadvect" to find dVf
253  label nSurfaceCells = 0;
254 
255  // Clear out the data for re-use and reset list containing information
256  // whether cells could possibly need bounding
257  clearIsoFaceData();
258 
259  // Get necessary references
260  const scalarField& phiIn = phi_.primitiveField();
261  const scalarField& magSfIn = mesh_.magSf().primitiveField();
262  scalarField& dVfIn = dVf_.primitiveFieldRef();
263 
264  // Get necessary mesh data
265  const cellList& cellFaces = mesh_.cells();
266  const labelList& own = mesh_.faceOwner();
267 
268 
269  // Storage for isoFace points. Only used if writeIsoFacesToFile_
270  DynamicList<List<point>> isoFacePts;
271  const DynamicField<label>& interfaceLabels = surf_->interfaceLabels();
272 
273  // Calculating isoface normal velocity
274  scalarField Un0(interfaceLabels.size());
275  forAll(Un0, i)
276  {
277  const label celli = interfaceLabels[i];
278  const point x0(surf_->centre()[celli]);
279  const vector n0(normalised(-surf_->normal()[celli]));
280  Un0[i] = UInterp.interpolate(x0, celli) & n0;
281  }
282 
283  // Taking acount of porosity if enabled
284  if (porosityEnabled_)
285  {
286  forAll(Un0, i)
287  {
288  const label celli = interfaceLabels[i];
289  Un0[i] /= porosityPtr_->primitiveField()[celli];
290  }
291  }
292 
293  // Loop through cells
294  forAll(interfaceLabels, i)
295  {
296  const label celli = interfaceLabels[i];
297  if (mag(surf_->normal()[celli]) != 0)
298  {
299 
300  // This is a surface cell, increment counter, append and mark cell
301  nSurfaceCells++;
302  surfCells_.append(celli);
303  const point x0(surf_->centre()[celli]);
304  const vector n0(normalised(-surf_->normal()[celli]));
305 
306  DebugInfo
307  << "\n------------ Cell " << celli << " with alpha1 = "
308  << alpha1In_[celli] << " and 1-alpha1 = "
309  << 1.0 - alpha1In_[celli] << " ------------"
310  << endl;
311 
312  // Estimate time integrated flux through each downwind face
313  // Note: looping over all cell faces - in reduced-D, some of
314  // these faces will be on empty patches
315  const cell& celliFaces = cellFaces[celli];
316  for (const label facei : celliFaces)
317  {
318  if (mesh_.isInternalFace(facei))
319  {
320  bool isDownwindFace = false;
321 
322  if (celli == own[facei])
323  {
324  if (phiIn[facei] >= 0)
325  {
326  isDownwindFace = true;
327  }
328  }
329  else
330  {
331  if (phiIn[facei] < 0)
332  {
333  isDownwindFace = true;
334  }
335  }
336 
337  if (isDownwindFace)
338  {
339  dVfIn[facei] = advectFace_.timeIntegratedFaceFlux
340  (
341  facei,
342  x0,
343  n0,
344  Un0[i],
345  dt,
346  phiIn[facei],
347  magSfIn[facei]
348  );
349  }
350 
351  }
352  else
353  {
354  bsFaces_.append(facei);
355  bsx0_.append(x0);
356  bsn0_.append(n0);
357  bsUn0_.append(Un0[i]);
358 
359  // Note: we must not check if the face is on the
360  // processor patch here.
361  }
362  }
363  }
364  }
365 
366  // Get references to boundary fields
367  const polyBoundaryMesh& boundaryMesh = mesh_.boundaryMesh();
368  const surfaceScalarField::Boundary& phib = phi_.boundaryField();
369  const surfaceScalarField::Boundary& magSfb = mesh_.magSf().boundaryField();
370  surfaceScalarField::Boundary& dVfb = dVf_.boundaryFieldRef();
371  const label nInternalFaces = mesh_.nInternalFaces();
372 
373  // Loop through boundary surface faces
374  forAll(bsFaces_, i)
375  {
376  // Get boundary face index (in the global list)
377  const label facei = bsFaces_[i];
378  const label patchi = boundaryMesh.patchID()[facei - nInternalFaces];
379  const label start = boundaryMesh[patchi].start();
380 
381  if (!phib[patchi].empty())
382  {
383  const label patchFacei = facei - start;
384  const scalar phiP = phib[patchi][patchFacei];
385 
386  if (phiP >= 0)
387  {
388  const scalar magSf = magSfb[patchi][patchFacei];
389 
390  dVfb[patchi][patchFacei] = advectFace_.timeIntegratedFaceFlux
391  (
392  facei,
393  bsx0_[i],
394  bsn0_[i],
395  bsUn0_[i],
396  dt,
397  phiP,
398  magSf
399  );
400 
401  // Check if the face is on processor patch and append it to
402  // the list if necessary
403  checkIfOnProcPatch(facei);
404  }
405  }
406  }
407 
408  // Synchronize processor patches
409  syncProcPatches(dVf_, phi_);
410 
411  writeIsoFaces(isoFacePts);
412 
413  DebugInfo << "Number of isoAdvector surface cells = "
414  << returnReduce(nSurfaceCells, sumOp<label>()) << endl;
415 }
416 
417 
418 void Foam::isoAdvection::setDownwindFaces
419 (
420  const label celli,
421  DynamicLabelList& downwindFaces
422 ) const
423 {
425 
426  // Get necessary mesh data and cell information
427  const labelList& own = mesh_.faceOwner();
428  const cellList& cells = mesh_.cells();
429  const cell& c = cells[celli];
430 
431  downwindFaces.clear();
432 
433  // Check all faces of the cell
434  for (const label facei: c)
435  {
436  // Get face and corresponding flux
437  const scalar phi = faceValue(phi_, facei);
438 
439  if (own[facei] == celli)
440  {
441  if (phi >= 0)
442  {
443  downwindFaces.append(facei);
444  }
445  }
446  else if (phi < 0)
447  {
448  downwindFaces.append(facei);
449  }
450  }
451 
452  downwindFaces.shrink();
453 }
454 
455 
456 Foam::scalar Foam::isoAdvection::netFlux
457 (
458  const surfaceScalarField& dVf,
459  const label celli
460 ) const
461 {
462  scalar dV = 0;
463 
464  // Get face indices
465  const cell& c = mesh_.cells()[celli];
466 
467  // Get mesh data
468  const labelList& own = mesh_.faceOwner();
469 
470  for (const label facei : c)
471  {
472  const scalar dVff = faceValue(dVf, facei);
473 
474  if (own[facei] == celli)
475  {
476  dV += dVff;
477  }
478  else
479  {
480  dV -= dVff;
481  }
482  }
483 
484  return dV;
485 }
486 
487 
488 Foam::DynamicList<Foam::label> Foam::isoAdvection::syncProcPatches
489 (
490  surfaceScalarField& dVf,
491  const surfaceScalarField& phi,
492  bool returnSyncedFaces
493 )
494 {
495  DynamicLabelList syncedFaces(0);
496  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
497 
498  if (Pstream::parRun())
499  {
500  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
501 
502  // Send
503  for (const label patchi : procPatchLabels_)
504  {
505  const processorPolyPatch& procPatch =
506  refCast<const processorPolyPatch>(patches[patchi]);
507 
508  UOPstream toNbr(procPatch.neighbProcNo(), pBufs);
509  const scalarField& pFlux = dVf.boundaryField()[patchi];
510 
511  const List<label>& surfCellFacesOnProcPatch =
512  surfaceCellFacesOnProcPatches_[patchi];
513 
514  const UIndirectList<scalar> dVfPatch
515  (
516  pFlux,
517  surfCellFacesOnProcPatch
518  );
519 
520  toNbr << surfCellFacesOnProcPatch << dVfPatch;
521  }
522 
523  pBufs.finishedSends();
524 
525 
526  // Receive and combine
527  for (const label patchi : procPatchLabels_)
528  {
529  const processorPolyPatch& procPatch =
530  refCast<const processorPolyPatch>(patches[patchi]);
531 
532  UIPstream fromNeighb(procPatch.neighbProcNo(), pBufs);
533  List<label> faceIDs;
534  List<scalar> nbrdVfs;
535 
536  fromNeighb >> faceIDs >> nbrdVfs;
537  if (returnSyncedFaces)
538  {
539  List<label> syncedFaceI(faceIDs);
540  for (label& faceI : syncedFaceI)
541  {
542  faceI += procPatch.start();
543  }
544  syncedFaces.append(syncedFaceI);
545  }
546 
547  if (debug)
548  {
549  Pout<< "Received at time = " << mesh_.time().value()
550  << ": surfCellFacesOnProcPatch = " << faceIDs << nl
551  << "Received at time = " << mesh_.time().value()
552  << ": dVfPatch = " << nbrdVfs << endl;
553  }
554 
555  // Combine fluxes
556  scalarField& localFlux = dVf.boundaryFieldRef()[patchi];
557 
558  forAll(faceIDs, i)
559  {
560  const label facei = faceIDs[i];
561  localFlux[facei] = - nbrdVfs[i];
562  if (debug && mag(localFlux[facei] + nbrdVfs[i]) > ROOTVSMALL)
563  {
564  Pout<< "localFlux[facei] = " << localFlux[facei]
565  << " and nbrdVfs[i] = " << nbrdVfs[i]
566  << " for facei = " << facei << endl;
567  }
568  }
569  }
570 
571  if (debug)
572  {
573  // Write out results for checking
574  forAll(procPatchLabels_, patchLabeli)
575  {
576  const label patchi = procPatchLabels_[patchLabeli];
577  const scalarField& localFlux = dVf.boundaryField()[patchi];
578  Pout<< "time = " << mesh_.time().value() << ": localFlux = "
579  << localFlux << endl;
580  }
581  }
582 
583  // Reinitialising list used for minimal parallel communication
584  forAll(surfaceCellFacesOnProcPatches_, patchi)
585  {
586  surfaceCellFacesOnProcPatches_[patchi].clear();
587  }
588  }
589 
590  return syncedFaces;
591 }
592 
593 
594 void Foam::isoAdvection::checkIfOnProcPatch(const label facei)
595 {
596  if (!mesh_.isInternalFace(facei))
597  {
598  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
599  const label patchi = pbm.patchID()[facei - mesh_.nInternalFaces()];
600 
601  if (isA<processorPolyPatch>(pbm[patchi]) && !pbm[patchi].empty())
602  {
603  const label patchFacei = pbm[patchi].whichFace(facei);
604  surfaceCellFacesOnProcPatches_[patchi].append(patchFacei);
605  }
606  }
607 }
608 
609 
610 void Foam::isoAdvection::applyBruteForceBounding()
611 {
612  addProfilingInFunction(geometricVoF);
613  bool alpha1Changed = false;
614 
615  const scalar snapAlphaTol = dict_.getOrDefault<scalar>("snapTol", 0);
616  if (snapAlphaTol > 0)
617  {
618  alpha1_ =
619  alpha1_
620  *pos0(alpha1_ - snapAlphaTol)
621  *neg0(alpha1_ - (1.0 - snapAlphaTol))
622  + pos0(alpha1_ - (1.0 - snapAlphaTol));
623 
624  alpha1Changed = true;
625  }
626 
627  if (dict_.getOrDefault("clip", true))
628  {
629  alpha1_ = min(scalar(1), max(scalar(0), alpha1_));
630  alpha1Changed = true;
631  }
632 
633  if (alpha1Changed)
634  {
635  alpha1_.correctBoundaryConditions();
636  }
637 }
638 
639 
641 {
642  if (!mesh_.time().writeTime()) return;
643 
644  if (dict_.getOrDefault("writeSurfCells", false))
645  {
646  cellSet cSet
647  (
648  IOobject
649  (
650  "surfCells",
651  mesh_.time().timeName(),
652  mesh_,
654  )
655  );
656 
657  cSet.insert(surfCells_);
658 
659  cSet.write();
660  }
661 }
662 
663 
665 (
666  const DynamicList<List<point>>& faces
667 ) const
668 {
669  if (!writeIsoFacesToFile_ || !mesh_.time().writeTime()) return;
670 
671  // Writing isofaces to obj file for inspection, e.g. in paraview
672  const fileName outputFile
673  (
674  mesh_.time().globalPath()
675  / "isoFaces"
676  / word::printf("isoFaces_%012d.obj", mesh_.time().timeIndex())
677  );
678 
679  if (Pstream::parRun())
680  {
681  // Collect points from all the processors
683  allProcFaces[Pstream::myProcNo()] = faces;
684  Pstream::gatherList(allProcFaces);
685 
686  if (Pstream::master())
687  {
688  mkDir(outputFile.path());
689  OBJstream os(outputFile);
690  Info<< nl << "isoAdvection: writing iso faces to file: "
691  << os.name() << nl << endl;
692 
693  face f;
694  forAll(allProcFaces, proci)
695  {
696  const DynamicList<List<point>>& procFacePts =
697  allProcFaces[proci];
698 
699  for (const List<point>& facePts : procFacePts)
700  {
701  if (facePts.size() != f.size())
702  {
703  f = face(identity(facePts.size()));
704  }
705 
706  os.write(f, facePts, false);
707  }
708  }
709  }
710  }
711  else
712  {
713  mkDir(outputFile.path());
714  OBJstream os(outputFile);
715  Info<< nl << "isoAdvection: writing iso faces to file: "
716  << os.name() << nl << endl;
717 
718  face f;
719  for (const List<point>& facePts : faces)
720  {
721  if (facePts.size() != f.size())
722  {
723  f = face(identity(facePts.size()));
724  }
725 
726  os.write(f, facePts, false);
727  }
728  }
729 }
730 
731 
732 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
volFields.H
profiling.H
Foam::objectRegistry::getObjectPtr
Type * getObjectPtr(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:423
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
meshTools.H
Foam::isoAdvection::writeSurfaceCells
void writeSurfaceCells() const
Return cellSet of surface cells.
Definition: isoAdvection.C:640
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::OBJstream
OFstream that keeps track of vertices.
Definition: OBJstream.H:58
Foam::OFstream::name
virtual const fileName & name() const
Read/write access to the name of the stream.
Definition: OSstream.H:107
Foam::fileName::path
static std::string path(const std::string &str)
Return directory path name (part before last /)
Definition: fileNameI.H:176
Foam::primitiveMesh::cellPoints
const labelListList & cellPoints() const
Definition: primitiveMeshCellPoints.C:34
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:55
Foam::GeometricField::primitiveField
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Definition: GeometricFieldI.H:53
Foam::neg0
dimensionedScalar neg0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:210
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
ref
rDeltaT ref()
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::OBJstream::write
virtual Ostream & write(const char c)
Write character.
Definition: OBJstream.C:78
syncTools.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::isoAdvection::writeIsoFaces
void writeIsoFaces(const DynamicList< List< point >> &isoFacePts) const
Write isoface points to .obj file.
Definition: isoAdvection.C:665
Foam::objectRegistry::foundObject
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
Definition: objectRegistryTemplates.C:379
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
interpolationCellPoint.H
Foam::polyBoundaryMesh::patchID
const labelList & patchID() const
Per boundary face label the patch index.
Definition: polyBoundaryMesh.C:456
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::syncTools::syncFaceList
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:396
Foam::regIOobject::write
virtual bool write(const bool valid=true) const
Write using setting from DB.
Definition: regIOobjectWrite.C:132
Foam::fvMesh::magSf
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Definition: fvMeshGeometry.C:330
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:388
upwind.H
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:187
Foam::cellList
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:47
Foam::primitiveMesh::cellCells
const labelListList & cellCells() const
Definition: primitiveMeshCellCells.C:102
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
fvcSurfaceIntegrate.H
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::FatalError
error FatalError
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
addProfilingInFunction
#define addProfilingInFunction(name)
Definition: profilingTrigger.H:124
Foam::PtrList::clear
void clear()
Clear the PtrList. Delete allocated entries and set size to zero.
Definition: PtrListI.H:97
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::primitiveMesh::cellVolumes
const scalarField & cellVolumes() const
Definition: primitiveMeshCellCentresAndVols.C:96
Foam::cellSet
A collection of cell labels.
Definition: cellSet.H:51
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
isoAdvection.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:487
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
volPointInterpolation.H
U
U
Definition: pEqn.H:72
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:52
Foam::UPstream::commsTypes::nonBlocking
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::nl
constexpr char nl
Definition: Ostream.H:404
f
labelList f(nPoints)
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::cutFace::debug
static int debug
Definition: cutFace.H:134
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
fvcGrad.H
Calculate the gradient of the given field.
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:77
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
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::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
bool
bool
Definition: EEqn.H:20
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
porosityProperties
IOdictionary porosityProperties(IOobject("porosityProperties", runTime.constant(), runTime, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE))
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::word::printf
static word printf(const char *fmt, const PrimitiveType &val)
Use a printf-style formatter for a primitive.
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::dimVol
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:64
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
cellSet.H
phib
phib
Definition: pEqn.H:189
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:507
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
boundary
faceListList boundary
Definition: createBlockMesh.H:4
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
OBJstream.H
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:445
Foam::primitiveMesh::faceAreas
const vectorField & faceAreas() const
Definition: primitiveMeshFaceCentresAndAreas.C:89