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-2022 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.primitiveFieldRef()),
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  porosityPtr_ = mesh_.getObjectPtr<volScalarField>("porosity");
160 
161  if (porosityPtr_)
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 
372  // Loop through boundary surface faces
373  forAll(bsFaces_, i)
374  {
375  // Get boundary face index (in the global list)
376  const label facei = bsFaces_[i];
377  const label patchi = boundaryMesh.patchID(facei);
378 
379  if (!phib[patchi].empty())
380  {
381  const label patchFacei = boundaryMesh[patchi].whichFace(facei);
382 
383  const scalar phiP = phib[patchi][patchFacei];
384 
385  if (phiP >= 0)
386  {
387  const scalar magSf = magSfb[patchi][patchFacei];
388 
389  dVfb[patchi][patchFacei] = advectFace_.timeIntegratedFaceFlux
390  (
391  facei,
392  bsx0_[i],
393  bsn0_[i],
394  bsUn0_[i],
395  dt,
396  phiP,
397  magSf
398  );
399 
400  // Check if the face is on processor patch and append it to
401  // the list if necessary
402  checkIfOnProcPatch(facei);
403  }
404  }
405  }
406 
407  // Synchronize processor patches
408  syncProcPatches(dVf_, phi_);
409 
410  writeIsoFaces(isoFacePts);
411 
412  DebugInfo << "Number of isoAdvector surface cells = "
413  << returnReduce(nSurfaceCells, sumOp<label>()) << endl;
414 }
415 
416 
417 void Foam::isoAdvection::setDownwindFaces
418 (
419  const label celli,
420  DynamicLabelList& downwindFaces
421 ) const
422 {
424 
425  // Get necessary mesh data and cell information
426  const labelList& own = mesh_.faceOwner();
427  const cellList& cells = mesh_.cells();
428  const cell& c = cells[celli];
429 
430  downwindFaces.clear();
431 
432  // Check all faces of the cell
433  for (const label facei: c)
434  {
435  // Get face and corresponding flux
436  const scalar phi = faceValue(phi_, facei);
437 
438  if (own[facei] == celli)
439  {
440  if (phi >= 0)
441  {
442  downwindFaces.append(facei);
443  }
444  }
445  else if (phi < 0)
446  {
447  downwindFaces.append(facei);
448  }
449  }
450 
451  downwindFaces.shrink();
452 }
453 
454 
455 Foam::scalar Foam::isoAdvection::netFlux
456 (
457  const surfaceScalarField& dVf,
458  const label celli
459 ) const
460 {
461  scalar dV = 0;
462 
463  // Get face indices
464  const cell& c = mesh_.cells()[celli];
465 
466  // Get mesh data
467  const labelList& own = mesh_.faceOwner();
468 
469  for (const label facei : c)
470  {
471  const scalar dVff = faceValue(dVf, facei);
472 
473  if (own[facei] == celli)
474  {
475  dV += dVff;
476  }
477  else
478  {
479  dV -= dVff;
480  }
481  }
482 
483  return dV;
484 }
485 
486 
487 Foam::DynamicList<Foam::label> Foam::isoAdvection::syncProcPatches
488 (
489  surfaceScalarField& dVf,
490  const surfaceScalarField& phi,
491  bool returnSyncedFaces
492 )
493 {
494  DynamicLabelList syncedFaces(0);
495  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
496 
497  if (Pstream::parRun())
498  {
499  DynamicList<label> neighProcs;
500  PstreamBuffers pBufs(UPstream::commsTypes::nonBlocking);
501 
502  // Send
503  for (const label patchi : procPatchLabels_)
504  {
505  const processorPolyPatch& procPatch =
506  refCast<const processorPolyPatch>(patches[patchi]);
507  const label nbrProci = procPatch.neighbProcNo();
508 
509  // Neighbour connectivity
510  neighProcs.push_uniq(nbrProci);
511 
512  const scalarField& pFlux = dVf.boundaryField()[patchi];
513  const List<label>& surfCellFacesOnProcPatch =
514  surfaceCellFacesOnProcPatches_[patchi];
515 
516  const UIndirectList<scalar> dVfPatch
517  (
518  pFlux,
519  surfCellFacesOnProcPatch
520  );
521 
522  UOPstream toNbr(nbrProci, pBufs);
523  toNbr << surfCellFacesOnProcPatch << dVfPatch;
524  }
525 
526  // Limited to involved neighbour procs
527  pBufs.finishedNeighbourSends(neighProcs);
528 
529 
530  // Receive and combine
531  for (const label patchi : procPatchLabels_)
532  {
533  const processorPolyPatch& procPatch =
534  refCast<const processorPolyPatch>(patches[patchi]);
535  const label nbrProci = procPatch.neighbProcNo();
536 
537  List<label> faceIDs;
538  List<scalar> nbrdVfs;
539 
540  {
541  UIPstream fromNbr(nbrProci, pBufs);
542  fromNbr >> faceIDs >> nbrdVfs;
543  }
544 
545  if (returnSyncedFaces)
546  {
547  List<label> syncedFaceI(faceIDs);
548  for (label& faceI : syncedFaceI)
549  {
550  faceI += procPatch.start();
551  }
552  syncedFaces.append(syncedFaceI);
553  }
554 
555  if (debug)
556  {
557  Pout<< "Received at time = " << mesh_.time().value()
558  << ": surfCellFacesOnProcPatch = " << faceIDs << nl
559  << "Received at time = " << mesh_.time().value()
560  << ": dVfPatch = " << nbrdVfs << endl;
561  }
562 
563  // Combine fluxes
564  scalarField& localFlux = dVf.boundaryFieldRef()[patchi];
565 
566  forAll(faceIDs, i)
567  {
568  const label facei = faceIDs[i];
569  localFlux[facei] = - nbrdVfs[i];
570  if (debug && mag(localFlux[facei] + nbrdVfs[i]) > ROOTVSMALL)
571  {
572  Pout<< "localFlux[facei] = " << localFlux[facei]
573  << " and nbrdVfs[i] = " << nbrdVfs[i]
574  << " for facei = " << facei << endl;
575  }
576  }
577  }
578 
579  if (debug)
580  {
581  // Write out results for checking
582  forAll(procPatchLabels_, patchLabeli)
583  {
584  const label patchi = procPatchLabels_[patchLabeli];
585  const scalarField& localFlux = dVf.boundaryField()[patchi];
586  Pout<< "time = " << mesh_.time().value() << ": localFlux = "
587  << localFlux << endl;
588  }
589  }
590 
591  // Reinitialising list used for minimal parallel communication
592  forAll(surfaceCellFacesOnProcPatches_, patchi)
593  {
594  surfaceCellFacesOnProcPatches_[patchi].clear();
595  }
596  }
597 
598  return syncedFaces;
599 }
600 
601 
602 void Foam::isoAdvection::checkIfOnProcPatch(const label facei)
603 {
604  if (!mesh_.isInternalFace(facei))
605  {
606  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
607  const label patchi = pbm.patchID(facei);
608 
609  if (isA<processorPolyPatch>(pbm[patchi]) && !pbm[patchi].empty())
610  {
611  const label patchFacei = pbm[patchi].whichFace(facei);
612  surfaceCellFacesOnProcPatches_[patchi].append(patchFacei);
613  }
614  }
615 }
616 
617 
618 void Foam::isoAdvection::applyBruteForceBounding()
619 {
620  addProfilingInFunction(geometricVoF);
621  bool alpha1Changed = false;
622 
623  const scalar snapAlphaTol = dict_.getOrDefault<scalar>("snapTol", 0);
624  if (snapAlphaTol > 0)
625  {
626  alpha1_ =
627  alpha1_
628  *pos0(alpha1_ - snapAlphaTol)
629  *neg0(alpha1_ - (1.0 - snapAlphaTol))
630  + pos0(alpha1_ - (1.0 - snapAlphaTol));
631 
632  alpha1Changed = true;
633  }
634 
635  if (dict_.getOrDefault("clip", true))
636  {
637  alpha1_.clamp_range(zero_one{});
638  alpha1Changed = true;
639  }
640 
641  if (alpha1Changed)
642  {
643  alpha1_.correctBoundaryConditions();
644  }
645 }
646 
647 
649 {
650  if (!mesh_.time().writeTime()) return;
651 
652  if (dict_.getOrDefault("writeSurfCells", false))
653  {
654  cellSet cSet
655  (
656  IOobject
657  (
658  "surfCells",
659  mesh_.time().timeName(),
660  mesh_,
662  )
663  );
664 
665  cSet.insert(surfCells_);
667  cSet.write();
668  }
669 }
670 
671 
673 (
674  const DynamicList<List<point>>& faces
675 ) const
676 {
677  if (!writeIsoFacesToFile_ || !mesh_.time().writeTime()) return;
678 
679  // Writing isofaces to obj file for inspection, e.g. in paraview
680  const fileName outputFile
681  (
682  mesh_.time().globalPath()
683  / "isoFaces"
684  / word::printf("isoFaces_%012d.obj", mesh_.time().timeIndex())
685  );
686 
687  if (Pstream::parRun())
688  {
689  // Collect points from all the processors
690  List<DynamicList<List<point>>> allProcFaces(Pstream::nProcs());
691  allProcFaces[Pstream::myProcNo()] = faces;
692  Pstream::gatherList(allProcFaces);
693 
694  if (Pstream::master())
695  {
696  mkDir(outputFile.path());
697  OBJstream os(outputFile);
698  Info<< nl << "isoAdvection: writing iso faces to file: "
699  << os.name() << nl << endl;
700 
701  for
702  (
703  const DynamicList<List<point>>& procFacePts
704  : allProcFaces
705  )
706  {
707  for (const List<point>& facePts : procFacePts)
708  {
709  os.writeFace(facePts, false);
710  }
711  }
712  }
713  }
714  else
715  {
716  mkDir(outputFile.path());
717  OBJstream os(outputFile);
718  Info<< nl << "isoAdvection: writing iso faces to file: "
719  << os.name() << nl << endl;
720 
721  for (const List<point>& facePts : faces)
722  {
723  os.writeFace(facePts, false);
724  }
725  }
726 }
727 
728 
729 // ************************************************************************* //
faceListList boundary
Ostream & writeFace(const UList< point > &points, const bool lines=true)
Write face loop points with lines/filled-polygon.
Definition: OBJstream.C:279
const polyBoundaryMesh & pbm
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
List< cell > cellList
List of cell.
Definition: cellListFwd.H:39
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Type gMin(const FieldField< Field, Type > &f)
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition: OSstream.H:128
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const labelList & patchID() const
Per boundary face label the patch index.
void writeIsoFaces(const DynamicList< List< point >> &isoFacePts) const
Write isoface points to .obj file.
Definition: isoAdvection.C:666
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:65
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:432
const cellList & cells() const
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
Type of boundary fields.
Ignore writing from objectRegistry::writeObject()
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
conserve primitiveFieldRef()+
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void writeSurfaceCells() const
Return cellSet of surface cells.
Definition: isoAdvection.C:641
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
word timeName
Definition: getTimeIndex.H:3
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1065
static void gatherList(const List< commsStruct > &comms, List< T > &values, const int tag, const label comm)
Gather data, but keep individual values separate. Uses the specified communication schedule...
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:614
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
Calculate the gradient of the given field.
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:137
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
#define DebugInFunction
Report an information message using Foam::Info.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
dimensionedScalar neg0(const dimensionedScalar &ds)
Reading is optional [identical to LAZY_READ].
Vector< scalar > vector
Definition: vector.H:57
const vectorField & cellCentres() const
static int debug
Definition: cutFace.H:143
#define DebugInfo
Report an information message using Foam::Info.
static word printf(const char *fmt, const PrimitiveType &val)
Use a printf-style formatter for a primitive.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
dimensionedScalar pos0(const dimensionedScalar &ds)
void clear()
Clear the patch list and all demand-driven data.
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:112
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
IOdictionary porosityProperties(IOobject("porosityProperties", runTime.constant(), runTime, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE))
const vectorField & faceCentres() const
void append(autoPtr< T > &ptr)
Move append an element to the end of the list.
Definition: PtrList.H:344
U
Definition: pEqn.H:72
vector point
Point is a vector.
Definition: point.H:37
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const vectorField & faceAreas() const
const dimensionedScalar c
Speed of light in a vacuum.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
const polyBoundaryMesh & patches
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
"nonBlocking" : (MPI_Isend, MPI_Irecv)
#define addProfilingInFunction(Name)
Define profiling trigger with specified name and description corresponding to the compiler-defined fu...
messageStream Info
Information stream (stdout output on master, null elsewhere)
phib
Definition: pEqn.H:189
List< label > labelList
A List of labels.
Definition: List.H:62
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const labelListList & cellCells() const
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
const labelListList & cellPoints() const
Namespace for OpenFOAM.
const scalarField & cellVolumes() const
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const volScalarField & alpha1