Loading...
Searching...
No Matches
extractEulerianParticles.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-2023 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29#include "regionSplit2D.H"
31#include "volFields.H"
32#include "surfaceFields.H"
33#include "surfaceInterpolate.H"
35#include "emptyPolyPatch.H"
36#include "coupledPolyPatch.H"
37#include "binned.H"
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
42namespace Foam
43{
44namespace functionObjects
45{
48 (
52 );
53}
54}
55
56// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
57
59{
61
63 if (zoneID_ == -1)
64 {
66 << "Unable to find faceZone " << faceZoneName_
67 << ". Available faceZones are: " << mesh_.faceZones().names()
68 << exit(FatalError);
69 }
70
71 const faceZone& fz = mesh_.faceZones()[zoneID_];
72 const label nFaces = fz.size();
73 const label allFaces = returnReduce(nFaces, sumOp<label>());
74
75 if (allFaces < nInjectorLocations_)
76 {
78 << "faceZone " << faceZoneName_
79 << ": Number of faceZone faces (" << allFaces
80 << ") is less than the number of requested locations ("
81 << nInjectorLocations_ << ")."
82 << exit(FatalError);
83 }
84
85 Info<< type() << " " << name() << " output:" << nl
86 << " faceZone : " << faceZoneName_ << nl
87 << " faces : " << allFaces << nl
88 << endl;
90 // Initialise old iteration blocked faces
91 // Note: for restart, this info needs to be written/read
92 regions0_.setSize(fz.size(), -1);
93}
94
95
97{
99
100 if (!nInjectorLocations_)
101 {
102 return;
103 }
104
105 const faceZone& fz = mesh_.faceZones()[zoneID_];
106
107 // Agglomerate faceZone faces into nInjectorLocations_ global locations
108 const indirectPrimitivePatch patch
109 (
110 IndirectList<face>(mesh_.faces(), fz),
111 mesh_.points()
112 );
113
114 const label nFaces = fz.size();
115 label nLocations = nInjectorLocations_;
116
117 if (Pstream::parRun())
118 {
119 label nGlobalFaces = returnReduce(nFaces, sumOp<label>());
120 scalar fraction = scalar(nFaces)/scalar(nGlobalFaces);
121 nLocations = ceil(fraction*nInjectorLocations_);
122 if (debug)
123 {
124 Pout<< "nFaces:" << nFaces
125 << ", nGlobalFaces:" << nGlobalFaces
126 << ", fraction:" << fraction
127 << ", nLocations:" << nLocations
128 << endl;
129 }
130 }
131
132 pairPatchAgglomeration ppa
133 (
134 patch.localFaces(),
135 patch.localPoints(),
136 10,
137 50,
138 nLocations,
139 labelMax,
140 180
141 );
142
143 ppa.agglomerate();
144
145 label nCoarseFaces = 0;
146 if (nFaces != 0)
147 {
148 fineToCoarseAddr_ = ppa.restrictTopBottomAddressing();
149 nCoarseFaces = max(fineToCoarseAddr_) + 1;
150 }
151
152 globalCoarseFaces_ = globalIndex(nCoarseFaces);
154 Info<< "Created " << returnReduce(nCoarseFaces, sumOp<label>())
155 << " coarse faces" << endl;
156}
157
158
161{
163
164 const auto& phi = mesh_.lookupObject<surfaceScalarField>(phiName_);
165
166 if (phi.dimensions() == dimMass/dimTime)
167 {
168 const auto& rho = mesh_.lookupObject<volScalarField>(rhoName_);
169
171 }
172
173 return phi;
174}
175
176
178(
179 const surfaceScalarField& alphaf,
180 const faceZone& fz,
181 boolList& blockedFaces
182)
183{
185
186 // Initialise storage for patch and patch-face indices where faceZone
187 // intersects mesh patch(es)
188 patchIDs_.setSize(fz.size(), -1);
189 patchFaceIDs_.setSize(fz.size(), -1);
190
191 label nBlockedFaces = 0;
192 forAll(fz, localFacei)
193 {
194 const label facei = fz[localFacei];
195
196 if (mesh_.isInternalFace(facei))
197 {
198 if (alphaf[facei] > alphaThreshold_)
199 {
200 blockedFaces[localFacei] = true;
201 }
202 }
203 else
204 {
205 label patchi = mesh_.boundaryMesh().whichPatch(facei);
206 label patchFacei = -1;
207
208 const polyPatch& pp = mesh_.boundaryMesh()[patchi];
209 const scalarField& alphafp = alphaf.boundaryField()[patchi];
210 const auto* cpp = isA<coupledPolyPatch>(pp);
211
212 if (cpp)
213 {
214 patchFacei = (cpp->owner() ? pp.whichFace(facei) : -1);
215 }
216 else if (!isA<emptyPolyPatch>(pp))
217 {
218 patchFacei = pp.whichFace(facei);
219 }
220
221 if (patchFacei == -1)
222 {
223 patchi = -1;
224 }
225 else if (alphafp[patchFacei] > alphaThreshold_)
226 {
227 blockedFaces[localFacei] = true;
228 }
229
230 patchIDs_[localFacei] = patchi;
231 patchFaceIDs_[localFacei] = patchFacei;
233 }
234
235 DebugInFunction << "Number of blocked faces: " << nBlockedFaces << endl;
236}
237
238
240(
241 const scalar time,
242 const label regioni
243)
244{
245 DebugInFunction << "collectParticle: " << regioni << endl;
246
247 const label particlei = regionToParticleMap_[regioni];
248 eulerianParticle p = particles_[particlei];
249
250 if (p.faceIHit != -1 && nInjectorLocations_)
251 {
252 // Use coarse face index for tag output
253 label coarseFacei = fineToCoarseAddr_[p.faceIHit];
254 p.faceIHit = globalCoarseFaces_.toGlobal(coarseFacei);
255 }
256
258
259 const scalar pDiameter = cbrt(6.0*p.V/constant::mathematical::pi);
260
261 if ((pDiameter > minDiameter_) && (pDiameter < maxDiameter_))
262 {
263 if (Pstream::master())
264 {
265 const scalar d = cbrt(6.0*p.V/constant::mathematical::pi);
266 const point position = p.VC/(p.V + ROOTVSMALL);
267 const vector U = p.VU/(p.V + ROOTVSMALL);
268 label tag = -1;
269 if (nInjectorLocations_)
270 {
271 tag = p.faceIHit;
272 }
273
274 injectedParticle* ip = new injectedParticle
275 (
276 mesh_,
277 position,
278 tag,
279 time,
280 d,
281 U,
282 false // not looking to set cell owner etc.
283 );
284
285 cloud_.addParticle(ip);
286
287 collectedVolume_ += p.V;
288 }
289
290 ++nCollectedParticles_;
291 }
292 else
293 {
294 // Discard particles over/under diameter threshold
296 discardedVolume_ += p.V;
297 }
298}
299
300
302(
303 const label nNewRegions,
304 const scalar time,
305 labelList& regionFaceIDs
306)
307{
309
310 // Determine mapping between old and new regions so that we can
311 // accumulate particle info
312 labelList oldToNewRegion(particles_.size(), -1);
313 labelList newToNewRegion(identity(nNewRegions));
314
315 forAll(regionFaceIDs, facei)
316 {
317 label newRegioni = regionFaceIDs[facei];
318 label oldRegioni = regions0_[facei];
319
320 if (newRegioni != -1 && oldRegioni != -1)
321 {
322 // If old region has split into multiple regions we need to
323 // renumber new regions to maintain connectivity with old regions
324 newToNewRegion[newRegioni] =
325 max(newRegioni, oldToNewRegion[oldRegioni]);
326 oldToNewRegion[oldRegioni] = newRegioni;
327 }
328 }
329
330 // Create map from new regions to slots in particles list
331 // - filter through new-to-new addressing to identify new particles
332 Pstream::listReduce(newToNewRegion, maxOp<label>());
333
334 label nParticle = -1;
335 labelHashSet newRegions;
336 Map<label> newRegionToParticleMap;
337 forAll(newToNewRegion, newRegioni0)
338 {
339 label newRegioni = newToNewRegion[newRegioni0];
340 if (newRegions.insert(newRegioni))
341 {
342 ++nParticle;
343 }
344
345 // New particle slot
346 newRegionToParticleMap.insert(newRegioni0, nParticle);
347 }
348
349 // Accumulate old region data or create a new particle if there is no
350 // mapping from the old-to-new region
351 Pstream::listReduce(oldToNewRegion, maxOp<label>());
352
353 List<eulerianParticle> newParticles(newRegionToParticleMap.size());
354 forAll(oldToNewRegion, oldRegioni)
355 {
356 label newRegioni = oldToNewRegion[oldRegioni];
357 if (newRegioni == -1)
358 {
359 // No mapping from old-to-new - collect new particle
361 << "Collecting particle from oldRegion:" << oldRegioni
362 << endl;
363
364 collectParticle(time, oldRegioni);
365 }
366 else
367 {
368 // Combine existing particle into new particle
369 label newParticlei = newRegionToParticleMap[newRegioni];
370 label oldParticlei = regionToParticleMap_[oldRegioni];
371
373 << "Combining newRegioni: " << newRegioni
374 << "(p:" << newParticlei << ") and "
375 << "oldRegioni: " << oldRegioni
376 << "(p:" << oldParticlei << ")"
377 << endl;
378
379 newParticles[newParticlei] =
381 (
382 newParticles[newParticlei],
383 particles_[oldParticlei]
384 );
385 }
386 }
387
388 // Reset the particles list and addressing for latest available info
389 particles_.transfer(newParticles);
390 regionToParticleMap_ = newRegionToParticleMap;
392 // Reset the region IDs for the next integration step
393 // - these become the oldRegioni's
394 regions0_ = regionFaceIDs;
395}
396
397
399(
400 const surfaceScalarField& alphaf,
401 const surfaceScalarField& phi,
402 const labelList& regionFaceIDs,
403 const faceZone& fz
404)
405{
407
408 const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
410
411 const scalar deltaT = mesh_.time().deltaTValue();
412 const pointField& faceCentres = mesh_.faceCentres();
413
414 forAll(regionFaceIDs, localFacei)
415 {
416 const label newRegioni = regionFaceIDs[localFacei];
417 if (newRegioni != -1)
418 {
419 const label particlei = regionToParticleMap_[newRegioni];
420 const label meshFacei = fz[localFacei];
421 eulerianParticle& p = particles_[particlei];
422
423 if (p.faceIHit < 0)
424 {
425 // New particle - does not exist in particles_ list
426 p.faceIHit = localFacei;
427 p.V = 0;
428 p.VC = vector::zero;
429 p.VU = vector::zero;
430 }
431
432 // Accumulate particle properties
433 scalar magPhii = mag(faceValue(phi, localFacei, meshFacei));
434 vector Ufi = faceValue(Uf, localFacei, meshFacei);
435 scalar dV = magPhii*deltaT;
436 p.V += dV;
437 p.VC += dV*faceCentres[meshFacei];
438 p.VU += dV*Ufi;
440 }
441}
442
443
444// * * * * * * * * * * * * * * * * Constructor * * * * * * * * * * * * * * * //
445
447(
448 const word& name,
449 const Time& runTime,
450 const dictionary& dict
451)
452:
453 fvMeshFunctionObject(name, runTime, dict),
454 writeFile(runTime, name),
455 cloud_(mesh_, "eulerianParticleCloud"),
456 faceZoneName_(),
457 zoneID_(-1),
458 patchIDs_(),
459 patchFaceIDs_(),
460 alphaName_("alpha"),
461 alphaThreshold_(0.1),
462 UName_("U"),
463 rhoName_("rho"),
464 phiName_("phi"),
465 nInjectorLocations_(0),
466 fineToCoarseAddr_(),
467 globalCoarseFaces_(),
468 regions0_(),
469 particles_(),
470 regionToParticleMap_(),
471 minDiameter_(ROOTVSMALL),
472 maxDiameter_(GREAT),
473 nCollectedParticles_(getProperty<label>("nCollectedParticles", 0)),
474 collectedVolume_(getProperty<scalar>("collectedVolume", 0)),
475 nDiscardedParticles_(getProperty<label>("nDiscardedParticles", 0)),
476 discardedVolume_(getProperty<scalar>("discardedVolume", 0))
477{
478 if (mesh_.nSolutionD() != 3)
479 {
480 FatalErrorInFunction
481 << name << " function object only applicable to 3-D cases"
482 << exit(FatalError);
483 }
485 read(dict);
486}
487
488
489// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
490
492(
493 const dictionary& dict
494)
495{
497
499 {
500 dict.readEntry("faceZone", faceZoneName_);
501 dict.readEntry("alpha", alphaName_);
502
503 dict.readIfPresent("alphaThreshold", alphaThreshold_);
504 dict.readIfPresent("U", UName_);
505 dict.readIfPresent("rho", rhoName_);
506 dict.readIfPresent("phi", phiName_);
507 dict.readIfPresent("nLocations", nInjectorLocations_);
508 dict.readIfPresent("minDiameter", minDiameter_);
509 dict.readIfPresent("maxDiameter", maxDiameter_);
510
511 checkFaceZone();
512
513 if (nInjectorLocations_)
514 {
515 initialiseBins();
516 }
517
518 return true;
519 }
520
521 return false;
522}
523
524
526{
528
529 Log << type() << " " << name() << " output:" << nl;
530
531 const volScalarField& alpha =
532 mesh_.lookupObject<volScalarField>(alphaName_);
533
534 auto talphaf = surfaceScalarField::New
535 (
539 );
540
541 const faceZone& fz = mesh_.faceZones()[zoneID_];
542 const indirectPrimitivePatch patch
543 (
544 IndirectList<face>(mesh_.faces(), fz),
545 mesh_.points()
546 );
547
548 // Set the blocked faces, i.e. where alpha > alpha threshold value
549 boolList blockedFaces(fz.size(), false);
550 setBlockedFaces(talphaf(), fz, blockedFaces);
551
552 // Split the faceZone according to the blockedFaces
553 // - Returns a list of (disconnected) region index per face zone face
554 regionSplit2D regionFaceIDs(mesh_, patch, blockedFaces);
555
556 // Global number of regions
557 const label nRegionsNew = regionFaceIDs.nRegions();
558
559 // Calculate the addressing between the old and new region information
560 // Also collects particles that have traversed the faceZone
561 // - Note: may also update regionFaceIDs
562 calculateAddressing
563 (
564 nRegionsNew,
565 mesh_.time().value(),
566 regionFaceIDs
567 );
568
569 // Process latest region information
570 tmp<surfaceScalarField> tphi = phiU();
571 accumulateParticleInfo(talphaf(), tphi(), regionFaceIDs, fz);
572
573 Log << " Collected particles : " << nCollectedParticles_ << nl
574 << " Collected volume : " << collectedVolume_ << nl
575 << " Discarded particles : " << nDiscardedParticles_ << nl
576 << " Discarded volume : " << discardedVolume_ << nl
577 << " Particles in progress : " << particles_.size() << nl
578 << endl;
579
580 return true;
581}
582
583
585{
587
588 cloud_.write();
589
590 setProperty("nCollectedParticles", nCollectedParticles_);
591 setProperty("collectedVolume", collectedVolume_);
592 setProperty("nDiscardedParticles", nDiscardedParticles_);
593 setProperty("discardedVolume", discardedVolume_);
594
595 return true;
596}
597
598
599// ************************************************************************* //
#define Log
Definition PDRblock.C:28
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvsPatchField< scalar >::calculatedType())
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
@ NO_REGISTER
Do not request registration (bool: false).
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
A List with indirect addressing.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
static const Form zero
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition ZoneMesh.C:757
wordList names() const
A list of the zone names.
Definition ZoneMesh.C:463
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Lightweight class to store particle data derived from VOF calculations, with special handling for inp...
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
virtual const word & type() const =0
Runtime type information.
Generates particle size information from Eulerian calculations, e.g. VoF.
virtual void accumulateParticleInfo(const surfaceScalarField &alphaf, const surfaceScalarField &phi, const labelList &regionFaceIDs, const faceZone &fz)
Process latest region information.
virtual tmp< surfaceScalarField > phiU() const
Return the volumetric flux.
label nDiscardedParticles_
Total number of discarded particles.
label nInjectorLocations_
Number of sample locations to generate.
scalar alphaThreshold_
Value of phase fraction used to identify particle boundaries.
virtual void setBlockedFaces(const surfaceScalarField &alphaf, const faceZone &fz, boolList &blockedFaces)
Set the blocked faces, i.e. where alpha > alpha threshold value.
word UName_
Name of the velocity field, default = U.
globalIndex globalCoarseFaces_
Global coarse face addressing.
label nCollectedParticles_
Total number of collected particles.
Type faceValue(const GeometricField< Type, fvsPatchField, surfaceMesh > &field, const label localFaceI, const label globalFaceI) const
virtual void collectParticle(const scalar time, const label regioni)
Collect particles that have passed through the faceZone.
labelList patchFaceIDs_
Patch face indices where faceZone face intersect patch.
extractEulerianParticles(const word &name, const Time &runTime, const dictionary &dict)
Construct from components.
injectedParticleCloud cloud_
Storage for collected particles.
Map< label > regionToParticleMap_
Map from region to index in particles_ list.
labelList regions0_
Region indices in faceZone faces from last iteration.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
word phiName_
Name of the flux field, default ="rho".
word rhoName_
Name of the density field, default = rho.
List< eulerianParticle > particles_
Particle properties (partial, being accumulated).
virtual void initialiseBins()
Initialise the particle collection bins.
virtual void checkFaceZone()
Check that the faceZone is valid.
labelList fineToCoarseAddr_
Agglomeration addressing from fine to coarse (local proc only).
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function-object results.
labelList patchIDs_
Patch indices where faceZone face intersect patch.
virtual void calculateAddressing(const label nRegionsNew, const scalar time, labelList &regionFaceIDs)
Calculate the addressing between regions between iterations Returns the number of active regions (par...
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
void setProperty(const word &entryName, const Type &value)
Add generic property.
Type getProperty(const word &entryName, const Type &defaultValue=Type(Zero)) const
Retrieve generic property.
const Time & time() const
Return time database.
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
Definition writeFile.C:200
virtual bool read(const dictionary &dict)
Read.
Definition writeFile.C:240
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
Primarily stores particle properties so that it can be injected at a later time. Note that this store...
Primitive patch pair agglomerate method.
const labelList & restrictTopBottomAddressing() const
Return restriction from top level to bottom level.
void agglomerate()
Agglomerate patch.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition polyMesh.H:671
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition polyMesh.C:878
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Splits a patch into regions based on a mask field. Result is a globally consistent label list of regi...
label nRegions() const noexcept
Return the global number of regions.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
volScalarField & p
engineTime & runTime
autoPtr< surfaceVectorField > Uf
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
constexpr scalar pi(M_PI)
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Namespace for OpenFOAM.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr label labelMax
Definition label.H:55
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.