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 Johan Roenby
11  Modified code Copyright (C) 2019-2020 DLR
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 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46  defineTypeNameAndDebug(isoAdvection, 0);
47 }
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
51 Foam::isoAdvection::isoAdvection
52 (
54  const surfaceScalarField& phi,
55  const volVectorField& U
56 )
57 :
58  mesh_(alpha1.mesh()),
59  dict_(mesh_.solverDict(alpha1.name())),
60  alpha1_(alpha1),
61  alpha1In_(alpha1.ref()),
62  phi_(phi),
63  U_(U),
64  dVf_
65  (
66  IOobject
67  (
68  "dVf_",
69  mesh_.time().timeName(),
70  mesh_,
71  IOobject::NO_READ,
72  IOobject::NO_WRITE
73  ),
74  mesh_,
76  ),
77  alphaPhi_
78  (
79  IOobject
80  (
81  "alphaPhi_",
82  mesh_.time().timeName(),
83  mesh_,
84  IOobject::NO_READ,
85  IOobject::NO_WRITE
86  ),
87  mesh_,
89  ),
90  advectionTime_(0),
91 
92  // Tolerances and solution controls
93  nAlphaBounds_(dict_.getOrDefault<label>("nAlphaBounds", 3)),
94  isoFaceTol_(dict_.getOrDefault<scalar>("isoFaceTol", 1e-8)),
95  surfCellTol_(dict_.getOrDefault<scalar>("surfCellTol", 1e-8)),
96  writeIsoFacesToFile_(dict_.getOrDefault("writeIsoFaces", false)),
97 
98  // Cell cutting data
99  surfCells_(label(0.2*mesh_.nCells())),
100  surf_(reconstructionSchemes::New(alpha1_, phi_, U_, dict_)),
101  advectFace_(alpha1.mesh(), alpha1),
102  bsFaces_(label(0.2*mesh_.nBoundaryFaces())),
103  bsx0_(bsFaces_.size()),
104  bsn0_(bsFaces_.size()),
105  bsUn0_(bsFaces_.size()),
106 
107  // Parallel run data
108  procPatchLabels_(mesh_.boundary().size()),
109  surfaceCellFacesOnProcPatches_(0)
110 {
112 
113  // Prepare lists used in parallel runs
114  if (Pstream::parRun())
115  {
116  // Force calculation of required demand driven data (else parallel
117  // communication may crash)
118  mesh_.cellCentres();
119  mesh_.cellVolumes();
120  mesh_.faceCentres();
121  mesh_.faceAreas();
122  mesh_.magSf();
123  mesh_.boundaryMesh().patchID();
124  mesh_.cellPoints();
125  mesh_.cellCells();
126  mesh_.cells();
127 
128  // Get boundary mesh and resize the list for parallel comms
129  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
130 
131  surfaceCellFacesOnProcPatches_.resize(patches.size());
132 
133  // Append all processor patch labels to the list
134  forAll(patches, patchi)
135  {
136  if
137  (
138  isA<processorPolyPatch>(patches[patchi])
139  && patches[patchi].size() > 0
140  )
141  {
142  procPatchLabels_.append(patchi);
143  }
144  }
145  }
146 }
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
151 void Foam::isoAdvection::timeIntegratedFlux()
152 {
153  // Get time step
154  const scalar dt = mesh_.time().deltaTValue();
155 
156  // Create object for interpolating velocity to isoface centres
157  interpolationCellPoint<vector> UInterp(U_);
158 
159  // For each downwind face of each surface cell we "isoadvect" to find dVf
160  label nSurfaceCells = 0;
161 
162  // Clear out the data for re-use and reset list containing information
163  // whether cells could possibly need bounding
164  clearIsoFaceData();
165 
166  // Get necessary references
167  const scalarField& phiIn = phi_.primitiveField();
168  const scalarField& magSfIn = mesh_.magSf().primitiveField();
169  scalarField& dVfIn = dVf_.primitiveFieldRef();
170 
171  // Get necessary mesh data
172  const cellList& cellFaces = mesh_.cells();
173  const labelList& own = mesh_.faceOwner();
174 
175 
176  // Storage for isoFace points. Only used if writeIsoFacesToFile_
177  DynamicList<List<point>> isoFacePts;
178  const DynamicField<label>& interfaceLabels = surf_->interfaceLabels();
179 
180  // Loop through cells
181  forAll(interfaceLabels, i)
182  {
183  const label celli = interfaceLabels[i];
184  if (mag(surf_->normal()[celli]) != 0)
185  {
186 
187  // This is a surface cell, increment counter, append and mark cell
188  nSurfaceCells++;
189  surfCells_.append(celli);
190 
191  DebugInfo
192  << "\n------------ Cell " << celli << " with alpha1 = "
193  << alpha1In_[celli] << " and 1-alpha1 = "
194  << 1.0 - alpha1In_[celli] << " ------------"
195  << endl;
196 
197  // Cell is cut
198  const point x0 = surf_->centre()[celli];
199  vector n0 = -surf_->normal()[celli];
200  n0 /= (mag(n0));
201 
202  // Get the speed of the isoface by interpolating velocity and
203  // dotting it with isoface unit normal
204  const scalar Un0 = UInterp.interpolate(x0, celli) & n0;
205 
206  DebugInfo
207  << "calcIsoFace gives initial surface: \nx0 = " << x0
208  << ", \nn0 = " << n0 << ", \nUn0 = "
209  << Un0 << endl;
210 
211  // Estimate time integrated flux through each downwind face
212  // Note: looping over all cell faces - in reduced-D, some of
213  // these faces will be on empty patches
214  const cell& celliFaces = cellFaces[celli];
215  forAll(celliFaces, fi)
216  {
217  const label facei = celliFaces[fi];
218 
219  if (mesh_.isInternalFace(facei))
220  {
221  bool isDownwindFace = false;
222 
223  if (celli == own[facei])
224  {
225  if (phiIn[facei] >= 0)
226  {
227  isDownwindFace = true;
228  }
229  }
230  else
231  {
232  if (phiIn[facei] < 0)
233  {
234  isDownwindFace = true;
235  }
236  }
237 
238  if (isDownwindFace)
239  {
240  dVfIn[facei] = advectFace_.timeIntegratedFaceFlux
241  (
242  facei,
243  x0,
244  n0,
245  Un0,
246  dt,
247  phiIn[facei],
248  magSfIn[facei]
249  );
250  }
251 
252  }
253  else
254  {
255  bsFaces_.append(facei);
256  bsx0_.append(x0);
257  bsn0_.append(n0);
258  bsUn0_.append(Un0);
259 
260  // Note: we must not check if the face is on the
261  // processor patch here.
262  }
263  }
264  }
265  }
266 
267  // Get references to boundary fields
268  const polyBoundaryMesh& boundaryMesh = mesh_.boundaryMesh();
269  const surfaceScalarField::Boundary& phib = phi_.boundaryField();
270  const surfaceScalarField::Boundary& magSfb = mesh_.magSf().boundaryField();
271  surfaceScalarField::Boundary& dVfb = dVf_.boundaryFieldRef();
272  const label nInternalFaces = mesh_.nInternalFaces();
273 
274  // Loop through boundary surface faces
275  forAll(bsFaces_, i)
276  {
277  // Get boundary face index (in the global list)
278  const label facei = bsFaces_[i];
279  const label patchi = boundaryMesh.patchID()[facei - nInternalFaces];
280  const label start = boundaryMesh[patchi].start();
281 
282  if (phib[patchi].size())
283  {
284  const label patchFacei = facei - start;
285  const scalar phiP = phib[patchi][patchFacei];
286 
287  if (phiP >= 0)
288  {
289  const scalar magSf = magSfb[patchi][patchFacei];
290 
291  dVfb[patchi][patchFacei] = advectFace_.timeIntegratedFaceFlux
292  (
293  facei,
294  bsx0_[i],
295  bsn0_[i],
296  bsUn0_[i],
297  dt,
298  phiP,
299  magSf
300  );
301 
302  // Check if the face is on processor patch and append it to
303  // the list if necessary
304  checkIfOnProcPatch(facei);
305  }
306  }
307  }
308 
309  // Synchronize processor patches
310  syncProcPatches(dVf_, phi_);
311 
312  writeIsoFaces(isoFacePts);
313 
314  DebugInfo << "Number of isoAdvector surface cells = "
315  << returnReduce(nSurfaceCells, sumOp<label>()) << endl;
316 }
317 
318 
319 void Foam::isoAdvection::setDownwindFaces
320 (
321  const label celli,
322  DynamicLabelList& downwindFaces
323 ) const
324 {
326 
327  // Get necessary mesh data and cell information
328  const labelList& own = mesh_.faceOwner();
329  const cellList& cells = mesh_.cells();
330  const cell& c = cells[celli];
331 
332  downwindFaces.clear();
333 
334  // Check all faces of the cell
335  forAll(c, fi)
336  {
337  // Get face and corresponding flux
338  const label facei = c[fi];
339  const scalar phi = faceValue(phi_, facei);
340 
341  if (own[facei] == celli)
342  {
343  if (phi >= 0)
344  {
345  downwindFaces.append(facei);
346  }
347  }
348  else if (phi < 0)
349  {
350  downwindFaces.append(facei);
351  }
352  }
353 
354  downwindFaces.shrink();
355 }
356 
357 
358 Foam::scalar Foam::isoAdvection::netFlux
359 (
360  const surfaceScalarField& dVf,
361  const label celli
362 ) const
363 {
364  scalar dV = 0;
365 
366  // Get face indices
367  const cell& c = mesh_.cells()[celli];
368 
369  // Get mesh data
370  const labelList& own = mesh_.faceOwner();
371 
372  forAll(c, fi)
373  {
374  const label facei = c[fi];
375  const scalar dVff = faceValue(dVf, facei);
376 
377  if (own[facei] == celli)
378  {
379  dV += dVff;
380  }
381  else
382  {
383  dV -= dVff;
384  }
385  }
386 
387  return dV;
388 }
389 
390 
391 Foam::DynamicList<Foam::label> Foam::isoAdvection::syncProcPatches
392 (
393  surfaceScalarField& dVf,
394  const surfaceScalarField& phi,
395  bool returnSyncedFaces
396 )
397 {
398  DynamicLabelList syncedFaces(0);
399  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
400 
401  if (Pstream::parRun())
402  {
403  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
404 
405  // Send
406  forAll(procPatchLabels_, i)
407  {
408  const label patchi = procPatchLabels_[i];
409 
410  const processorPolyPatch& procPatch =
411  refCast<const processorPolyPatch>(patches[patchi]);
412 
413  UOPstream toNbr(procPatch.neighbProcNo(), pBufs);
414  const scalarField& pFlux = dVf.boundaryField()[patchi];
415 
416  const List<label>& surfCellFacesOnProcPatch =
417  surfaceCellFacesOnProcPatches_[patchi];
418 
419  const UIndirectList<scalar> dVfPatch
420  (
421  pFlux,
422  surfCellFacesOnProcPatch
423  );
424 
425  toNbr << surfCellFacesOnProcPatch << dVfPatch;
426  }
427 
428  pBufs.finishedSends();
429 
430 
431  // Receive and combine
432  forAll(procPatchLabels_, patchLabeli)
433  {
434  const label patchi = procPatchLabels_[patchLabeli];
435 
436  const processorPolyPatch& procPatch =
437  refCast<const processorPolyPatch>(patches[patchi]);
438 
439  UIPstream fromNeighb(procPatch.neighbProcNo(), pBufs);
440  List<label> faceIDs;
441  List<scalar> nbrdVfs;
442 
443  fromNeighb >> faceIDs >> nbrdVfs;
444  if (returnSyncedFaces)
445  {
446  List<label> syncedFaceI(faceIDs);
447  for (label& faceI : syncedFaceI)
448  {
449  faceI += procPatch.start();
450  }
451  syncedFaces.append(syncedFaceI);
452  }
453 
454  if (debug)
455  {
456  Pout<< "Received at time = " << mesh_.time().value()
457  << ": surfCellFacesOnProcPatch = " << faceIDs << nl
458  << "Received at time = " << mesh_.time().value()
459  << ": dVfPatch = " << nbrdVfs << endl;
460  }
461 
462  // Combine fluxes
463  scalarField& localFlux = dVf.boundaryFieldRef()[patchi];
464 
465  forAll(faceIDs, i)
466  {
467  const label facei = faceIDs[i];
468  localFlux[facei] = - nbrdVfs[i];
469  if (debug && mag(localFlux[facei] + nbrdVfs[i]) > 10*SMALL)
470  {
471  Pout<< "localFlux[facei] = " << localFlux[facei]
472  << " and nbrdVfs[i] = " << nbrdVfs[i]
473  << " for facei = " << facei << endl;
474  }
475  }
476  }
477 
478  if (debug)
479  {
480  // Write out results for checking
481  forAll(procPatchLabels_, patchLabeli)
482  {
483  const label patchi = procPatchLabels_[patchLabeli];
484  const scalarField& localFlux = dVf.boundaryField()[patchi];
485  Pout<< "time = " << mesh_.time().value() << ": localFlux = "
486  << localFlux << endl;
487  }
488  }
489 
490  // Reinitialising list used for minimal parallel communication
491  forAll(surfaceCellFacesOnProcPatches_, patchi)
492  {
493  surfaceCellFacesOnProcPatches_[patchi].clear();
494  }
495  }
496 
497  return syncedFaces;
498 }
499 
500 
501 void Foam::isoAdvection::checkIfOnProcPatch(const label facei)
502 {
503  if (!mesh_.isInternalFace(facei))
504  {
505  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
506  const label patchi = pbm.patchID()[facei - mesh_.nInternalFaces()];
507 
508  if (isA<processorPolyPatch>(pbm[patchi]) && pbm[patchi].size())
509  {
510  const label patchFacei = pbm[patchi].whichFace(facei);
511  surfaceCellFacesOnProcPatches_[patchi].append(patchFacei);
512  }
513  }
514 }
515 
516 
517 
518 
520 {
521  bool alpha1Changed = false;
522 
523  const scalar snapAlphaTol = dict_.getOrDefault<scalar>("snapTol", 0);
524  if (snapAlphaTol > 0)
525  {
526  alpha1_ =
527  alpha1_
528  *pos0(alpha1_ - snapAlphaTol)
529  *neg0(alpha1_ - (1.0 - snapAlphaTol))
530  + pos0(alpha1_ - (1.0 - snapAlphaTol));
531 
532  alpha1Changed = true;
533  }
534 
535  if (dict_.getOrDefault("clip", true))
536  {
537  alpha1_ = min(scalar(1), max(scalar(0), alpha1_));
538  alpha1Changed = true;
539  }
540 
541  if (alpha1Changed)
542  {
543  alpha1_.correctBoundaryConditions();
544  }
545 }
546 
547 
549 {
550  if (!mesh_.time().writeTime()) return;
551 
552  if (dict_.getOrDefault("writeSurfCells", false))
553  {
554  cellSet cSet
555  (
556  IOobject
557  (
558  "surfCells",
559  mesh_.time().timeName(),
560  mesh_,
562  )
563  );
564 
565  cSet.insert(surfCells_);
566 
567  cSet.write();
568  }
569 }
570 
572 (
573  const DynamicList<List<point>>& faces
574 ) const
575 {
576  if (!writeIsoFacesToFile_ || !mesh_.time().writeTime()) return;
577 
578  // Writing isofaces to obj file for inspection, e.g. in paraview
579  const fileName outputFile
580  (
581  mesh_.time().globalPath()
582  / "isoFaces"
583  / word::printf("isoFaces_%012d.obj", mesh_.time().timeIndex())
584  );
585 
586  if (Pstream::parRun())
587  {
588  // Collect points from all the processors
590  allProcFaces[Pstream::myProcNo()] = faces;
591  Pstream::gatherList(allProcFaces);
592 
593  if (Pstream::master())
594  {
595  mkDir(outputFile.path());
596  OBJstream os(outputFile);
597  Info<< nl << "isoAdvection: writing iso faces to file: "
598  << os.name() << nl << endl;
599 
600  face f;
601  forAll(allProcFaces, proci)
602  {
603  const DynamicList<List<point>>& procFacePts =
604  allProcFaces[proci];
605 
606  forAll(procFacePts, i)
607  {
608  const List<point>& facePts = procFacePts[i];
609 
610  if (facePts.size() != f.size())
611  {
612  f = face(identity(facePts.size()));
613  }
614 
615  os.write(f, facePts, false);
616  }
617  }
618  }
619  }
620  else
621  {
622  mkDir(outputFile.path());
623  OBJstream os(outputFile);
624  Info<< nl << "isoAdvection: writing iso faces to file: "
625  << os.name() << nl << endl;
626 
627  face f;
628  forAll(faces, i)
629  {
630  const List<point>& facePts = faces[i];
631 
632  if (facePts.size() != f.size())
633  {
634  f = face(identity(facePts.size()));
635  }
636 
637  os.write(f, facePts, false);
638  }
639  }
640 }
641 
642 
643 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
volFields.H
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:104
meshTools.H
Foam::UPstream::commsTypes::nonBlocking
Foam::isoAdvection::writeSurfaceCells
void writeSurfaceCells() const
Return cellSet of surface cells.
Definition: isoAdvection.C:548
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
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::fileName::path
static std::string path(const std::string &str)
Return directory path name (part before last /)
Definition: fileNameI.H:186
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::UPstream::parRun
static bool & parRun()
Test if this a parallel run, or allow modify access.
Definition: UPstream.H:434
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:458
Foam::isoAdvection::applyBruteForceBounding
void applyBruteForceBounding()
Apply the bounding based on user inputs.
Definition: isoAdvection.C:519
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:350
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::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:572
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:452
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::regIOobject::write
virtual bool write(const bool valid=true) const
Write using setting from DB.
Definition: regIOobjectWrite.C:165
Foam::fvMesh::magSf
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Definition: fvMeshGeometry.C:330
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:365
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:474
upwind.H
Foam::Vector::centre
const Vector< Cmpt > & centre(const Foam::List< Vector< Cmpt >> &) const
Return *this (used for point which is a typedef to Vector<scalar>.
Definition: VectorI.H:114
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:43
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
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
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::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
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
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:359
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:464
Foam::nl
constexpr char nl
Definition: Ostream.H:385
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::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:115
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:181
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::word::printf
static word printf(const char *fmt, const PrimitiveType &val)
Use a printf-style formatter for a primitive.
Foam::dimVol
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:65
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::point
vector point
Point is a vector.
Definition: point.H:43
Foam::IOobject::NO_READ
Definition: IOobject.H:123
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
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:446
Foam::primitiveMesh::faceAreas
const vectorField & faceAreas() const
Definition: primitiveMeshFaceCentresAndAreas.C:89