multiphaseInterSystem.H
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) 2017-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::multiphaseInterSystem
28
29Description
30
31SourceFiles
32 multiphaseInterSystem.C
33
34\*---------------------------------------------------------------------------*/
35
36#ifndef multiphaseInterSystem_H
37#define multiphaseInterSystem_H
38
39#include "basicThermo.H"
40
41#include "phaseModel.H"
42#include "phasePair.H"
43#include "orderedPhasePair.H"
44
45#include "volFields.H"
46#include "surfaceFields.H"
47#include "fvMatricesFwd.H"
49#include "localMin.H"
50
52
53// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54
55namespace Foam
56{
57
58// Forward Declarations
59namespace multiphaseInter
60{
61class surfaceTensionModel;
62}
63class porousModel;
64
65/*---------------------------------------------------------------------------*\
66 Class multiphaseInterSystem Declaration
67\*---------------------------------------------------------------------------*/
70:
71 public basicThermo,
73{
74public:
75
76 // Public typedefs
77
78 typedef
80 <
82 >
84
87
90
91
92protected:
93
94 // Protected typedefs
95
96 typedef
98
99
100 typedef
102 <
108
109
110 typedef
112 <
118
119
120
121 // Protected data
122
123 //- Reference to the mesh
124 const fvMesh& mesh_;
125
126 //- Dynamic viscocity
128
129 //- Phase names
131
132 //- Mixture total volumetric flux
134
135 //- Mixture total mass flux
137
138 //- Phase models
140
141 //- Phase pairs
143
144 //- Total ordered phase pairs in the system
146
147 //- Turbulent Prandt number
149
150 //- Turbulence model
152
153
154 // Sub Models
155
156 //- Surface tension models
158
159 //- Interface porous models
161
162
163 // Protected member functions
164
165
166 //- Calculate and return the laminar viscosity
167 void calcMu();
168
169 //- Generate the phases
171 (
172 const wordList& names
173 ) const;
174
175 //- Generate the mixture flux
177 (
179 ) const;
180
181 //- Generate pairs
182 void generatePairs(const dictTable& modelDicts);
183
184 //- Generate pair table
185 void generatePairsTable();
186
187 //- Generate pairs and sub-model tables using pair keys
188 template<class modelType>
189 void createSubModels
190 (
191 const dictTable& modelDicts,
193 <
197 >& models
198 );
199
200 //- Generate pairs and sub-model tables using mesh
201 template<class modelType>
202 void createSubModels
203 (
204 const dictTable& modelDicts,
205 const fvMesh& mesh,
207 <
211 >& models
212 );
213
214 //- Generate pairs and sub-model tables
215 template<class modelType>
217 (
218 const word& modelName,
220 <
224 >& models
225 );
226
227
228 //- Generate pairs and per-phase sub-model tables with mesh ref
229 template<class modelType>
231 (
232 const word& modelName,
233 const fvMesh& mesh,
235 <
239 >& models
240 );
241
242
243 //- Generate pairs and per-phase sub-model tables
244 template<class modelType>
246 (
247 const word& modelName,
249 <
253 >& models
254 );
255
256
257public:
258
259 //- Runtime type information
260 TypeName("multiphaseInterSystem");
261
262 //- Default name of the phase properties dictionary
263 static const word phasePropertiesName;
264
265
266 // Constructors
267
268 //- Construct from fvMesh
270
271
272 //- Destructor
273 virtual ~multiphaseInterSystem();
274
275
276 // Member Functions
277
278 // Energy related thermo functionaliy functions
279
280 //- Return access to the internal energy field [J/Kg]
281 // \note this mixture thermo is prepared to work with T
282 virtual volScalarField& he()
283 {
285 return const_cast<volScalarField&>(volScalarField::null());
286 }
287
288 //- Return access to the internal energy field [J/Kg]
289 // \note this mixture thermo is prepared to work with T
290 virtual const volScalarField& he() const
291 {
293 return volScalarField::null();
294 }
295
296 //- Enthalpy/Internal energy
297 // for given pressure and temperature [J/kg]
298 virtual tmp<volScalarField> he
299 (
300 const volScalarField& p,
301 const volScalarField& T
302 ) const;
303
304 //- Enthalpy/Internal energy for cell-set [J/kg]
305 virtual tmp<scalarField> he
306 (
307 const scalarField& p,
308 const scalarField& T,
309 const labelList& cells
310 ) const;
311
312 //- Enthalpy/Internal energy for patch [J/kg]
313 virtual tmp<scalarField> he
314 (
315 const scalarField& p,
316 const scalarField& T,
317 const label patchi
318 ) const;
319
320 //- Chemical enthalpy of the mixture [J/kg]
321 virtual tmp<volScalarField> hc() const;
322
323 //- Temperature from enthalpy/internal energy for cell-set
324 virtual tmp<scalarField> THE
325 (
326 const scalarField& h,
327 const scalarField& p,
328 const scalarField& T0,
329 const labelList& cells
330 ) const;
331
332 //- Temperature from enthalpy/internal energy for patch
333 virtual tmp<scalarField> THE
334 (
335 const scalarField& h,
336 const scalarField& p,
337 const scalarField& T0,
338 const label patchi
339 ) const;
340
341
342 // Thermo
343
344 //- Return the mixture density
345 virtual tmp<volScalarField> rho() const;
346
347 //- Return the mixture density on a patch
348 virtual tmp<scalarField> rho(const label patchi) const;
349
350 //- Return Cp of the mixture
351 virtual tmp<volScalarField> Cp() const;
352
353 //- Heat capacity at constant pressure for patch [J/kg/K]
354 virtual tmp<scalarField> Cp
355 (
356 const scalarField& p,
357 const scalarField& T,
358 const label patchi
359 ) const;
360
361 //- Heat capacity using pressure and temperature
362 virtual tmp<scalarField> Cp
363 (
364 const scalarField& p,
365 const scalarField& T,
366 const labelList& cells
367 ) const
368 {
370 return tmp<scalarField>::New(p);
371 }
372
373 //- Return Cv of the mixture
374 virtual tmp<volScalarField> Cv() const;
375
376 //- Heat capacity at constant volume for patch [J/kg/K]
377 virtual tmp<scalarField> Cv
378 (
379 const scalarField& p,
380 const scalarField& T,
381 const label patchI
382 ) const;
383
384 //- Density from pressure and temperature
386 (
387 const scalarField& p,
388 const scalarField& T,
389 const labelList& cells
390 ) const;
391
392 //- Gamma = Cp/Cv []
393 virtual tmp<volScalarField> gamma() const;
394
395 //- Gamma = Cp/Cv for patch []
396 virtual tmp<scalarField> gamma
397 (
398 const scalarField& p,
399 const scalarField& T,
400 const label patchi
401 ) const;
402
403 //- Heat capacity at constant pressure/volume [J/kg/K]
404 virtual tmp<volScalarField> Cpv() const;
405
406 //- Heat capacity at constant pressure/volume for patch [J/kg/K]
407 virtual tmp<scalarField> Cpv
408 (
409 const scalarField& p,
410 const scalarField& T,
411 const label patchi
412 ) const;
413
414 //- Heat capacity ratio []
415 virtual tmp<volScalarField> CpByCpv() const;
416
417 //- Heat capacity ratio for patch []
419 (
420 const scalarField& p,
421 const scalarField& T,
422 const label patchi
423 ) const;
424
425 //- Molecular weight [kg/kmol] of the mixture
426 virtual tmp<volScalarField> W() const;
427
428
429 // Transport
430
431 //- Thermal diffusivity for temperature of mixture [J/m/s/K]
432 virtual tmp<volScalarField> kappa() const;
433
434 //- Thermal diffusivity for temperature
435 // of mixture for patch [J/m/s/K]
436 virtual tmp<scalarField> kappa
437 (
438 const label patchi
439 ) const;
440
441 //- Thermal diffusivity for energy of mixture [kg/m/s]
442 virtual tmp<volScalarField> alphahe() const;
443
444 //- Thermal diffusivity for energy of mixture for patch [kg/m/s]
445 virtual tmp<scalarField> alphahe(const label patchi) const;
446
447 //- Effective thermal diffusivity for temperature
448 // of mixture [J/m/s/K]
450 (
451 const volScalarField& kappat
452 ) const;
453
454 //- Effective thermal diffusivity for temperature
455 // of mixture for patch [J/m/s/K]
457 (
458 const scalarField& alphat,
459 const label patchi
460 ) const;
461
462 //- Effective thermal diffusivity of mixture [kg/m/s]
464 (
465 const volScalarField& alphat
466 ) const;
467
468 //- Effective thermal diffusivity of mixture for patch [kg/m/s]
470 (
471 const scalarField& alphat,
472 const label patchi
473 ) const;
474
475 //- Return Prandt number
476 const dimensionedScalar& Prt() const;
477
478
479 // Access to transport state variables
480
481 //- Dynamic viscosity of mixture [kg/m/s]
482 virtual tmp<volScalarField> mu() const;
483
484 //- Dynamic viscosity of mixture for patch [kg/m/s]
485 virtual tmp<scalarField> mu(const label patchi) const;
486
487 //- Kinematic viscosity of mixture [m^2/s]
488 virtual tmp<volScalarField> nu() const;
489
490 //- Kinematic viscosity of mixture for patch [m^2/s]
491 virtual tmp<scalarField> nu(const label patchi) const;
492
493
494 // Turbulence
495
496 //- Set turbulence model
498 {
499 turb_ = &turb;
500 }
501
502 //- Return pointer to turbulence model
504 {
505 return turb_;
506 }
507
508 //- Return the turbulent dynamic viscosity
509 tmp<volScalarField> mut() const;
510
511 //- Return the effective dynamic viscosity
513
514 //- Return the turbulent kinematic viscosity
515 tmp<volScalarField> nut() const;
516
517 //- Return the effective kinematic viscosity
519
520 //- Effective thermal turbulent diffusivity for temperature
521 // of mixture [J/m/s/K]
523
524 //- Effective thermal turbulent diffusivity for temperature
525 // of mixture for patch [J/m/s/K]
526 tmp<scalarField> kappaEff(const label patchi) const;
527
528 //- Effective thermal turbulent diffusivity of mixture [kg/m/s]
530
531 //- Effective thermal turbulent diffusivity of mixture [kg/m/s]
532 tmp<scalarField> alphaEff(const label patchi) const;
533
534
535 // Phase fluxes
536
537 //- Constant access to the total flux
538 const surfaceScalarField& phi() const;
539
540 //- Access to the total mixture flux
542
543 //- Constant access to the mixture mass flux
544 const surfaceScalarField& rhoPhi() const;
545
546 //- Access to the total mixture mass flux
548
549 //- Mixture U
550 tmp<volVectorField> U() const;
551
552
553 // Surface tension
554
555 //- Calculate surface tension of the mixture
557
558 //- Return the surface tension coefficient
560 (
561 const phasePairKey& key
562 ) const;
563
564
565 //- Return coefficients (1/rho)
566 virtual tmp<volScalarField> coeffs(const word& key) const;
567
568
569 // Interface porous between solid/fluid phases
570
571 //- Add interface porosity on phasePair
573
574
575 // Inter-Phase mass and heat transfe
576
577 //- Return interfacial source mass rate per phase pair
578 virtual tmp<volScalarField> dmdt(const phasePairKey& key) const = 0;
579
580 //- Return the heat transfer matrices
582 (
583 const volScalarField& T
584 ) = 0;
585
586 //- Return the volumetric rate transfer matrix
588 (
589 const volScalarField& p
590 ) = 0;
591
592 //- Calculate mass transfer for alpha's
593 virtual void alphaTransfer(SuSpTable& Su, SuSpTable& Sp) = 0;
594
595 //- Calculate mass transfer for species
596 virtual void massSpeciesTransfer
597 (
601 const word speciesName
602 ) = 0;
603
604 //- Add volume change in pEq
605 virtual bool includeVolChange() = 0;
606
607
608 // Solve phases and correct models
609
610 //- Solve for the phase transport equations
611 virtual void solve() = 0;
612
613 //- Correct the mixture thermos
614 virtual void correct();
615
616 //- Correct mass sources
617 virtual void correctMassSources(const volScalarField& T) = 0;
618
619 //- Return the name of the thermo physics
620 virtual word thermoName() const
621 {
623 return word();
624 }
625
626 //- Correct the turbulence
627 // \note Each phase could have its own turbulence
628 virtual void correctTurbulence();
629
630 //- Read base phaseProperties dictionary
631 virtual bool read();
632
633
634 // Access to phases models
635
636 //- Constant access the total phase pairs
637 const phasePairTable& totalPhasePairs() const;
638
639 //- Non-constant access the total phase pairs
641
642 //- Constant access the phases
643 const phaseModelTable& phases() const;
644
645 //- Access the phases
647
648 //- Access a sub model between a phase pair
649 template <class modelType>
650 const modelType& lookupSubModel(const phasePair& key) const;
651
652 //- Access a sub model between two phases
653 template <class modelType>
654 const modelType& lookupSubModel
655 (
656 const multiphaseInter::phaseModel& from,
658 ) const;
659
660
661 // Query phase thermo information
662
663 //- Return true if the equation of state is incompressible for all
664 // phases
665 virtual bool incompressible() const;
666
667 //- Return true if a phase is incompressible
668 virtual bool incompressible(const word) const;
669
670 //- Return true if the equation of state is isochoric for all phasses
671 // i.e. rho = const
672 virtual bool isochoric() const;
673
674 //- Return mesh
675 const fvMesh& mesh() const;
676
677
678 // Help functions for the interfaces
679
680 //- Interface normal surface vector
682 (
683 const volScalarField& alpha1,
685 ) const;
686
687 //- Interface normal volField vector
689 (
690 const volScalarField& alpha1,
692 ) const;
693
694 //- Interface normal surface vector
696 (
697 const volScalarField& alpha1,
699 ) const;
700
701 //- Interface curvature
703 (
704 const volScalarField& alpha1,
706 ) const;
707
708 //- Near Interface of alpha1 and alpha2
710 (
711 const volScalarField& alpha1,
713 ) const;
714
715 //- Near Interface of alpha'n
717};
718
719
720// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
721
722} // End namespace Foam
723
724// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
725
726#ifdef NoRepository
728#endif
729
730// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
731
732#endif
733
734// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
compressible::turbulenceModel & turb
const volScalarField & alpha1
const volScalarField & alpha2
static const GeometricField< scalar, fvPatchField, volMesh > & null()
Return a null geometric field.
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:66
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:622
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:610
Base-class for all transport models used by the compressible turbulence models.
Abstract base class for turbulence models (RAS, LES and laminar).
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
virtual tmp< fvScalarMatrix > volTransfer(const volScalarField &p)=0
Return the volumetric rate transfer matrix.
void createSubModels(const dictTable &modelDicts, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables using pair keys.
void generatePairsTable()
Generate pair table.
tmp< surfaceScalarField > surfaceTensionForce() const
Calculate surface tension of the mixture.
tmp< volScalarField > nearInterface() const
Near Interface of alpha'n.
virtual ~multiphaseInterSystem()
Destructor.
compressibleTurbulenceModel * turb_
Turbulence model.
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol] of the mixture.
dimensionedScalar Prt_
Turbulent Prandt number.
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const =0
Return interfacial source mass rate per phase pair.
static const word phasePropertiesName
Default name of the phase properties dictionary.
const dimensionedScalar & Prt() const
Return Prandt number.
virtual void alphaTransfer(SuSpTable &Su, SuSpTable &Sp)=0
Calculate mass transfer for alpha's.
virtual tmp< volScalarField > Cp() const
Return Cp of the mixture.
void generatePairs(const dictTable &modelDicts)
Generate pairs.
tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
virtual word thermoName() const
Return the name of the thermo physics.
const fvMesh & mesh_
Reference to the mesh.
virtual tmp< volScalarField > surfaceTensionCoeff(const phasePairKey &key) const
Return the surface tension coefficient.
HashTable< autoPtr< multiphaseInter::surfaceTensionModel >, phasePairKey, phasePairKey::hash > surfaceTensionModelTable
const surfaceScalarField & phi() const
Constant access to the total flux.
virtual void correct()
Correct the mixture thermos.
tmp< volVectorField > U() const
Mixture U.
const phaseModelTable & phases() const
Constant access the phases.
TypeName("multiphaseInterSystem")
Runtime type information.
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
virtual tmp< volScalarField > mu() const
Dynamic viscosity of mixture [kg/m/s].
compressibleTurbulenceModel * turbulence() const
Return pointer to turbulence model.
const modelType & lookupSubModel(const phasePair &key) const
Access a sub model between a phase pair.
tmp< volScalarField > mut() const
Return the turbulent dynamic viscosity.
interfacePorousModelTable interfacePorousModelTable_
Interface porous models.
virtual tmp< fvScalarMatrix > heatTransfer(const volScalarField &T)=0
Return the heat transfer matrices.
virtual tmp< scalarField > Cp(const scalarField &p, const scalarField &T, const labelList &cells) const
Heat capacity using pressure and temperature.
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
virtual void correctMassSources(const volScalarField &T)=0
Correct mass sources.
const surfaceScalarField & rhoPhi() const
Constant access to the mixture mass flux.
surfaceTensionModelTable surfaceTensionModels_
Surface tension models.
const phasePairTable & totalPhasePairs() const
Constant access the total phase pairs.
virtual bool includeVolChange()=0
Add volume change in pEq.
virtual void massSpeciesTransfer(const multiphaseInter::phaseModel &phase, volScalarField::Internal &Su, volScalarField::Internal &Sp, const word speciesName)=0
Calculate mass transfer for species.
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
virtual void correctTurbulence()
Correct the turbulence.
virtual tmp< volScalarField > coeffs(const word &key) const
Return coefficients (1/rho)
tmp< volScalarField > alphaEff() const
Effective thermal turbulent diffusivity of mixture [kg/m/s].
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
HashTable< autoPtr< multiphaseInter::phaseModel > > generatePhaseModels(const wordList &names) const
Generate the phases.
phasePairTable phasePairs_
Phase pairs.
surfaceScalarField phi_
Mixture total volumetric flux.
HashTable< autoPtr< multiphaseInter::phaseModel > > phaseModelTable
HashTable< autoPtr< porousModel >, phasePairKey, phasePairKey::hash > interfacePorousModelTable
HashTable< volScalarField::Internal > SuSpTable
virtual tmp< volScalarField > hc() const
Chemical enthalpy of the mixture [J/kg].
tmp< volScalarField > nuEff() const
Return the effective kinematic viscosity.
virtual void solve()=0
Solve for the phase transport equations.
void calcMu()
Calculate and return the laminar viscosity.
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
virtual volScalarField & he()
Return access to the internal energy field [J/Kg].
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [J/m/s/K].
virtual tmp< scalarField > THE(const scalarField &h, const scalarField &p, const scalarField &T0, const labelList &cells) const
Temperature from enthalpy/internal energy for cell-set.
volScalarField mu_
Dynamic viscocity.
void addInterfacePorosity(fvVectorMatrix &UEqn)
Add interface porosity on phasePair.
const fvMesh & mesh() const
Return mesh.
tmp< surfaceScalarField > nHatf(const volScalarField &alpha1, const volScalarField &alpha2) const
Interface normal surface vector.
virtual tmp< volScalarField > rho() const
Return the mixture density.
tmp< volScalarField > nut() const
Return the turbulent kinematic viscosity.
surfaceScalarField rhoPhi_
Mixture total mass flux.
tmp< surfaceVectorField > nHatfv(const volScalarField &alpha1, const volScalarField &alpha2) const
Interface normal surface vector.
virtual tmp< scalarField > rhoEoS(const scalarField &p, const scalarField &T, const labelList &cells) const
Density from pressure and temperature.
tmp< volScalarField > kappaEff() const
Effective thermal turbulent diffusivity for temperature.
virtual tmp< volScalarField > Cv() const
Return Cv of the mixture.
tmp< volVectorField > nVolHatfv(const volScalarField &alpha1, const volScalarField &alpha2) const
Interface normal volField vector.
void setTurbulence(compressibleTurbulenceModel &turb)
Set turbulence model.
phaseModelTable phaseModels_
Phase models.
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
virtual bool read()
Read base phaseProperties dictionary.
virtual bool isochoric() const
Return true if the equation of state is isochoric for all phasses.
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
virtual bool incompressible() const
Return true if the equation of state is incompressible for all.
tmp< surfaceScalarField > generatePhi(const HashTable< autoPtr< multiphaseInter::phaseModel > > &phaseModels) const
Generate the mixture flux.
phasePairTable totalPhasePairs_
Total ordered phase pairs in the system.
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:68
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:56
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:57
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
fvVectorMatrix & UEqn
Definition: UEqn.H:13
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
Forward declarations of fvMatrix specializations.
const cellShapeList & cells
zeroField Su
Definition: alphaSuSp.H:1
zeroField Sp
Definition: alphaSuSp.H:2
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: List.H:66
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
volScalarField & h
scalar T0
Definition: createFields.H:22
Hashing functor for phasePairKey.
Definition: phasePairKey.H:123
Foam::surfaceFields.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73