Loading...
Searching...
No Matches
reactingOneDim.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2021 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 "reactingOneDim.H"
31#include "fvm.H"
32#include "fvcDiv.H"
33#include "fvcVolumeIntegrate.H"
34#include "fvcLaplacian.H"
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39namespace Foam
40{
41namespace regionModels
42{
45
46// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
47
49
52
53// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
54
55void reactingOneDim::readReactingOneDimControls()
56{
57 const dictionary& solution = this->solution().subDict("SIMPLE");
58 solution.readEntry("nNonOrthCorr", nNonOrthCorr_);
59 time().controlDict().readEntry("maxDi", maxDiff_);
60 coeffs().readEntry("minimumDelta", minimumDelta_);
61 gasHSource_ = coeffs().getOrDefault("gasHSource", false);
62 coeffs().readEntry("qrHSource", qrHSource_);
64 coeffs().getOrDefault("useChemistrySolvers", true);
65}
66
67
69{
71 {
72 readReactingOneDimControls();
73 return true;
74 }
75
76 return false;
77}
78
79
81{
83 {
84 readReactingOneDimControls();
85 return true;
86 }
87
88 return false;
89}
90
91
93{
94 // Update local qr from coupled qr field
95 qr_ == Zero;
96
97 // Retrieve field from coupled region using mapped boundary conditions
99
101
103 {
104 const label patchi = intCoupledPatchIDs_[i];
105
106 // qr is positive going in the solid
107 // If the surface is emitting the radiative flux is set to zero
108 qrBf[patchi] = max(qrBf[patchi], scalar(0));
109 }
110
111 const vectorField& cellC = regionMesh().cellCentres();
112
114
115 // Propagate qr through 1-D regions
116 label localPyrolysisFacei = 0;
118 {
119 const label patchi = intCoupledPatchIDs_[i];
120
121 const scalarField& qrp = qr_.boundaryField()[patchi];
122 const vectorField& Cf = regionMesh().Cf().boundaryField()[patchi];
123
124 forAll(qrp, facei)
125 {
126 const scalar qr0 = qrp[facei];
127 point Cf0 = Cf[facei];
128 const labelList& cells = boundaryFaceCells_[localPyrolysisFacei++];
129 scalar kappaInt = 0.0;
130 forAll(cells, k)
131 {
132 const label celli = cells[k];
133 const point& Cf1 = cellC[celli];
134 const scalar delta = mag(Cf1 - Cf0);
135 kappaInt += kappa()[celli]*delta;
136 qr_[celli] = qr0*exp(-kappaInt);
137 Cf0 = Cf1;
138 }
139 }
140 }
141}
142
143
145{
146 phiHsGas_ == Zero;
147 phiGas_ == Zero;
148
149 const speciesTable& gasTable = solidChemistry_->gasTable();
150
151 forAll(gasTable, gasI)
152 {
153 tmp<volScalarField> tHsiGas =
154 solidChemistry_->gasHs(solidThermo_->p(), solidThermo_->T(), gasI);
155
156 const volScalarField& HsiGas = tHsiGas();
157
158 const volScalarField::Internal& RRiGas = solidChemistry_->RRg(gasI);
159
161
162 label totalFaceId = 0;
164 {
165 const label patchi = intCoupledPatchIDs_[i];
166
167 scalarField& phiGasp = phiGasBf[patchi];
168 const scalarField& cellVol = regionMesh().V();
169
170 forAll(phiGasp, facei)
171 {
172 const labelList& cells = boundaryFaceCells_[totalFaceId++];
173 scalar massInt = 0.0;
175 {
176 const label celli = cells[k];
177 massInt += RRiGas[celli]*cellVol[celli];
178 phiHsGas_[celli] += massInt*HsiGas[celli];
179 }
180
181 phiGasp[facei] += massInt;
182
183 if (debug)
184 {
185 Info<< " Gas : " << gasTable[gasI]
186 << " on patch : " << patchi
187 << " mass produced at face(local) : "
188 << facei
189 << " is : " << massInt
190 << " [kg/s] " << endl;
192 }
193 }
194 }
195}
196
197
199{
200 if (qrHSource_)
201 {
203 }
204
205 updatePhiGas();
206}
207
208
209void reactingOneDim::updateMesh(const scalarField& deltaV)
210{
211 Info<< "Initial/final volumes = " << gSum(deltaV) << endl;
212
213 // Move the mesh
214 const labelList moveMap = moveMesh(deltaV, minimumDelta_);
215
216 // Flag any cells that have not moved as non-reacting
217 forAll(moveMap, i)
218 {
219 if (moveMap[i] == 1)
221 solidChemistry_->setCellReacting(i, false);
222 }
223 }
224}
225
226
228{
230
231 if (!moveMesh_)
232 {
233 fvScalarMatrix rhoEqn
234 (
235 fvm::ddt(rho_) == -solidChemistry_->RRg()
236 );
237
238 rhoEqn.solve();
239 }
240 else
241 {
242 const scalarField deltaV
243 (
245 );
246
247 updateMesh(deltaV);
248 }
249}
250
251
253{
255
256 volScalarField Yt(0.0*Ys_[0]);
257
258 for (label i=0; i<Ys_.size()-1; i++)
259 {
260 volScalarField& Yi = Ys_[i];
261
262 fvScalarMatrix YiEqn
263 (
264 fvm::ddt(rho_, Yi) == solidChemistry_->RRs(i)
265 );
266
267 if (regionMesh().moving())
268 {
269 surfaceScalarField phiYiRhoMesh
270 (
272 );
273
274 YiEqn -= fvc::div(phiYiRhoMesh);
275
276 }
277
278 YiEqn.solve(regionMesh().solver("Yi"));
279 Yi.clamp_min(0);
280 Yt += Yi;
282
283 Ys_[Ys_.size() - 1] = 1.0 - Yt;
284
285}
286
287
289{
291
293
294 fvScalarMatrix hEqn
295 (
299 - fvc::laplacian(kappa(), T())
300 ==
302 + solidChemistry_->RRsHs()
303 );
304
305/*
306 NOTE: gas Hs is included in hEqn
307
308 if (gasHSource_)
309 {
310 const surfaceScalarField phiGas(fvc::interpolate(phiHsGas_));
311 hEqn += fvc::div(phiGas);
312 }
313*/
314
315 if (qrHSource_)
316 {
318 hEqn += fvc::div(phiqr);
319 }
320
321/*
322 NOTE: The moving mesh option is only correct for reaction such as
323 Solid -> Gas, thus the ddt term is compensated exactly by chemistrySh and
324 the mesh flux is not necessary.
325
326 if (regionMesh().moving())
327 {
328 surfaceScalarField phihMesh
329 (
330 fvc::interpolate(rho_*h_)*regionMesh().phi()
331 );
332
333 hEqn -= fvc::div(phihMesh);
335*/
336 hEqn.relax();
337 hEqn.solve();
338}
339
340
342{
343 /*
344 totalGasMassFlux_ = 0;
345 forAll(intCoupledPatchIDs_, i)
346 {
347 const label patchi = intCoupledPatchIDs_[i];
348 totalGasMassFlux_ += gSum(phiGas_.boundaryField()[patchi]);
349 }
350 */
351
352 if (infoOutput_)
353 {
355
360 }
361}
362
363
364// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
365
366reactingOneDim::reactingOneDim
367(
368 const word& modelType,
369 const fvMesh& mesh,
370 const word& regionType
371)
372:
373 pyrolysisModel(modelType, mesh, regionType),
374 solidThermo_(solidReactionThermo::New(regionMesh())),
375 solidChemistry_(basicSolidChemistryModel::New(solidThermo_())),
376 radiation_(radiation::radiationModel::New(solidThermo_->T())),
377 rho_
378 (
380 (
381 "rho",
382 regionMesh().time().timeName(),
383 regionMesh(),
384 IOobject::NO_READ,
385 IOobject::AUTO_WRITE
386 ),
387 solidThermo_->rho()
388 ),
389 Ys_(solidThermo_->composition().Y()),
390 h_(solidThermo_->he()),
391 nNonOrthCorr_(-1),
392 maxDiff_(10),
393 minimumDelta_(1e-4),
394
395 phiGas_
396 (
398 (
399 "phiGas",
400 time().timeName(),
401 regionMesh(),
402 IOobject::READ_IF_PRESENT,
403 IOobject::AUTO_WRITE
404 ),
405 regionMesh(),
407 ),
408
409 phiHsGas_
410 (
412 (
413 "phiHsGas",
414 time().timeName(),
415 regionMesh(),
416 IOobject::NO_READ,
417 IOobject::NO_WRITE
418 ),
419 regionMesh(),
421 ),
422
423 chemistryQdot_
424 (
426 (
427 "chemistryQdot",
428 time().timeName(),
429 regionMesh(),
430 IOobject::NO_READ,
431 IOobject::AUTO_WRITE
432 ),
433 regionMesh(),
435 ),
436
437 qr_
438 (
440 (
441 "qr",
442 time().timeName(),
443 regionMesh(),
444 IOobject::MUST_READ,
445 IOobject::AUTO_WRITE
446 ),
447 regionMesh()
448 ),
449
450 lostSolidMass_(dimensionedScalar(dimMass, Zero)),
451 addedGasMass_(dimensionedScalar(dimMass, Zero)),
452 totalGasMassFlux_(0.0),
453 totalHeatRR_(dimensionedScalar(dimEnergy/dimTime, Zero)),
454 gasHSource_(false),
455 qrHSource_(false),
456 useChemistrySolvers_(true)
457{
459 {
460 read();
461 }
462}
463
464
465reactingOneDim::reactingOneDim
466(
467 const word& modelType,
468 const fvMesh& mesh,
469 const dictionary& dict,
470 const word& regionType
471)
472:
473 pyrolysisModel(modelType, mesh, dict, regionType),
474 solidThermo_(solidReactionThermo::New(regionMesh())),
475 solidChemistry_(basicSolidChemistryModel::New(solidThermo_())),
476 radiation_(radiation::radiationModel::New(solidThermo_->T())),
477 rho_
478 (
480 (
481 "rho",
482 regionMesh().time().timeName(),
483 regionMesh(),
484 IOobject::NO_READ,
485 IOobject::AUTO_WRITE
486 ),
487 solidThermo_->rho()
488 ),
489 Ys_(solidThermo_->composition().Y()),
490 h_(solidThermo_->he()),
491 nNonOrthCorr_(-1),
492 maxDiff_(10),
493 minimumDelta_(1e-4),
494
495 phiGas_
496 (
498 (
499 "phiGas",
500 time().timeName(),
501 regionMesh(),
502 IOobject::NO_READ,
503 IOobject::AUTO_WRITE
504 ),
505 regionMesh(),
507 ),
508
509 phiHsGas_
510 (
512 (
513 "phiHsGas",
514 time().timeName(),
515 regionMesh(),
516 IOobject::NO_READ,
517 IOobject::NO_WRITE
518 ),
519 regionMesh(),
521 ),
522
523 chemistryQdot_
524 (
526 (
527 "chemistryQdot",
528 time().timeName(),
529 regionMesh(),
530 IOobject::NO_READ,
531 IOobject::AUTO_WRITE
532 ),
533 regionMesh(),
535 ),
536
537 qr_
538 (
540 (
541 "qr",
542 time().timeName(),
543 regionMesh(),
544 IOobject::MUST_READ,
545 IOobject::AUTO_WRITE
546 ),
547 regionMesh()
548 ),
549
550 lostSolidMass_(dimensionedScalar(dimMass, Zero)),
551 addedGasMass_(dimensionedScalar(dimMass, Zero)),
552 totalGasMassFlux_(0.0),
553 totalHeatRR_(dimensionedScalar(dimEnergy/dimTime, Zero)),
554 gasHSource_(false),
555 qrHSource_(false),
556 useChemistrySolvers_(true)
557{
558 if (active_)
559 {
561 }
562}
563
564
565// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
568{}
569
570
571// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
572
573scalar reactingOneDim::addMassSources(const label patchi, const label facei)
574{
575 label index = 0;
577 {
578 if (primaryPatchIDs_[i] == patchi)
579 {
580 index = i;
581 break;
582 }
583 }
584
585 const label localPatchId = intCoupledPatchIDs_[index];
586
587 const scalar massAdded = phiGas_.boundaryField()[localPatchId][facei];
588
589 if (debug)
590 {
591 Info<< "\nPyrolysis region: " << type() << "added mass : "
592 << massAdded << endl;
593 }
594
595 return massAdded;
596}
597
598
600{
601 scalar DiNum = -GREAT;
602
603 if (regionMesh().nInternalFaces() > 0)
604 {
605 surfaceScalarField KrhoCpbyDelta
606 (
610 );
611
612 DiNum = max(KrhoCpbyDelta.primitiveField())*time().deltaTValue();
613 }
614
616}
617
619scalar reactingOneDim::maxDiff() const
620{
621 return maxDiff_;
622}
623
626{
627 return rho_;
628}
629
632{
633 return solidThermo_->T();
634}
635
638{
639 return solidThermo_->Cp();
640}
641
644{
645 return radiation_->absorptionEmission().a();
646}
647
650{
651 return solidThermo_->kappa();
652}
653
656{
657 return phiGas_;
658}
659
662{
664}
665
666
668{
669 Info<< "\nEvolving pyrolysis in region: " << regionMesh().name() << endl;
670
672 {
673 solidChemistry_->solve(time().deltaTValue());
674 }
675 else
676 {
677 solidChemistry_->calculate();
678 }
679
681
683
684 updateFields();
685
687
688 for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
689 {
690 solveEnergy();
691 }
692
694
695 solidThermo_->correct();
696
697 auto limits = gMinMax(solidThermo_->T().primitiveField());
698
699 Info<< "pyrolysis min/max(T) = "
700 << limits.min() << ", " << limits.max() << endl;
701}
702
703
705{
706 Info<< "\nPyrolysis in region: " << regionMesh().name() << endl;
707
708 Info<< indent << "Total gas mass produced [kg] = "
709 << addedGasMass_.value() << nl
710 << indent << "Total solid mass lost [kg] = "
711 << lostSolidMass_.value() << nl
712 //<< indent << "Total pyrolysis gases [kg/s] = "
713 //<< totalGasMassFlux_ << nl
714 << indent << "Total heat release rate [J/s] = "
715 << totalHeatRR_.value() << nl;
716}
717
718
719// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
720
721} // End namespace Foam
722} // End namespace regionModels
723} // End namespace pyrolysisModels
724
725// ************************************************************************* //
label k
scalar delta
volScalarField Yt(0.0 *Y[0])
volScalarField & he
Definition YEEqn.H:52
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
DimensionedField< scalar, volMesh > Internal
void clamp_min(const Type &lower)
Impose lower (floor) clamp on the field values (in-place).
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
void correctBoundaryConditions()
Correct boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
scalar deltaTValue() const noexcept
Return time step value.
Definition TimeStateI.H:49
dimensionedScalar deltaT() const
Return time step.
Definition TimeStateI.H:61
const dictionary & controlDict() const noexcept
Return read access to the controlDict dictionary.
Definition Time.H:539
Chemistry model for solid thermodynamics.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
const Type & value() const noexcept
Return const reference to value.
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition fvMatrix.C:1094
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const word & name() const
Return reference to name.
Definition fvMesh.H:387
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
const vectorField & cellCentres() const
Top level model for radiation modelling.
static autoPtr< pyrolysisModel > New(const fvMesh &mesh, const word &regionType="pyrolysis")
Return a reference to the selected pyrolysis model.
virtual bool read()
Read control parameters.
Reacting, 1-D pyrolysis model NOTE: The moveMesh option can only be applied to solid reaction such as...
void solveSpeciesMass()
Solve solid species mass conservation.
bool qrHSource_
Add in depth radiation source term.
virtual scalar solidRegionDiffNo() const
Mean diffusion number of the solid region.
void solveContinuity()
Solve continuity equation.
const volScalarField & rho() const
Fields.
virtual tmp< volScalarField > kappa() const
Return the region thermal conductivity [W/m/k].
dimensionedScalar lostSolidMass_
Cumulative lost mass of the condensed phase [kg].
void updateqr()
Update radiative flux in pyrolysis region.
scalar minimumDelta_
Minimum delta for combustion.
bool useChemistrySolvers_
Use chemistry solvers (ode or sequential).
autoPtr< solidReactionThermo > solidThermo_
Reference to solid thermo.
PtrList< volScalarField > & Ys_
List of solid components.
volScalarField phiHsGas_
Sensible enthalpy gas flux [J/m2/s].
virtual const volScalarField & T() const
Return const temperature [K].
dimensionedScalar addedGasMass_
Cumulative mass generation of the gas phase [kg].
virtual tmp< volScalarField > kappaRad() const
Return the region absorptivity [1/m].
volScalarField qr_
Coupled region radiative heat flux [W/m2].
virtual const tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
virtual scalar addMassSources(const label patchi, const label facei)
External hook to add mass to the primary region.
void updatePhiGas()
Update enthalpy flux for pyrolysis gases.
autoPtr< basicSolidChemistryModel > solidChemistry_
Reference to the solid chemistry model.
void updateMesh(const scalarField &mass0)
Update/move mesh based on change in mass.
virtual void preEvolveRegion()
Pre-evolve region.
autoPtr< radiation::radiationModel > radiation_
Pointer to radiation model.
dimensionedScalar totalHeatRR_
Total heat release rate [J/s].
surfaceScalarField phiGas_
Total gas mass flux to the primary region [kg/m2/s].
virtual void info()
Provide some feedback.
virtual const surfaceScalarField & phiGas() const
Return the total gas mass flux to primary region [kg/m2/s].
virtual scalar maxDiff() const
Return max diffusivity allowed in the solid.
bool gasHSource_
Add gas enthalpy source term.
label nNonOrthCorr_
Number of non-orthogonal correctors.
bool read()
Read control parameters from dictionary.
volScalarField chemistryQdot_
Heat release rate [J/s/m3].
virtual void evolveRegion()
Evolve the pyrolysis equations.
scalar totalGasMassFlux_
Total mass gas flux at the pyrolysing walls [kg/s].
tmp< labelField > moveMesh(const scalarList &deltaV, const scalar minDelta=0.0)
Move mesh points according to change in cell volumes.
const surfaceScalarField & nMagSf() const
Return the face area magnitudes / [m2].
Switch moveMesh_
Flag to allow mesh movement.
labelListList boundaryFaceCells_
Global cell IDs.
const Time & time_
Reference to the time database.
Definition regionModel.H:97
const dictionary & solution() const
Return the solution dictionary.
Switch infoOutput_
Active information output.
labelList primaryPatchIDs_
List of patch IDs on the primary region coupled to this region.
const Time & time() const noexcept
Return the reference to the time database.
const fvMesh & regionMesh() const
Return the region mesh database.
virtual void preEvolveRegion()
Pre-evolve region.
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
const dictionary & coeffs() const noexcept
Return the model coefficients dictionary.
Foam::solidReactionThermo.
Selector class for relaxation factors, solver type and solution.
Definition solution.H:95
Base solver class.
Definition solver.H:48
virtual const surfaceScalarField & deltaCoeffs() const
Return reference to cell-centre difference coefficients.
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
basicSpecieMixture & composition
PtrList< volScalarField > & Y
scalar DiNum
auto limits
Definition setRDeltaT.H:186
dynamicFvMesh & mesh
Calculate the divergence of the given field.
Calculate the laplacian of the given field.
Volume integrate volField creating a volField.
const cellShapeList & cells
word timeName
Definition getTimeIndex.H:3
#define DebugInFunction
Report an information message using Foam::Info.
Namespace for handling debugging switches.
Definition debug.C:45
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.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
Namespace for radiation modelling.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
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.
Type gSum(const FieldField< Field, Type > &f)
dimensionedScalar exp(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimEnergy
fvMatrix< scalar > fvScalarMatrix
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
hashedWordList speciesTable
A table of species as a hashedWordList.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
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
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition stdFoam.H:315