Loading...
Searching...
No Matches
ParticleCollector.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) 2012-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "ParticleCollector.H"
30#include "Pstream.H"
31#include "surfaceWriter.H"
32#include "unitConversion.H"
33#include "Random.H"
34#include "triangle.H"
35#include "cloud.H"
36
37// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38
39template<class CloudType>
40void Foam::ParticleCollector<CloudType>::makeLogFile
41(
42 const faceList& faces,
43 const Field<point>& points,
44 const Field<scalar>& area
45)
46{
47 // Create the output file if not already created
48 if (logToFile_)
49 {
50 DebugInfo<< "Creating output file" << endl;
51
52 if (Pstream::master())
53 {
54 // Create directory if does not exist
55 mkDir(this->writeTimeDir());
56
57 // Open new file at start up
58 outputFilePtr_.reset
59 (
60 new OFstream(this->writeTimeDir()/(type() + ".dat"))
61 );
62
63 outputFilePtr_()
64 << "# Source : " << type() << nl
65 << "# Bins : " << faces.size() << nl
66 << "# Total area : " << sum(area) << nl;
67
68 outputFilePtr_()
69 << "# Geometry :" << nl
70 << '#'
71 << tab << "Bin"
72 << tab << "(Centre_x Centre_y Centre_z)"
73 << tab << "Area"
74 << nl;
75
76 forAll(faces, i)
77 {
78 outputFilePtr_()
79 << '#'
80 << tab << i
81 << tab << faces[i].centre(points)
82 << tab << area[i]
83 << nl;
84 }
85
86 outputFilePtr_()
87 << '#' << nl
88 << "# Output format:" << nl;
89
90 forAll(faces, i)
91 {
92 word id = Foam::name(i);
93 word binId = "bin_" + id;
94
95 outputFilePtr_()
96 << '#'
97 << tab << "Time"
98 << tab << binId
99 << tab << "mass[" << id << "]"
100 << tab << "massFlowRate[" << id << "]"
101 << endl;
102 }
103 }
104 }
105}
106
107
108template<class CloudType>
109void Foam::ParticleCollector<CloudType>::initPolygons
110(
111 const List<Field<point>>& polygons
112)
113{
114 mode_ = mtPolygon;
115
116 label nPoints = 0;
117 forAll(polygons, polyI)
118 {
119 label np = polygons[polyI].size();
120 if (np < 3)
121 {
122 FatalIOErrorInFunction(this->coeffDict())
123 << "polygons must consist of at least 3 points"
124 << exit(FatalIOError);
125 }
126
127 nPoints += np;
128 }
129
130 label pointOffset = 0;
131 points_.setSize(nPoints);
132 faces_.setSize(polygons.size());
133 faceTris_.setSize(polygons.size());
134 area_.setSize(polygons.size());
135 forAll(faces_, facei)
136 {
137 const Field<point>& polyPoints = polygons[facei];
138 face f(identity(polyPoints.size(), pointOffset));
139 UIndirectList<point>(points_, f) = polyPoints;
140 area_[facei] = f.mag(points_);
141
142 DynamicList<face> tris;
143 f.triangles(points_, tris);
144 faceTris_[facei].transfer(tris);
145
146 faces_[facei].transfer(f);
147
148 pointOffset += polyPoints.size();
149 }
150}
151
152
153template<class CloudType>
154void Foam::ParticleCollector<CloudType>::initConcentricCircles()
155{
156 mode_ = mtConcentricCircle;
157
158 vector origin(this->coeffDict().lookup("origin"));
159
160 this->coeffDict().readEntry("radius", radius_);
161 this->coeffDict().readEntry("nSector", nSector_);
162
163 label nS = nSector_;
164
165 vector refDir;
166 if (nSector_ > 1)
167 {
168 this->coeffDict().readEntry("refDir", refDir);
169
170 refDir -= normal_[0]*(normal_[0] & refDir);
171 refDir.normalise();
172 }
173 else
174 {
175 // Set 4 quadrants for single sector cases
176 nS = 4;
177
178 vector tangent = Zero;
179 scalar magTangent = 0.0;
180
181 Random& rnd = this->owner().rndGen();
182 while (magTangent < SMALL)
183 {
184 vector v = rnd.sample01<vector>();
185
186 tangent = v - (v & normal_[0])*normal_[0];
187 magTangent = mag(tangent);
188 }
189
190 refDir = tangent/magTangent;
191 }
192
193 scalar dTheta = 5.0;
194 scalar dThetaSector = 360.0/scalar(nS);
195 label intervalPerSector = max(1, ceil(dThetaSector/dTheta));
196 dTheta = dThetaSector/scalar(intervalPerSector);
197
198 label nPointPerSector = intervalPerSector + 1;
199
200 label nPointPerRadius = nS*(nPointPerSector - 1);
201 label nPoint = radius_.size()*nPointPerRadius;
202 label nFace = radius_.size()*nS;
203
204 // Add origin
205 nPoint++;
206
207 points_.setSize(nPoint);
208 faces_.setSize(nFace);
209 area_.setSize(nFace);
210
211 coordSys_ = coordSystem::cylindrical(origin, normal_[0], refDir);
212
213 List<label> ptIDs(identity(nPointPerRadius));
214
215 points_[0] = origin;
216
217 // Points
218 forAll(radius_, radI)
219 {
220 label pointOffset = radI*nPointPerRadius + 1;
221
222 for (label i = 0; i < nPointPerRadius; i++)
223 {
224 label pI = i + pointOffset;
225 point pCyl(radius_[radI], degToRad(i*dTheta), 0.0);
226 points_[pI] = coordSys_.globalPosition(pCyl);
227 }
228 }
229
230 // Faces
231 DynamicList<label> facePts(2*nPointPerSector);
232 forAll(radius_, radI)
233 {
234 if (radI == 0)
235 {
236 for (label secI = 0; secI < nS; secI++)
237 {
238 facePts.clear();
239
240 // Append origin point
241 facePts.append(0);
242
243 for (label ptI = 0; ptI < nPointPerSector; ptI++)
244 {
245 label i = ptI + secI*(nPointPerSector - 1);
246 label id = ptIDs.fcIndex(i - 1) + 1;
247 facePts.append(id);
248 }
249
250 label facei = secI + radI*nS;
251
252 faces_[facei] = face(facePts);
253 area_[facei] = faces_[facei].mag(points_);
254 }
255 }
256 else
257 {
258 for (label secI = 0; secI < nS; secI++)
259 {
260 facePts.clear();
261
262 label offset = (radI - 1)*nPointPerRadius + 1;
263
264 for (label ptI = 0; ptI < nPointPerSector; ptI++)
265 {
266 label i = ptI + secI*(nPointPerSector - 1);
267 label id = offset + ptIDs.fcIndex(i - 1);
268 facePts.append(id);
269 }
270 for (label ptI = nPointPerSector-1; ptI >= 0; ptI--)
271 {
272 label i = ptI + secI*(nPointPerSector - 1);
273 label id = offset + nPointPerRadius + ptIDs.fcIndex(i - 1);
274 facePts.append(id);
275 }
276
277 label facei = secI + radI*nS;
278
279 faces_[facei] = face(facePts);
280 area_[facei] = faces_[facei].mag(points_);
281 }
282 }
283 }
284}
285
286
287template<class CloudType>
288void Foam::ParticleCollector<CloudType>::collectParcelPolygon
289(
290 const point& p1,
291 const point& p2
292) const
293{
294 forAll(faces_, facei)
295 {
296 const label facePoint0 = faces_[facei][0];
297
298 const point& pf = points_[facePoint0];
299
300 const scalar d1 = normal_[facei] & (p1 - pf);
301 const scalar d2 = normal_[facei] & (p2 - pf);
302
303 if (sign(d1) == sign(d2))
304 {
305 // Did not cross polygon plane
306 continue;
307 }
308
309 // Intersection point
310 const point pIntersect = p1 + (d1/(d1 - d2))*(p2 - p1);
311
312 // Identify if point is within the bounds of the face. Create triangles
313 // between the intersection point and each edge of the face. If all the
314 // triangle normals point in the same direction as the face normal, then
315 // the particle is within the face. Note that testing for pointHits on
316 // the face's decomposed triangles does not work due to ambiguity along
317 // the diagonals.
318 const face& f = faces_[facei];
319 const vector areaNorm = f.areaNormal(points_);
320 bool inside = true;
321 for (label i = 0; i < f.size(); ++i)
322 {
323 const label j = f.fcIndex(i);
324 const triPointRef t(pIntersect, points_[f[i]], points_[f[j]]);
325
326 if ((areaNorm & t.areaNormal()) < 0)
327 {
328 inside = false;
329 break;
330 }
331 }
332
333 // Add to the list of hits
334 if (inside)
335 {
336 hitFaceIDs_.append(facei);
337 }
338 }
339}
340
341
342template<class CloudType>
343void Foam::ParticleCollector<CloudType>::collectParcelConcentricCircles
344(
345 const point& p1,
346 const point& p2
347) const
348{
349 label secI = -1;
350
351 const scalar d1 = normal_[0] & (p1 - coordSys_.origin());
352 const scalar d2 = normal_[0] & (p2 - coordSys_.origin());
353
354 if (sign(d1) == sign(d2))
355 {
356 // Did not cross plane
357 return;
358 }
359
360 // Intersection point in cylindrical coordinate system
361 const point pCyl = coordSys_.localPosition(p1 + (d1/(d1 - d2))*(p2 - p1));
362
363 scalar r = pCyl[0];
364
365 if (r < radius_.last())
366 {
367 label radI = 0;
368 while (r > radius_[radI])
369 {
370 radI++;
371 }
372
373 if (nSector_ == 1)
374 {
375 secI = 4*radI;
376 }
377 else
378 {
379 scalar theta = pCyl[1] + constant::mathematical::pi;
380
381 secI =
382 nSector_*radI
383 + floor
384 (
385 scalar(nSector_)*theta/constant::mathematical::twoPi
386 );
387 }
388 }
389
390 if (secI != -1)
391 {
392 hitFaceIDs_.append(secI);
394}
395
396
397// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
398
399template<class CloudType>
401{
402 const fvMesh& mesh = this->owner().mesh();
403 const Time& time = mesh.time();
404 scalar timeNew = time.value();
405 scalar timeElapsed = timeNew - timeOld_;
406
407 totalTime_ += timeElapsed;
408
409 const scalar alpha = (totalTime_ - timeElapsed)/totalTime_;
410 const scalar beta = timeElapsed/totalTime_;
411
412 forAll(faces_, facei)
413 {
414 massFlowRate_[facei] =
415 alpha*massFlowRate_[facei] + beta*mass_[facei]/timeElapsed;
416 massTotal_[facei] += mass_[facei];
417 }
418
419 Log_<< type() << " output:" << nl;
420
421 Field<scalar> faceMassTotal(mass_.size(), Zero);
422 this->getModelProperty("massTotal", faceMassTotal);
423
424 Field<scalar> faceMassFlowRate(massFlowRate_.size(), Zero);
425 this->getModelProperty("massFlowRate", faceMassFlowRate);
426
427
428 scalar sumTotalMass = 0.0;
429 scalar sumAverageMFR = 0.0;
430 forAll(faces_, facei)
431 {
432 faceMassTotal[facei] +=
433 returnReduce(massTotal_[facei], sumOp<scalar>());
434
435 faceMassFlowRate[facei] +=
436 returnReduce(massFlowRate_[facei], sumOp<scalar>());
437
438 sumTotalMass += faceMassTotal[facei];
439 sumAverageMFR += faceMassFlowRate[facei];
440
441 if (outputFilePtr_)
442 {
443 outputFilePtr_()
444 << time.timeName()
445 << tab << facei
446 << tab << faceMassTotal[facei]
447 << tab << faceMassFlowRate[facei]
448 << endl;
449 }
450 }
451
452 Log_<< " sum(total mass) = " << sumTotalMass << nl
453 << " sum(average mass flow rate) = " << sumAverageMFR << nl
454 << endl;
455
456
457 if (Pstream::master() && surfaceFormat_ != "none")
458 {
459 auto writer = surfaceWriter::New
460 (
461 surfaceFormat_,
462 surfaceWriter::formatOptions(this->coeffDict(), surfaceFormat_)
463 );
464
465 if (debug)
466 {
467 writer->verbose(true);
468 }
469
470 writer->open
471 (
472 points_,
473 faces_,
474 (this->writeTimeDir() / "collector"),
475 false // serial - already merged
476 );
477
478 writer->nFields(2); // Legacy VTK
479 writer->write("massFlowRate", faceMassFlowRate);
480 writer->write("massTotal", faceMassTotal);
481 }
482
483
484 if (resetOnWrite_)
485 {
486 Field<scalar> dummy(faceMassTotal.size(), Zero);
487 this->setModelProperty("massTotal", dummy);
488 this->setModelProperty("massFlowRate", dummy);
489
490 timeOld_ = timeNew;
491 totalTime_ = 0.0;
492 }
493 else
494 {
495 this->setModelProperty("massTotal", faceMassTotal);
496 this->setModelProperty("massFlowRate", faceMassFlowRate);
497 }
498
499 forAll(faces_, facei)
500 {
501 mass_[facei] = 0.0;
502 massTotal_[facei] = 0.0;
503 massFlowRate_[facei] = 0.0;
505}
506
507
508// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
509
510template<class CloudType>
512(
513 const dictionary& dict,
514 CloudType& owner,
515 const word& modelName
516)
517:
518 CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
519 mode_(mtUnknown),
520 parcelType_(this->coeffDict().getOrDefault("parcelType", -1)),
521 removeCollected_(this->coeffDict().getBool("removeCollected")),
522 resetOnWrite_(this->coeffDict().getBool("resetOnWrite")),
523 logToFile_(this->coeffDict().getBool("log")),
524 points_(),
525 faces_(),
526 faceTris_(),
527 nSector_(0),
528 radius_(),
529 coordSys_(),
530 normal_(),
531 negateParcelsOppositeNormal_
532 (
533 this->coeffDict().getBool("negateParcelsOppositeNormal")
534 ),
535 surfaceFormat_(this->coeffDict().getWord("surfaceFormat")),
536 totalTime_(0.0),
537 mass_(),
538 massTotal_(),
539 massFlowRate_(),
540 outputFilePtr_(),
541 timeOld_(owner.mesh().time().value()),
542 hitFaceIDs_()
543{
544 normal_ /= mag(normal_);
545
546 word mode(this->coeffDict().lookup("mode"));
547 if (mode == "polygon")
548 {
549 List<Field<point>> polygons(this->coeffDict().lookup("polygons"));
550
551 initPolygons(polygons);
552
553 vector n0(this->coeffDict().lookup("normal"));
554 normal_ = vectorField(faces_.size(), n0);
555 }
556 else if (mode == "polygonWithNormal")
557 {
558 List<Tuple2<Field<point>, vector>> polygonAndNormal
559 (
560 this->coeffDict().lookup("polygons")
561 );
562
563 List<Field<point>> polygons(polygonAndNormal.size());
564 normal_.setSize(polygonAndNormal.size());
565
566 forAll(polygons, polyI)
567 {
568 polygons[polyI] = polygonAndNormal[polyI].first();
569 normal_[polyI] = normalised(polygonAndNormal[polyI].second());
570 }
571
572 initPolygons(polygons);
573 }
574 else if (mode == "concentricCircle")
575 {
576 vector n0(this->coeffDict().lookup("normal"));
577 normal_ = vectorField(1, n0);
578
579 initConcentricCircles();
580 }
581 else
582 {
584 << "Unknown mode " << mode << ". Available options are "
585 << "polygon, polygonWithNormal and concentricCircle"
586 << exit(FatalIOError);
587 }
588
589 mass_.setSize(faces_.size(), 0.0);
590 massTotal_.setSize(faces_.size(), 0.0);
591 massFlowRate_.setSize(faces_.size(), 0.0);
592
593 makeLogFile(faces_, points_, area_);
594}
595
596
597template<class CloudType>
599(
601)
602:
604 mode_(pc.mode_),
605 parcelType_(pc.parcelType_),
606 removeCollected_(pc.removeCollected_),
607 resetOnWrite_(pc.resetOnWrite_),
608 logToFile_(pc.logToFile_),
609 points_(pc.points_),
610 faces_(pc.faces_),
611 faceTris_(pc.faceTris_),
612 nSector_(pc.nSector_),
613 radius_(pc.radius_),
614 coordSys_(pc.coordSys_),
615 area_(pc.area_),
616 normal_(pc.normal_),
617 negateParcelsOppositeNormal_(pc.negateParcelsOppositeNormal_),
618 surfaceFormat_(pc.surfaceFormat_),
619 totalTime_(pc.totalTime_),
620 mass_(pc.mass_),
621 massTotal_(pc.massTotal_),
622 massFlowRate_(pc.massFlowRate_),
623 outputFilePtr_(),
624 timeOld_(0.0),
625 hitFaceIDs_()
626{}
627
628
629// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
630
631template<class CloudType>
633(
634 parcelType& p,
635 const scalar dt,
636 const point& position0,
637 const typename parcelType::trackingData& td
638)
639{
640 bool keepParticle = true;
641
642 if ((parcelType_ != -1) && (parcelType_ != p.typeId()))
643 {
644 return keepParticle;
645 }
646
647 hitFaceIDs_.clear();
648
649 switch (mode_)
650 {
651 case mtPolygon:
652 {
653 collectParcelPolygon(position0, p.position());
654 break;
655 }
656 case mtConcentricCircle:
657 {
658 collectParcelConcentricCircles(position0, p.position());
659 break;
660 }
661 default:
662 {}
663 }
664
665 forAll(hitFaceIDs_, i)
666 {
667 label facei = hitFaceIDs_[i];
668 scalar m = p.nParticle()*p.mass();
669
670 if (negateParcelsOppositeNormal_)
671 {
672 scalar Unormal = 0;
673 vector Uhat = p.U();
674 switch (mode_)
675 {
676 case mtPolygon:
677 {
678 Unormal = Uhat & normal_[facei];
679 break;
680 }
681 case mtConcentricCircle:
682 {
683 Unormal = Uhat & normal_[0];
684 break;
685 }
686 default:
687 {}
688 }
689
690 Uhat /= mag(Uhat) + ROOTVSMALL;
691
692 if (Unormal < 0)
693 {
694 m = -m;
695 }
696 }
697
698 // Add mass contribution
699 mass_[facei] += m;
700
701 if (nSector_ == 1)
702 {
703 mass_[facei + 1] += m;
704 mass_[facei + 2] += m;
705 mass_[facei + 3] += m;
706 }
707
708 if (removeCollected_)
709 {
710 keepParticle = false;
711 }
712 }
713
714 return keepParticle;
715}
716
717
718// ************************************************************************* //
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edgesCentres")))
Templated cloud function object base class.
fileName writeTimeDir() const
Return the output time path.
CloudFunctionObject(CloudType &owner)
Construct null from owner.
const CloudType & owner() const
Return const access to the owner cloud.
particle::trackingData trackingData
Definition DSMCParcel.H:155
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
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
Function object to collect the parcel mass- and mass flow rate over a set of polygons....
ParticleCollector(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
void write()
Write post-processing info.
virtual bool postMove(parcelType &p, const scalar dt, const point &position0, const typename parcelType::trackingData &td)
Post-move hook.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
T & first()
Access first element of the list, position [0].
Definition UList.H:957
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const Type & value() const noexcept
Return const reference to value.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Lookup type of boundary radiation properties.
Definition lookup.H:60
bool getModelProperty(const word &entryName, Type &value) const
Retrieve generic property from the sub-model.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
const word & modelName() const
Return const access to the name of the sub-model.
void setModelProperty(const word &entryName, const Type &value)
Add generic property to the sub-model.
static autoPtr< surfaceWriter > New(const word &writeType)
Select construct a surfaceWriter.
static dictionary formatOptions(const dictionary &dict, const word &formatName, const word &entryName="formatOptions")
Same as fileFormats::getFormatOptions.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
const pointField & points
label nPoints
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
return returnReduce(nRefine-oldNRefine, sumOp< label >())
#define Log_
Report write to Foam::Info if the class log switch is true.
#define DebugInfo
Report an information message using Foam::Info.
Namespace for handling debugging switches.
Definition debug.C:45
const wordList area
Standard area field types (scalar, vector, tensor, etc).
type
Types of root.
Definition Roots.H:53
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
DSMCCloud< dsmcParcel > CloudType
dimensionedScalar sign(const dimensionedScalar &ds)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition POSIX.C:775
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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...
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
vector point
Point is a vector.
Definition point.H:37
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
labelList f(nPoints)
volScalarField & alpha
dictionary dict
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.
mkDir(pdfPath)