PhaseChangeModel.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::PhaseChangeModel
29 
30 Group
31  grpLagrangianIntermediatePhaseChangeSubModels
32 
33 Description
34  Templated phase change model class
35 
36 SourceFiles
37  PhaseChangeModel.C
38  PhaseChangeModelNew.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef PhaseChangeModel_H
43 #define PhaseChangeModel_H
44 
45 #include "IOdictionary.H"
46 #include "autoPtr.H"
47 #include "runTimeSelectionTables.H"
48 #include "CloudSubModelBase.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
55 /*---------------------------------------------------------------------------*\
56  Class PhaseChangeModel Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 template<class CloudType>
60 class PhaseChangeModel
61 :
62  public CloudSubModelBase<CloudType>
63 {
64 public:
65 
66  // Public Enumerations
67 
68  //- Enthalpy transfer type
70  {
73  };
74 
75  //- Name representations of enthalpy transfer types
77 
78 
79 protected:
80 
81  // Protected data
82 
83  //- Enthalpy transfer type enumeration
85 
86 
87  // Counters
88 
89  //- Mass of lagrangian phase converted
90  scalar dMass_;
91 
92 
93  // Protected Member Functions
94 
95  //- Convert word to enthalpy transfer type
97 
98  //- Sherwood number
99  scalar Sh() const;
100 
101 
102 public:
103 
104  //- Runtime type information
105  TypeName("phaseChangeModel");
106 
107  //- Declare runtime constructor selection table
109  (
110  autoPtr,
112  dictionary,
113  (
114  const dictionary& dict,
116  ),
117  (dict, owner)
118  );
119 
120 
121  // Constructors
122 
123  //- Construct null from owner
125 
126  //- Construct from dictionary
128  (
129  const dictionary& dict,
130  CloudType& owner,
131  const word& type
132  );
133 
134  //- Construct copy
136 
137  //- Construct and return a clone
138  virtual autoPtr<PhaseChangeModel<CloudType>> clone() const = 0;
139 
140 
141  //- Destructor
142  virtual ~PhaseChangeModel() = default;
143 
144 
145  //- Selector
147  (
148  const dictionary& dict,
150  );
151 
152 
153  // Access
154 
155  //- Return the enthalpy transfer type enumeration
156  const enthalpyTransferType& enthalpyTransfer() const;
157 
158 
159  // Member Functions
160 
161  //- Update model
162  virtual void calculate
163  (
164  const scalar dt,
165  const label celli,
166  const scalar Re,
167  const scalar Pr,
168  const scalar d,
169  const scalar nu,
170  const scalar rho,
171  const scalar T,
172  const scalar Ts,
173  const scalar pc,
174  const scalar Tc,
175  const scalarField& X,
176  const scalarField& solMass,
177  const scalarField& liqMass,
178  scalarField& dMassPC
179  ) const = 0;
180 
181  //- Return the enthalpy per unit mass
182  virtual scalar dh
183  (
184  const label idc,
185  const label idl,
186  const scalar p,
187  const scalar T
188  ) const;
189 
190  //- Return vapourisation temperature
191  virtual scalar Tvap(const scalarField& X) const;
192 
193  //- Return maximum/limiting temperature
194  virtual scalar TMax(const scalar p, const scalarField& X) const;
195 
196  //- Add to phase change mass
197  void addToPhaseChangeMass(const scalar dMass);
198 
199 
200  // I-O
201 
202  //- Write injection info
203  virtual void info();
204 };
205 
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 } // End namespace Foam
210 
211 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212 
213 #define makePhaseChangeModel(CloudType) \
214  \
215  typedef Foam::CloudType::reactingCloudType reactingCloudType; \
216  defineNamedTemplateTypeNameAndDebug \
217  ( \
218  Foam::PhaseChangeModel<reactingCloudType>, \
219  0 \
220  ); \
221  namespace Foam \
222  { \
223  defineTemplateRunTimeSelectionTable \
224  ( \
225  PhaseChangeModel<reactingCloudType>, \
226  dictionary \
227  ); \
228  }
229 
230 
231 #define makePhaseChangeModelType(SS, CloudType) \
232  \
233  typedef Foam::CloudType::reactingCloudType reactingCloudType; \
234  defineNamedTemplateTypeNameAndDebug(Foam::SS<reactingCloudType>, 0); \
235  \
236  Foam::PhaseChangeModel<reactingCloudType>:: \
237  adddictionaryConstructorToTable<Foam::SS<reactingCloudType>> \
238  add##SS##CloudType##reactingCloudType##ConstructorToTable_;
239 
240 
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242 
243 #ifdef NoRepository
244  #include "PhaseChangeModel.C"
245 #endif
246 
247 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248 
249 #endif
251 // ************************************************************************* //
virtual void info()
Write injection info.
virtual scalar TMax(const scalar p, const scalarField &X) const
Return maximum/limiting temperature.
dimensionedScalar Pr("Pr", dimless, laminarTransport)
static const wordList enthalpyTransferTypeNames
Name representations of enthalpy transfer types.
Templated phase change model class.
Definition: ReactingCloud.H:57
static autoPtr< PhaseChangeModel< CloudType > > New(const dictionary &dict, CloudType &owner)
Selector.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
void addToPhaseChangeMass(const scalar dMass)
Add to phase change mass.
const enthalpyTransferType & enthalpyTransfer() const
Return the enthalpy transfer type enumeration.
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:104
PhaseChangeModel(CloudType &owner)
Construct null from owner.
const CloudType & owner() const
Return const access to the owner cloud.
virtual autoPtr< PhaseChangeModel< CloudType > > clone() const =0
Construct and return a clone.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
scalar dMass_
Mass of lagrangian phase converted.
declareRunTimeSelectionTable(autoPtr, PhaseChangeModel, dictionary,(const dictionary &dict, CloudType &owner),(dict, owner))
Declare runtime constructor selection table.
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual void calculate(const scalar dt, const label celli, const scalar Re, const scalar Pr, const scalar d, const scalar nu, const scalar rho, const scalar T, const scalar Ts, const scalar pc, const scalar Tc, const scalarField &X, const scalarField &solMass, const scalarField &liqMass, scalarField &dMassPC) const =0
Update model.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
enthalpyTransferType enthalpyTransfer_
Enthalpy transfer type enumeration.
scalarField Re(const UList< complex > &cmplx)
Extract real component.
Definition: complexField.C:207
virtual scalar dh(const label idc, const label idl, const scalar p, const scalar T) const
Return the enthalpy per unit mass.
TypeName("phaseChangeModel")
Runtime type information.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Macros to ease declaration of run-time selection tables.
volScalarField & p
virtual scalar Tvap(const scalarField &X) const
Return vapourisation temperature.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
volScalarField & nu
Namespace for OpenFOAM.
scalar Sh() const
Sherwood number.
virtual ~PhaseChangeModel()=default
Destructor.
enthalpyTransferType wordToEnthalpyTransfer(const word &etName) const
Convert word to enthalpy transfer type.