heatTransferModel.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) 2011-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 Class
25  Foam::heatTransferModel
26 
27 Description
28  Model for heat transfer between phases
29 
30 SourceFiles
31  heatTransferModel.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef heatTransferModel_H
36 #define heatTransferModel_H
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 #include "volFields.H"
41 #include "dictionary.H"
42 #include "runTimeSelectionTables.H"
44 #include "SidedInterfacialModel.H"
45 
46 namespace Foam
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Class heatTransferModel Declaration
51 \*---------------------------------------------------------------------------*/
52 
54 :
55  public regIOobject
56 {
57 protected:
58 
59  // Protected data
60 
61  //- Residual phase fraction
63 
64 
65 public:
66 
67  //- Runtime type information
68  TypeName("heatTransferModel");
69 
70 
71  // Declare runtime construction
72 
74  (
75  autoPtr,
77  dictionary,
78  (
79  const dictionary& dict,
80  const phaseInterface& interface,
81  const bool registerObject
82  ),
83  (dict, interface, registerObject)
84  );
85 
86 
87  // Static Data Members
88 
89  //- Coefficient dimensions
90  static const dimensionSet dimK;
91 
92  //- Does this model require correcting on fixed flux boundaries?
93  static const bool correctFixedFluxBCs = false;
94 
95 
96  // Constructors
97 
98  //- Construct from a dictionary and an interface
100  (
101  const dictionary& dict,
102  const phaseInterface& interface,
103  const bool registerObject
104  );
105 
106 
107  //- Destructor
108  virtual ~heatTransferModel();
109 
110 
111  // Selectors
112 
114  (
115  const dictionary& dict,
116  const phaseInterface& interface,
117  const bool outer=true,
118  const bool registerObject=true
119  );
120 
121 
122  // Member Functions
123 
124  //- The heat transfer function K used in the enthalpy equation
125  // ddt(alpha1*rho1*ha) + ... = ... K*(Ta - Tb)
126  // ddt(alpha2*rho2*hb) + ... = ... K*(Tb - Ta)
127  tmp<volScalarField> K() const;
128 
129  //- The heat transfer function K used in the enthalpy equation
130  // ddt(alpha1*rho1*ha) + ... = ... K*(Ta - Tb)
131  // ddt(alpha2*rho2*hb) + ... = ... K*(Tb - Ta)
132  // with a specified residual volume fraction
133  virtual tmp<volScalarField> K(const scalar residualAlpha) const = 0;
134 
135  //- Dummy write for regIOobject
136  bool writeData(Ostream& os) const;
137 };
138 
139 
140 /*---------------------------------------------------------------------------*\
141  Class blendedHeatTransferModel Declaration
142 \*---------------------------------------------------------------------------*/
143 
145 :
146  public BlendedInterfacialModel<heatTransferModel>
147 {
148 public:
149 
150  // Constructors
151 
152  //- Inherit base class constructors
153  using
156 
157 
158  // Selectors
159 
161  (
162  const dictionary& dict,
164  );
165 
166 
167  // Member Functions
168 
169  //- Return the heatTransfer coefficient K
170  tmp<volScalarField> K() const;
171 
172  //- Return the heatTransfer coefficient K
173  tmp<volScalarField> K(const scalar residualAlpha) const;
174 };
175 
176 
177 /*---------------------------------------------------------------------------*\
178  Class sidedBlendedHeatTransferModel Declaration
179 \*---------------------------------------------------------------------------*/
180 
182 :
183  public SidedInterfacialModel<blendedHeatTransferModel>
184 {
185 public:
186 
187  // Constructors
188 
189  //- Inherit base class constructors
190  using
193 
194 
195  // Selectors
196 
198  (
199  const dictionary& dict,
201  );
202 };
203 
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 } // End namespace Foam
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 #endif
212 
213 // ************************************************************************* //
const phaseInterface & interface() const
Access the interface.
BlendedInterfacialModel(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
bool & registerObject()
Register object created from this IOobject with registry if true.
Definition: IOobject.H:346
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
const phaseInterface & interface() const
Access the interface.
SidedInterfacialModel(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
tmp< volScalarField > K() const
Return the heatTransfer coefficient K.
static autoPtr< blendedHeatTransferModel > New(const dictionary &dict, const phaseInterface &interface)
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
Model for heat transfer between phases.
virtual ~heatTransferModel()
Destructor.
bool writeData(Ostream &os) const
Dummy write for regIOobject.
static const bool correctFixedFluxBCs
Does this model require correcting on fixed flux boundaries?
const dimensionedScalar residualAlpha_
Residual phase fraction.
declareRunTimeSelectionTable(autoPtr, heatTransferModel, dictionary,(const dictionary &dict, const phaseInterface &interface, const bool registerObject),(dict, interface, registerObject))
heatTransferModel(const dictionary &dict, const phaseInterface &interface, const bool registerObject)
Construct from a dictionary and an interface.
tmp< volScalarField > K() const
The heat transfer function K used in the enthalpy equation.
TypeName("heatTransferModel")
Runtime type information.
static const dimensionSet dimK
Coefficient dimensions.
static autoPtr< heatTransferModel > New(const dictionary &dict, const phaseInterface &interface, const bool outer=true, const bool registerObject=true)
Class to represent an interface between phases. Derivations can further specify the configuration of ...
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
static autoPtr< sidedBlendedHeatTransferModel > New(const dictionary &dict, const phaseInterface &interface)
A class for managing temporary objects.
Definition: tmp.H:55
Namespace for OpenFOAM.
void outer(FieldField< Field1, typename outerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Macros to ease declaration of run-time selection tables.
dictionary dict