massDiffusionLimitedPhaseChange.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2024-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "fvmSup.H"
28 #include "multiphaseEuler.H"
32 
33 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace fv
38 {
41  (
42  fvModel,
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 void Foam::fv::massDiffusionLimitedPhaseChange::readCoeffs
53 (
54  const dictionary& dict
55 )
56 {
57  const phaseInterface interface(phase1_, phase2_);
58 
59  const dictionary& interfaceCompositionDict =
60  dict.subDict("interfaceComposition");
61 
62  checkInterfacialModelsDict<sidedInterfaceCompositionModel>
63  (
64  fluid_,
65  interfaceCompositionDict
66  );
67 
68  interfaceCompositionModel_.reset
69  (
71  (
72  interfaceCompositionDict,
73  interface
74  ).ptr()
75  );
76 
77  const dictionary& diffusiveMassTransferDict =
78  dict.subDict("diffusiveMassTransfer");
79 
80  checkBlendedInterfacialModelsDict<blendedSidedDiffusiveMassTransferModel>
81  (
82  fluid_,
83  diffusiveMassTransferDict
84  );
85 
86  diffusiveMassTransferModel_.reset
87  (
89  (
90  diffusiveMassTransferDict,
91  interface,
92  blendingDict<blendedSidedDiffusiveMassTransferModel>
93  (
94  fluid_,
95  diffusiveMassTransferDict
96  )
97  ).ptr()
98  );
99 
100  nIter_ = dict.lookupOrDefault<label>("nIter", 1);
101 }
102 
103 
104 Foam::wordList Foam::fv::massDiffusionLimitedPhaseChange::getSpecies() const
105 {
106  wordHashSet species;
107 
108  forAll(phaseNames(), i)
109  {
110  const phaseModel& phase = i ? phase2_ : phase1_;
111 
112  if (!interfaceCompositionModel_->haveModelInThe(phase)) continue;
113 
114  species.insert(interfaceCompositionModel_->modelInThe(phase).species());
115  }
116 
117  return species.toc();
118 }
119 
120 
121 void Foam::fv::massDiffusionLimitedPhaseChange::correctMDot() const
122 {
124 
125  forAll(phaseNames(), i)
126  {
127  const phaseModel& phase = i ? phase2_ : phase1_;
128  const label s = i ? +1 : -1;
129 
130  if (!interfaceCompositionModel_->haveModelInThe(phase)) continue;
131 
132  const interfaceCompositionModel& phaseIcm =
133  interfaceCompositionModel_->modelInThe(phase);
134 
135  forAll(phaseIcm.species(), phaseIcmSpeciei)
136  {
137  const word& specieName = phaseIcm.species()[phaseIcmSpeciei];
138  const label speciei = phaseIcm.thermo().species()[specieName];
139 
140  mDot_ +=
141  s
142  *(
143  mDotSus_[i][phaseIcmSpeciei]
144  + mDotSps_[i][phaseIcmSpeciei]*phaseIcm.thermo().Y()[speciei]
145  );
146  }
147  }
148 }
149 
150 
151 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
152 
154 (
155  const word& name,
156  const word& modelType,
157  const fvMesh& mesh,
158  const dictionary& dict
159 )
160 :
161  phaseChange(name, modelType, mesh, dict, NullObjectRef<wordList>()),
162  solver_(mesh().lookupObject<solvers::multiphaseEuler>(solver::typeName)),
163  fluid_(solver_.fluid),
164  phase1_(fluid_.phases()[phaseNames().first()]),
165  phase2_(fluid_.phases()[phaseNames().second()]),
166  interfaceCompositionModel_(nullptr),
167  diffusiveMassTransferModel_(nullptr),
168  nIter_(1),
169  Ts_
170  (
171  IOobject
172  (
173  name + ":Ts",
174  mesh.time().name(),
175  mesh,
176  IOobject::READ_IF_PRESENT,
177  IOobject::AUTO_WRITE
178  ),
179  (phase1_.thermo().T()() + phase2_.thermo().T()())/2
180  ),
181  mDotSus_(),
182  mDotSps_(),
183  pressureEquationIndex_(-1),
184  mDot_
185  (
186  IOobject
187  (
188  name + ":mDot",
189  mesh.time().name(),
190  mesh,
191  IOobject::READ_IF_PRESENT,
192  IOobject::AUTO_WRITE
193  ),
194  mesh,
196  )
197 {
198  readCoeffs(coeffs(dict));
199  setSpecies(name, modelType, getSpecies());
200 
201  // Allocate the coefficient fields
202  forAll(phaseNames(), i)
203  {
204  const phaseModel& phase = i ? phase2_ : phase1_;
205 
206  if (!interfaceCompositionModel_->haveModelInThe(phase)) continue;
207 
208  const interfaceCompositionModel& phaseIcm =
209  interfaceCompositionModel_->modelInThe(phase);
210 
211  mDotSus_[i].resize(phaseIcm.species().size());
212 
213  mDotSps_[i].resize(phaseIcm.species().size());
214 
215  forAll(phaseIcm.species(), phaseIcmSpeciei)
216  {
217  const word& specieName = phaseIcm.species()[phaseIcmSpeciei];
218 
219  mDotSus_[i].set
220  (
221  phaseIcmSpeciei,
223  (
224  IOobject
225  (
226  name + ":mDot" + specieName.capitalise() + "Su",
227  mesh.time().name(),
228  mesh
229  ),
230  mesh,
232  )
233  );
234 
235  mDotSps_[i].set
236  (
237  phaseIcmSpeciei,
239  (
240  IOobject
241  (
242  name + ":mDot" + specieName.capitalise() + "Sp",
243  mesh.time().name(),
244  mesh
245  ),
246  mesh,
248  )
249  );
250  }
251  }
252 }
253 
254 
255 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
256 
259 {
260  return Ts_;
261 }
262 
263 
266 {
267  // Heat transfer coefficients
268  const Pair<tmp<volScalarField>> Hs =
269  solver_.heatTransfer.Hs(phase1_, phase2_);
270  const volScalarField::Internal& H1 = Hs.first();
271  const volScalarField::Internal& H2 = Hs.second();
272 
273  return H2/(H1 + H2);
274 }
275 
276 
279 {
280  return mDot_;
281 }
282 
283 
286 {
287  const word& specieName = species()[mDoti];
288 
291  (
292  name() + ":mDot" + specieName.capitalise(),
293  this->mesh(),
295  );
296 
297  forAll(phaseNames(), i)
298  {
299  const phaseModel& phase = i ? phase2_ : phase1_;
300  const label s = i ? +1 : -1;
301 
302  if (!interfaceCompositionModel_->haveModelInThe(phase)) continue;
303 
304  const interfaceCompositionModel& phaseIcm =
305  interfaceCompositionModel_->modelInThe(phase);
306 
307  if (!phaseIcm.species().found(specieName)) continue;
308 
309  const label phaseIcmSpeciei = phaseIcm.species()[specieName];
310  const label speciei = phaseIcm.thermo().species()[specieName];
311 
312  tResult.ref() +=
313  s
314  *(
315  mDotSus_[i][phaseIcmSpeciei]
316  + mDotSps_[i][phaseIcmSpeciei]*phaseIcm.thermo().Y()[speciei]
317  );
318  }
319 
320  return tResult;
321 }
322 
323 
325 (
326  const volScalarField& alpha,
327  const volScalarField& rho,
328  fvMatrix<scalar>& eqn
329 ) const
330 {
331  // Pressure equation (i.e., continuity, linearised in the pressure)
332  if
333  (
334  (&alpha == &phase1_ || &alpha == &phase2_)
335  && (&rho == &phase1_.rho() || &rho == &phase2_.rho())
336  && &eqn.psi() == &solver_.p_rgh
337  )
338  {
339  // Ensure that the source is up-to date if this is the first call in
340  // this outer corrector
341  if (pressureEquationIndex_ == 0) correctMDot();
342  pressureEquationIndex_ ++;
343  }
344 
346 }
347 
348 
350 (
351  const volScalarField& alpha,
352  const volScalarField& rho,
353  const volScalarField& heOrYi,
354  fvMatrix<scalar>& eqn
355 ) const
356 {
357  const label index = this->index(phaseNames(), alpha.group());
358 
359  const ThermoRefPair<multicomponentThermo>& mcThermos =
360  thermos().thermos<multicomponentThermo>();
361 
362  const word specieName = heOrYi.member();
363 
364  // Mass fraction equation
365  if (mcThermos.valid()[index] && mcThermos[index].containsSpecie(specieName))
366  {
367  // A non-transferring specie. Do not add a source.
368  if (!species().found(specieName)) return;
369 
370  // A transferring specie. Add a (potentially) linearised source.
371  forAll(phaseNames(), i)
372  {
373  const phaseModel& phase = i ? phase2_ : phase1_;
374 
375  if (!interfaceCompositionModel_->haveModelInThe(phase)) continue;
376 
377  const interfaceCompositionModel& phaseIcm =
378  interfaceCompositionModel_->modelInThe(phase);
379 
380  if (!phaseIcm.species().found(specieName)) continue;
381 
382  const label phaseIcmSpeciei = phaseIcm.species()[specieName];
383  const label speciei = phaseIcm.thermo().species()[specieName];
384 
385  if (i == index)
386  {
387  eqn +=
388  mDotSus_[i][phaseIcmSpeciei]
389  + fvm::Sp(mDotSps_[i][phaseIcmSpeciei], eqn.psi());
390  }
391  else
392  {
393  eqn -=
394  mDotSus_[i][phaseIcmSpeciei]
395  + mDotSps_[i][phaseIcmSpeciei]*phaseIcm.thermo().Y()[speciei];
396  }
397  }
398 
399  return;
400  }
401 
402  // Something else. Fall back.
403  phaseChange::addSup(alpha, rho, heOrYi, eqn);
404 }
405 
406 
408 {
409  Info<< type() << ": " << name() << endl << incrIndent;
410 
411  const volScalarField::Internal& T1 = phase1_.thermo().T();
412  const volScalarField::Internal& T2 = phase2_.thermo().T();
413 
414  // Heat transfer coefficients
415  const Pair<tmp<volScalarField>> Hs =
416  solver_.heatTransfer.Hs(phase1_, phase2_, scalar(0));
417  const volScalarField::Internal& H1 = Hs.first();
418  const volScalarField::Internal& H2 = Hs.second();
419 
420  // Stabilisation heat transfer coefficient
421  static const dimensionedScalar HSmall
422  (
423  "small",
425  small
426  );
427 
428  // Mass transfer coefficients
430  forAll(phaseNames(), i)
431  {
432  const phaseModel& phase = i ? phase2_ : phase1_;
433 
434  if (!interfaceCompositionModel_->haveModelInThe(phase)) continue;
435 
436  KPtrs[i].set(diffusiveMassTransferModel_->KinThe(phase).ptr());
437  }
438 
439  // Mass diffusivities
441  forAll(phaseNames(), i)
442  {
443  const phaseModel& phase = i ? phase2_ : phase1_;
444 
445  if (!interfaceCompositionModel_->haveModelInThe(phase)) continue;
446 
447  interfaceCompositionModel& phaseIcm =
448  interfaceCompositionModel_->modelInThe(phase);
449 
450  Ds[i].setSize(phaseIcm.species().size());
451 
452  forAll(phaseIcm.species(), phaseIcmSpeciei)
453  {
454  const word& specieName = phaseIcm.species()[phaseIcmSpeciei];
455 
456  Ds[i].set
457  (
458  phaseIcmSpeciei,
459  phaseIcm.D(specieName).ptr()
460  );
461  }
462  }
463 
464  // Iterative solution of the interface heat-mass transfer balance
465  for (label iteri = 0; iteri < nIter_; ++ iteri)
466  {
469  (
470  name() + ":mDotL",
471  this->mesh(),
473  );
474  tmp<volScalarField::Internal> mDotLPrime =
476  (
477  name() + ":mDotLPrime",
478  this->mesh(),
480  );
481 
482  // Add latent heats from forward and backward models
483  forAll(phaseNames(), i)
484  {
485  const phaseModel& phase = i ? phase2_ : phase1_;
486  const label s = i ? +1 : -1;
487 
488  if (!interfaceCompositionModel_->haveModelInThe(phase)) continue;
489 
490  interfaceCompositionModel& phaseIcm =
491  interfaceCompositionModel_->modelInThe(phase);
492 
493  forAll(phaseIcm.species(), phaseIcmSpeciei)
494  {
495  const word& specieName = phaseIcm.species()[phaseIcmSpeciei];
496  const label speciei = phaseIcm.thermo().species()[specieName];
497  const label mDoti = species()[specieName];
498 
499  const volScalarField::Internal Yf
500  (
501  phaseIcm.Yf(specieName, vifToVf(Ts_))
502  );
503  const volScalarField::Internal YfPrime
504  (
505  phaseIcm.YfPrime(specieName, vifToVf(Ts_))
506  );
507 
508  const volScalarField::Internal rhoKDL
509  (
510  phase.rho()()
511  *KPtrs[i]()
512  *Ds[i][phaseIcmSpeciei]
513  *s*L(Ts_, mDoti)
514  );
515 
516  mDotL.ref() += rhoKDL*(Yf - phaseIcm.thermo().Y()[speciei]);
517  mDotLPrime.ref() += rhoKDL*YfPrime;
518  }
519  }
520 
521  // Update the interface temperature by applying one step of newton's
522  // method to the interface relation
523  Ts_ -=
524  (H1*(Ts_ - T1) + H2*(Ts_ - T2) + mDotL)
525  /(max(H1 + H2 + mDotLPrime, HSmall));
526 
527  Info<< indent << "min/mean/max Ts"
528  << " = " << gMin(Ts_.primitiveField())
529  << '/' << gAverage(Ts_.primitiveField())
530  << '/' << gMax(Ts_.primitiveField())
531  << endl;
532 
533  // Update the interface compositions
534  forAll(phaseNames(), i)
535  {
536  const phaseModel& phase = i ? phase2_ : phase1_;
537 
538  if (!interfaceCompositionModel_->haveModelInThe(phase)) continue;
539 
540  interfaceCompositionModel& phaseIcm =
541  interfaceCompositionModel_->modelInThe(phase);
542 
543  phaseIcm.update(vifToVf(Ts_));
544  }
545  }
546 
547  // Build the coefficients of the phase change rate
548  forAll(phaseNames(), i)
549  {
550  const phaseModel& phase = i ? phase2_ : phase1_;
551 
552  if (!interfaceCompositionModel_->haveModelInThe(phase)) continue;
553 
554  interfaceCompositionModel& phaseIcm =
555  interfaceCompositionModel_->modelInThe(phase);
556 
557  forAll(phaseIcm.species(), phaseIcmSpeciei)
558  {
559  const word& specieName = phaseIcm.species()[phaseIcmSpeciei];
560 
561  const volScalarField::Internal KD
562  (
563  KPtrs[i]()*Ds[i][phaseIcmSpeciei]
564  );
565 
566  const volScalarField::Internal Yf
567  (
568  phaseIcm.Yf(specieName, vifToVf(Ts_))
569  );
570 
571  mDotSus_[i][phaseIcmSpeciei] = phase.rho()()*KD*Yf;
572  mDotSps_[i][phaseIcmSpeciei] = -phase.rho()()*KD;
573  }
574  }
575 
576  // Reset the p_rgh equation solution counter
577  pressureEquationIndex_ = 0;
578 
579  // Correct the total phase change rate
580  correctMDot();
581 
582  Info<< indent << "min/mean/max mDot"
583  << " = " << gMin(mDot_.primitiveField())
584  << '/' << gAverage(mDot_.primitiveField())
585  << '/' << gMax(mDot_.primitiveField())
586  << endl;
587 
588  Info<< decrIndent;
589 }
590 
591 
593 {
594  if (phaseChange::read(dict))
595  {
596  readCoeffs(coeffs(dict));
597  reSetSpecies(getSpecies());
598  return true;
599  }
600  else
601  {
602  return false;
603  }
604 }
605 
606 
607 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
tmp< DimensionedField< Type, GeoMesh, Field > > T() const
Return the field transpose (only defined for second rank tensors)
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word member(const word &name)
Return member (name without the extension)
Definition: IOobject.C:146
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An ordered pair of two objects of type <Type> with first() and second() elements.
Definition: Pair.H:66
const Type & second() const
Return second.
Definition: PairI.H:121
const Type & first() const
Return first.
Definition: PairI.H:107
Class containing a pair of thermo references. Handles down-casting to more specific thermo types by c...
Definition: ThermoRefPair.H:51
const Pair< bool > & valid() const
Access the validity flags.
ThermoRefPair< OtherThermoType > thermos() const
Cast to a different thermo type, with error checking.
static autoPtr< blendedSidedDiffusiveMassTransferModel > New(const dictionary &dict, const phaseInterface &interface, const dictionary &blendingDict)
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const word & name() const
Return const reference to name.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
Finite volume model abstract base class.
Definition: fvModel.H:60
static const dictionary & coeffs(const word &modelType, const dictionary &)
Return the coefficients sub-dictionary for a given model type.
Definition: fvModelI.H:31
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:57
Model for mass-diffusion rate limited phase change between two phases.
virtual tmp< DimensionedField< scalar, volMesh > > Lfraction() const
Return the fraction of the latent heat that is transferred into.
virtual tmp< DimensionedField< scalar, volMesh > > mDot() const
Return the mass transfer rate.
virtual bool read(const dictionary &dict)
Read source dictionary.
void addSup(const volScalarField &alpha, const volScalarField &rho, fvMatrix< scalar > &eqn) const
Override the pressure equation to trigger correction of the.
virtual tmp< DimensionedField< scalar, volMesh > > Tchange() const
Return the temperature at which the phases are considered to be.
massDiffusionLimitedPhaseChange(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
const Pair< word > & phaseNames() const
Return the names of the phases.
Definition: massTransferI.H:52
virtual void addSup(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
Definition: massTransfer.C:225
Base class for phase change models.
Definition: phaseChange.H:61
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: phaseChange.C:624
void addSup(const volScalarField &alpha, const volScalarField &rho, const volScalarField &heOrYi, fvMatrix< scalar > &eqn) const
Override the energy equation to add the phase change heat, or.
Definition: phaseChange.C:551
void setSpecies(const word &name, const word &modelType, const wordList &species)
Set the names of the transferring species.
Definition: phaseChange.C:152
bool found(const word &) const
Does the list contain the specified name.
static const dimensionSet dimK
Coefficient dimensions.
Generic base class for interface composition models. These models describe the composition in phase 1...
tmp< volScalarField > D(const word &speciesName) const
Mass diffusivity.
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const =0
Interface mass fraction.
virtual void update(const volScalarField &Tf)=0
Update the composition.
const rhoFluidMulticomponentThermo & thermo() const
Return the thermo.
const hashedWordList & species() const
Return the transferring species names.
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const =0
The interface mass fraction derivative w.r.t. temperature.
Base-class for multi-component thermodynamic properties.
virtual const speciesTable & species() const =0
The table of species.
virtual PtrList< volScalarField > & Y()=0
Access the mass-fraction fields.
virtual const volScalarField & rho() const =0
Return the density field.
static autoPtr< sidedInterfaceCompositionModel > New(const dictionary &dict, const phaseInterface &interface)
Abstract base class for run-time selectable region solvers.
Definition: solver.H:56
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
word capitalise() const
Return the word with the first letter capitalised.
Definition: wordI.H:131
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Calculate the matrix for implicit and explicit sources.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:242
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimEnergy
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
messageStream Info
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
const dimensionSet dimTemperature
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:235
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
const dimensionSet dimTime
const dimensionSet dimDensity
const dimensionSet dimVolume
const T & NullObjectRef()
Return const reference to the nullObject of type T.
Definition: nullObjectI.H:27
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Type gAverage(const FieldField< Field, Type > &f)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:210
Type gMin(const FieldField< Field, Type > &f)
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:228
Type gMax(const FieldField< Field, Type > &f)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31