InterfaceCompositionPhaseChangePhaseSystem.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) 2015-2023 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 
28 #include "heatTransferModel.H"
30 #include "fvmSup.H"
31 
32 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
33 
34 template<class BasePhaseSystem>
37 {
39  (
40  interfaceCompositionModelTable,
41  interfaceCompositionModels_,
42  interfaceCompositionModelIter
43  )
44  {
45  const sidedInterfaceCompositionModel& compositionModel =
46  interfaceCompositionModelIter();
47 
48  const phaseInterface& interface = compositionModel.interface();
49 
50  *dmdtfs_[interface] = Zero;
51 
52  forAllConstIter(phaseInterface, interface, interfaceIter)
53  {
54  const phaseModel& phase = interfaceIter();
55 
56  if (!compositionModel.haveModelInThe(phase)) continue;
57 
59  (
60  hashedWordList,
61  compositionModel.modelInThe(phase).species(),
62  specieIter
63  )
64  {
65  const word& specie = *specieIter;
66 
67  *dmdtfs_[interface] +=
68  (interfaceIter.index() == 0 ? +1 : -1)
69  *(
70  *(*dmidtfSus_[interface])[specie]
71  + *(*dmidtfSps_[interface])[specie]*phase.Y(specie)
72  );
73  }
74  }
75  }
76 }
77 
78 
79 template<class BasePhaseSystem>
82 totalDmidtfs() const
83 {
84  autoPtr<phaseSystem::dmidtfTable> totalDmidtfsPtr
85  (
86  new phaseSystem::dmidtfTable
87  );
88  phaseSystem::dmidtfTable& totalDmidtfs = totalDmidtfsPtr();
89 
91  (
92  interfaceCompositionModelTable,
93  interfaceCompositionModels_,
94  interfaceCompositionModelIter
95  )
96  {
97  const sidedInterfaceCompositionModel& compositionModel =
98  interfaceCompositionModelIter();
99 
100  const phaseInterface& interface = compositionModel.interface();
101 
102  if (!totalDmidtfs.found(interface))
103  {
104  totalDmidtfs.insert(interface, new HashPtrTable<volScalarField>());
105  }
106 
107  forAllConstIter(phaseInterface, interface, interfaceIter)
108  {
109  const phaseModel& phase = interfaceIter();
110 
111  if (!compositionModel.haveModelInThe(phase)) continue;
112 
114  (
115  hashedWordList,
116  compositionModel.modelInThe(phase).species(),
117  specieIter
118  )
119  {
120  const word& specie = *specieIter;
121 
122  tmp<volScalarField> dmidtf
123  (
124  (interfaceIter.index() == 0 ? +1 : -1)
125  *(
126  *(*dmidtfSus_[interface])[specie]
127  + *(*dmidtfSps_[interface])[specie]*phase.Y(specie)
128  )
129  );
130 
131  if (totalDmidtfs[interface]->found(specie))
132  {
133  *(*totalDmidtfs[interface])[specie] += dmidtf;
134  }
135  else
136  {
137  totalDmidtfs[interface]->insert(specie, dmidtf.ptr());
138  }
139  }
140  }
141  }
142 
143  return totalDmidtfsPtr;
144 }
145 
146 
147 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
148 
149 template<class BasePhaseSystem>
152 (
153  const fvMesh& mesh
154 )
155 :
156  BasePhaseSystem(mesh),
157  nInterfaceCorrectors_
158  (
159  this->template lookupOrDefault<label>("nInterfaceCorrectors", 1)
160  )
161 {
162  this->generateInterfacialModels(interfaceCompositionModels_);
163  this->generateInterfacialModels(diffusiveMassTransferModels_);
164 
165  // Check that models have been specified in the correct combinations
167  (
169  interfaceCompositionModels_,
170  interfaceCompositionModelIter
171  )
172  {
173  const sidedInterfaceCompositionModel& sidedCompositionModel =
174  interfaceCompositionModelIter();
175 
176  const phaseInterface& interface = sidedCompositionModel.interface();
177  const phaseModel& phase1 = interface.phase1();
178  const phaseModel& phase2 = interface.phase2();
179 
180  this->template
181  validateMassTransfer<interfaceCompositionModel>(interface);
182 
183  if (!this->diffusiveMassTransferModels_.found(interface))
184  {
186  << "A diffusive mass transfer model for the " << interface
187  << " interface is not specified. This is required by the "
188  << "corresponding interface composition model."
189  << exit(FatalError);
190  }
191 
192  forAllConstIter(phaseInterface, interface, interfaceIter)
193  {
194  if
195  (
196  sidedCompositionModel.haveModelInThe(interfaceIter())
197  && !diffusiveMassTransferModels_[interface]
198  ->haveModelInThe(interfaceIter())
199  )
200  {
202  << "A diffusive mass transfer model for the "
203  << interfaceIter().name() << " side of the "
204  << interface.name() << " interface is not "
205  << "specified. This is required by the corresponding "
206  << "interface composition model."
207  << exit(FatalError);
208  }
209  }
210 
211  if
212  (
213  !this->heatTransferModels_.found(interface)
214  || !this->heatTransferModels_[interface]->haveModelInThe(phase1)
215  || !this->heatTransferModels_[interface]->haveModelInThe(phase2)
216  )
217  {
219  << "A heat transfer model for both sides of the "
220  << interface.name()
221  << " interface is not specified. This is required by the "
222  << "corresponding interface composition model"
223  << exit(FatalError);
224  }
225  }
226 
227  // Generate mass transfer fields, initially assumed to be zero
229  (
231  interfaceCompositionModels_,
232  interfaceCompositionModelIter
233  )
234  {
235  const sidedInterfaceCompositionModel& sidedCompositionModel =
236  interfaceCompositionModelIter();
237 
238  const phaseInterface& interface = sidedCompositionModel.interface();
239 
240  dmdtfs_.insert
241  (
242  interface,
243  new volScalarField
244  (
245  IOobject
246  (
248  (
249  "interfaceCompositionPhaseChange:dmdtf",
250  interface.name()
251  ),
252  this->mesh().time().name(),
253  this->mesh(),
256  ),
257  this->mesh(),
259  )
260  );
261 
262  dmidtfSus_.insert(interface, new HashPtrTable<volScalarField>());
263 
264  dmidtfSps_.insert(interface, new HashPtrTable<volScalarField>());
265 
266  Tfs_.insert
267  (
268  interface,
269  new volScalarField
270  (
271  IOobject
272  (
274  (
275  "interfaceCompositionPhaseChange:Tf",
276  interface.name()
277  ),
278  this->mesh().time().name(),
279  this->mesh(),
282  ),
283  (
284  interface.phase1().thermo().T()
285  + interface.phase2().thermo().T()
286  )/2
287  )
288  );
289 
290  forAllConstIter(phaseInterface, interface, interfaceIter)
291  {
292  const phaseModel& phase = interfaceIter();
293 
294  if (!sidedCompositionModel.haveModelInThe(phase)) continue;
295 
296  const interfaceCompositionModel& compositionModel =
297  sidedCompositionModel.modelInThe(phase);
298 
300  (
302  compositionModel.species(),
303  specieIter
304  )
305  {
306  const word& specie = *specieIter;
307 
308  dmidtfSus_[interface]->insert
309  (
310  specie,
311  new volScalarField
312  (
313  IOobject
314  (
316  (
318  (
319  "interfaceCompositionPhaseChange:dmidtfSu",
320  specie
321  ),
322  interface.name()
323  ),
324  this->mesh().time().name(),
325  this->mesh()
326  ),
327  this->mesh(),
329  )
330  );
331 
332  dmidtfSps_[interface]->insert
333  (
334  specie,
335  new volScalarField
336  (
337  IOobject
338  (
340  (
342  (
343  "interfaceCompositionPhaseChange:dmidtfSp",
344  specie
345  ),
346  interface.name()
347  ),
348  this->mesh().time().name(),
349  this->mesh()
350  ),
351  this->mesh(),
353  )
354  );
355  }
356  }
357  }
358 }
359 
360 
361 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
362 
363 template<class BasePhaseSystem>
366 {}
367 
368 
369 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
370 
371 template<class BasePhaseSystem>
374 (
375  const phaseInterfaceKey& key
376 ) const
377 {
378  tmp<volScalarField> tDmdtf = BasePhaseSystem::dmdtf(key);
379 
380  if (dmdtfs_.found(key))
381  {
382  tDmdtf.ref() += *dmdtfs_[key];
383  }
384 
385  return tDmdtf;
386 }
387 
388 
389 template<class BasePhaseSystem>
392 {
393  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
394 
395  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs_, dmdtfIter)
396  {
397  const phaseInterface interface(*this, dmdtfIter.key());
398 
399  addField(interface.phase1(), "dmdt", *dmdtfIter(), dmdts);
400  addField(interface.phase2(), "dmdt", - *dmdtfIter(), dmdts);
401  }
402 
403  return dmdts;
404 }
405 
406 
407 template<class BasePhaseSystem>
411 {
413  BasePhaseSystem::momentumTransfer();
414 
415  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
416 
417  this->addDmdtUfs(dmdtfs_, eqns);
418 
419  return eqnsPtr;
420 }
421 
422 
423 template<class BasePhaseSystem>
427 {
429  BasePhaseSystem::momentumTransferf();
430 
431  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
432 
433  this->addDmdtUfs(dmdtfs_, eqns);
434 
435  return eqnsPtr;
436 }
437 
438 
439 template<class BasePhaseSystem>
442 heatTransfer() const
443 {
445  BasePhaseSystem::heatTransfer();
446 
447  phaseSystem::heatTransferTable& eqns = eqnsPtr();
448 
449  this->addDmidtHefs
450  (
451  totalDmidtfs(),
452  Tfs_,
453  latentHeatScheme::symmetric,
454  latentHeatTransfer::mass,
455  eqns
456  );
457 
458  return eqnsPtr;
459 }
460 
461 
462 template<class BasePhaseSystem>
465 specieTransfer() const
466 {
468  BasePhaseSystem::specieTransfer();
469 
470  phaseSystem::specieTransferTable& eqns = eqnsPtr();
471 
472  // Explicit
473  /*
474  this->addDmidtYf(totalDmidtfs(), eqns);
475  */
476 
477  // Semi-implicit
479  (
481  interfaceCompositionModels_,
482  interfaceCompositionModelIter
483  )
484  {
485  const sidedInterfaceCompositionModel& compositionModel =
486  interfaceCompositionModelIter();
487 
488  const phaseInterface& interface = compositionModel.interface();
489 
490  forAllConstIter(phaseInterface, interface, interfaceIter)
491  {
492  const phaseModel& phase = interfaceIter();
493  const phaseModel& otherPhase = interfaceIter.otherPhase();
494 
495  if (!compositionModel.haveModelInThe(phase)) continue;
496 
498  (
500  compositionModel.modelInThe(phase).species(),
501  specieIter
502  )
503  {
504  const word& specie = *specieIter;
505 
506  // Implicit transport through this phase
507  *eqns[phase.Y(specie).name()] +=
508  *(*dmidtfSus_[interface])[specie]
509  + fvm::Sp(*(*dmidtfSps_[interface])[specie], phase.Y(specie));
510 
511  // Explicit transport out of the other phase
512  if (eqns.found(IOobject::groupName(specie, otherPhase.name())))
513  {
514  *eqns[otherPhase.Y(specie).name()] -=
515  *(*dmidtfSus_[interface])[specie]
516  + *(*dmidtfSps_[interface])[specie]*phase.Y(specie);
517  }
518  }
519  }
520  }
521 
522  return eqnsPtr;
523 }
524 
525 
526 template<class BasePhaseSystem>
528 correct()
529 {
531 
532  // Sum up the contribution from each interface composition model
534  (
536  interfaceCompositionModels_,
537  interfaceCompositionModelIter
538  )
539  {
540  const sidedInterfaceCompositionModel& compositionModel =
541  interfaceCompositionModelIter();
542 
543  const phaseInterface& interface = compositionModel.interface();
544 
545  const volScalarField& Tf(*this->Tfs_[interface]);
546 
547  forAllConstIter(phaseInterface, interface, interfaceIter)
548  {
549  const phaseModel& phase = interfaceIter();
550 
551  if (!compositionModel.haveModelInThe(phase)) continue;
552 
553  const volScalarField K
554  (
555  diffusiveMassTransferModels_[interface]->modelInThe(phase).K()
556  );
557 
559  (
561  compositionModel.modelInThe(phase).species(),
562  specieIter
563  )
564  {
565  const word& specie = *specieIter;
566 
567  const volScalarField KD
568  (
569  K*compositionModel.modelInThe(phase).D(specie)
570  );
571  const volScalarField Yf
572  (
573  compositionModel.modelInThe(phase).Yf(specie, Tf)
574  );
575 
576  *(*dmidtfSus_[interface])[specie] = phase.rho()*KD*Yf;
577  *(*dmidtfSps_[interface])[specie] = - phase.rho()*KD;
578  }
579  }
580  }
581 
582  correctDmdtfs();
583 }
584 
585 
586 template<class BasePhaseSystem>
589 {
590  BasePhaseSystem::correctSpecies();
591 
592  correctDmdtfs();
593 }
594 
595 
596 template<class BasePhaseSystem>
599 {
600  // This loop solves for the interface temperatures, Tf, and updates the
601  // interface composition models.
602  //
603  // In the presence of thermally-coupled mass transfer, the relation between
604  // heat transfers across the interface between phases 1 and 2 is:
605  //
606  // Q1 + Q2 = mDot*L
607  // H1*(Tf - T1) + H2*(Tf - T1) = K*rho*(Yfi - Yi)*Li
608  //
609  // Where Q1 and Q2 are the net transfer into phases 1 and 2 respectively,
610  // H1 and H2 are the heat transfer coefficients on either side, Tf is the
611  // temperature at the interface, mDot is the mass transfer rate from phase
612  // 2 to phase 1, and L is the latent heat of phase 2 minus phase 1, K is
613  // the diffusive mass transfer coefficient, Yfi - Yi is the concentration
614  // difference of a transferring specie between the interface and the bulk
615  // driving the transfer, Li is the latent heat change of the specie, and
616  // rho is the density in the phase in which the diffusive mass transfer is
617  // being represented.
618  //
619  // Yfi is likely to be a strong non-linear (typically exponential) function
620  // of Tf, so the solution for the temperature is newton-accelerated.
621 
622  forAllIter
623  (
625  interfaceCompositionModels_,
626  interfaceCompositionModelIter
627  )
628  {
629  sidedInterfaceCompositionModel& compositionModel =
630  *interfaceCompositionModelIter();
631 
632  const phaseInterface& interface = compositionModel.interface();
633  const phaseModel& phase1 = interface.phase1();
634  const phaseModel& phase2 = interface.phase2();
635 
637  this->heatTransferModels_[interface];
638 
641  diffusiveMassTransferModels_[interface];
642 
643  const volScalarField H1(heatTransferModel.modelInThe(phase1).K());
644  const volScalarField H2(heatTransferModel.modelInThe(phase2).K());
645  const dimensionedScalar HSmall("small", heatTransferModel::dimK, small);
646 
647  volScalarField& Tf = *this->Tfs_[interface];
648 
649  for (label i = 0; i < nInterfaceCorrectors_; ++ i)
650  {
651  tmp<volScalarField> dmdtLf =
653  (
654  IOobject::groupName("dmdtLf", interface.name()),
655  this->mesh(),
657  );
658  tmp<volScalarField> dmdtLfPrime =
660  (
661  IOobject::groupName("dmdtLfPrime", interface.name()),
662  this->mesh(),
663  dimensionedScalar(dmdtLf().dimensions()/dimTemperature, 0)
664  );
665 
666  // Add latent heats from forward and backward models
667  forAllConstIter(phaseInterface, interface, interfaceIter)
668  {
669  const phaseModel& phase = interfaceIter();
670 
671  if (!compositionModel.haveModelInThe(phase)) continue;
672 
673  const label sign = interfaceIter.index() == 0 ? 1 : -1;
674 
676  (
678  compositionModel.modelInThe(phase).species(),
679  specieIter
680  )
681  {
682  const word& specie = *specieIter;
683 
684  const volScalarField dY
685  (
686  compositionModel.modelInThe(phase).dY(specie, Tf)
687  );
688 
689  const volScalarField dYfPrime
690  (
691  compositionModel.modelInThe(phase).dYfPrime(specie, Tf)
692  );
693 
694  const volScalarField rhoKDL
695  (
696  phase.rho()
697  *diffusiveMassTransferModel.modelInThe(phase).K()
698  *compositionModel.modelInThe(phase).D(specie)
699  *this->Li
700  (
701  interface,
702  specie,
703  dY,
704  Tf,
705  latentHeatScheme::symmetric
706  )
707  );
708 
709  dmdtLf.ref() += sign*rhoKDL*dY;
710  dmdtLfPrime.ref() += sign*rhoKDL*dYfPrime;
711  }
712  }
713 
714  // Update the interface temperature by applying one step of newton's
715  // method to the interface relation
716  Tf -=
717  (
718  H1*(Tf - interface.phase1().thermo().T())
719  + H2*(Tf - interface.phase2().thermo().T())
720  - dmdtLf
721  )
722  /(
723  max(H1 + H2 - dmdtLfPrime, HSmall)
724  );
725 
727 
728  Info<< "Tf." << interface.name()
729  << ": min = " << min(Tf.primitiveField())
730  << ", mean = " << average(Tf.primitiveField())
731  << ", max = " << max(Tf.primitiveField())
732  << endl;
733 
734  // Update the interface compositions
735  forAllConstIter(phaseInterface, interface, interfaceIter)
736  {
737  const phaseModel& phase = interfaceIter();
738 
739  if (!compositionModel.haveModelInThe(phase)) continue;
740 
741  compositionModel.modelInThe(phase).update(Tf);
742  }
743  }
744  }
745 }
746 
747 
748 template<class BasePhaseSystem>
750 {
751  if (BasePhaseSystem::read())
752  {
753  bool readOK = true;
754 
755  // Models ...
756 
757  return readOK;
758  }
759  else
760  {
761  return false;
762  }
763 }
764 
765 
766 // ************************************************************************* //
bool found
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Generic GeometricField class.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
void correctBoundaryConditions()
Correct boundary field.
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:67
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word groupName(Name name, const word &group)
Class to provide interfacial heat and mass transfer between a number of phases according to a interfa...
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
virtual void correctSpecies()
Correct the mass transfer rates for the new species mass fractions.
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
virtual void correct()
Correct the fluid properties other than those listed below.
virtual void correctInterfaceThermo()
Correct the interface temperatures.
virtual tmp< volScalarField > dmdtf(const phaseInterfaceKey &key) const
Return the mass transfer rate for an interface.
virtual autoPtr< phaseSystem::specieTransferTable > specieTransfer() const
Return the specie transfer matrices.
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
Return the momentum transfer matrices for the face-based algorithm.
virtual bool read()
Read base phaseProperties dictionary.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
const phaseInterface & interface() const
Access the interface.
bool haveModelInThe(const phaseModel &phase) const
Does a model exist in the given phase?
const ModelType & modelInThe(const phaseModel &phase) const
Access the model within the given phase.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Model for diffusive mass transfer coefficients between two phases.
virtual tmp< volScalarField > K() const =0
The implicit mass transfer coefficient.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
A wordList with hashed indices for faster lookup by name.
Model for heat transfer between phases.
tmp< volScalarField > K() const
The heat transfer function K used in the enthalpy equation.
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.
tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
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.
tmp< volScalarField > dYfPrime(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
const hashedWordList & species() const
Return the transferring species names.
Word-pair based class used for keying interface models in hash tables.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
const phaseModel & phase1() const
Return phase 1.
const phaseModel & phase2() const
Return phase 2.
virtual const PtrList< volScalarField > & Y() const =0
Return the species mass fractions.
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:109
virtual const volScalarField & rho() const =0
Return the density field.
Base class of the thermophysical property types.
Definition: specie.H:66
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:181
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Calculate the matrix for implicit and explicit sources.
K
Definition: pEqn.H:75
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
void addField(const Group &group, const word &name, tmp< GeoField > field, PtrList< GeoField > &fieldList)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
dimensionedScalar sign(const dimensionedScalar &ds)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
const dimensionSet dimTemperature
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
const dimensionSet dimDensity
const dimensionSet dimVolume
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError