thermophysicalTransportModel.H
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) 2020-2021 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 Class
25  Foam::thermophysicalTransportModel
26 
27 Description
28  Abstract base class for thermophysical transport models
29  (RAS, LES and laminar).
30 
31 SourceFiles
32  thermophysicalTransportModel.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef thermophysicalTransportModel_H
37 #define thermophysicalTransportModel_H
38 
40 #include "fluidThermo.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 /*---------------------------------------------------------------------------*\
48  Class thermophysicalTransportModel Declaration
49 \*---------------------------------------------------------------------------*/
50 
52 :
53  public IOdictionary
54 {
55 protected:
56 
57  // Protected data
58 
60 
61 
62 public:
63 
64  //- Runtime type information
65  TypeName("thermophysicalTransport");
66 
67 
68  // Constructors
69 
70  //- Construct from compressibleMomentumTransportModel
72  (
74  );
75 
76  //- Disallow default bitwise copy construction
78  (
80  ) = delete;
81 
82 
83  //- Destructor
85  {}
86 
87 
88  // Member Functions
89 
90  //- Read model coefficients if they have changed
91  virtual bool read() = 0;
92 
94  {
96  }
97 
98  //- Access function to incompressible transport model
99  virtual const fluidThermo& thermo() const = 0;
100 
101  //- Const access to the coefficients dictionary
102  virtual const dictionary& coeffDict() const = 0;
103 
104  //- Effective thermal turbulent diffusivity for temperature
105  // of mixture [W/m/K]
106  virtual tmp<volScalarField> kappaEff() const = 0;
107 
108  //- Effective thermal turbulent diffusivity for temperature
109  // of mixture for patch [W/m/K]
110  virtual tmp<scalarField> kappaEff(const label patchi) const = 0;
111 
112  //- Effective thermal turbulent diffusivity of mixture [kg/m/s]
113  virtual tmp<volScalarField> alphaEff() const = 0;
114 
115  //- Effective thermal turbulent diffusivity of mixture
116  // for patch [kg/m/s]
117  virtual tmp<scalarField> alphaEff(const label patchi) const = 0;
118 
119  //- Effective mass diffusion coefficient
120  // for a given specie mass-fraction [kg/m/s]
121  virtual tmp<volScalarField> DEff(const volScalarField& Yi) const = 0;
122 
123  //- Effective mass diffusion coefficient
124  // for a given specie mass-fraction for patch [kg/m/s]
125  virtual tmp<scalarField> DEff
126  (
127  const volScalarField& Yi,
128  const label patchi
129  ) const = 0;
130 
131  //- Return the heat flux [W/m^2]
132  virtual tmp<surfaceScalarField> q() const = 0;
133 
134  //- Return the source term for the energy equation
135  virtual tmp<fvScalarMatrix> divq(volScalarField& he) const = 0;
136 
137  //- Return the specie flux for the given specie mass-fraction [kg/m^2/s]
138  virtual tmp<surfaceScalarField> j(const volScalarField& Yi) const = 0;
139 
140  //- Return the source term for the given specie mass-fraction equation
141  virtual tmp<fvScalarMatrix> divj(volScalarField& Yi) const = 0;
142 
143  //- Solve the thermophysical transport model equations
144  // and correct the transport coefficients
145  virtual void correct() = 0;
146 
147 
148  // Member Operators
149 
150  //- Disallow default bitwise assignment
151  void operator=(const thermophysicalTransportModel&) = delete;
152 };
153 
154 
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 
157 } // End namespace Foam
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 
161 #endif
162 
163 // ************************************************************************* //
void operator=(const thermophysicalTransportModel &)=delete
Disallow default bitwise assignment.
TypeName("thermophysicalTransport")
Runtime type information.
virtual tmp< surfaceScalarField > j(const volScalarField &Yi) const =0
Return the specie flux for the given specie mass-fraction [kg/m^2/s].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual const fluidThermo & thermo() const =0
Access function to incompressible transport model.
const compressibleMomentumTransportModel & momentumTransport() const
thermophysicalTransportModel(const compressibleMomentumTransportModel &momentumTransport)
Construct from compressibleMomentumTransportModel.
virtual tmp< surfaceScalarField > q() const =0
Return the heat flux [W/m^2].
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const =0
Effective mass diffusion coefficient.
virtual bool read()=0
Read model coefficients if they have changed.
virtual void correct()=0
Solve the thermophysical transport model equations.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const =0
Return the source term for the energy equation.
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const =0
Return the source term for the given specie mass-fraction equation.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:53
thermo he()
Abstract base class for thermophysical transport models (RAS, LES and laminar).
label patchi
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent diffusivity for temperature.
virtual tmp< volScalarField > alphaEff() const =0
Effective thermal turbulent diffusivity of mixture [kg/m/s].
A class for managing temporary objects.
Definition: PtrList.H:53
virtual const dictionary & coeffDict() const =0
Const access to the coefficients dictionary.
Base class for single-phase compressible turbulence models.
Namespace for OpenFOAM.
const compressibleMomentumTransportModel & momentumTransportModel_