PhaseChangeModel.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-2022 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::PhaseChangeModel
26 
27 Description
28  Templated phase change model class
29 
30 SourceFiles
31  PhaseChangeModel.C
32  PhaseChangeModelNew.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef PhaseChangeModel_H
37 #define PhaseChangeModel_H
38 
39 #include "IOdictionary.H"
40 #include "autoPtr.H"
41 #include "runTimeSelectionTables.H"
42 #include "CloudSubModelBase.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Class PhaseChangeModel Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 template<class CloudType>
54 class PhaseChangeModel
55 :
56  public CloudSubModelBase<CloudType>
57 {
58 public:
59 
60  // Public enumerations
61 
62  //- Enthalpy transfer type
64  {
67  };
68 
69  //- Name representations of enthalpy transfer types
71 
72 
73 protected:
74 
75  // Protected data
76 
77  //- Enthalpy transfer type enumeration
79 
80 
81  // Counters
82 
83  //- Mass of lagrangian phase converted
84  scalar dMass_;
85 
86 
87  // Protected Member Functions
88 
89  //- Convert word to enthalpy transfer type
91 
92  //- Sherwood number
93  scalar Sh() const;
94 
95 
96 public:
97 
98  //- Runtime type information
99  TypeName("phaseChangeModel");
100 
101  //- Declare runtime constructor selection table
103  (
104  autoPtr,
106  dictionary,
107  (
108  const dictionary& dict,
110  ),
111  (dict, owner)
112  );
113 
114 
115  // Constructors
116 
117  //- Construct null from owner
119 
120  //- Construct from dictionary
122  (
123  const dictionary& dict,
124  CloudType& owner,
125  const word& type
126  );
127 
128  //- Construct copy
130 
131  //- Construct and return a clone
132  virtual autoPtr<PhaseChangeModel<CloudType>> clone() const = 0;
133 
134 
135  //- Destructor
136  virtual ~PhaseChangeModel();
137 
138 
139  //- Selector
141  (
142  const dictionary& dict,
144  );
145 
146 
147  // Access
148 
149  //- Return the enthalpy transfer type enumeration
150  const enthalpyTransferType& enthalpyTransfer() const;
151 
152 
153  // Member Functions
154 
155  //- Update model
156  virtual void calculate
157  (
158  const typename CloudType::parcelType& p,
159  const typename CloudType::parcelType::trackingData& td,
160  const scalar dt,
161  const scalar Re,
162  const scalar Pr,
163  const scalar d,
164  const scalar nu,
165  const scalar T,
166  const scalar Ts,
167  const scalar pc,
168  const scalar Tc,
169  const scalarField& X,
170  scalarField& dMassPC
171  ) const = 0;
172 
173  //- Return the enthalpy per unit mass
174  virtual scalar dh
175  (
176  const label idc,
177  const label idl,
178  const scalar p,
179  const scalar T
180  ) const;
181 
182  //- Return vapourisation temperature
183  virtual scalar Tvap(const scalarField& X) const;
184 
185  //- Return maximum/limiting temperature
186  virtual scalar TMax(const scalar p, const scalarField& X) const;
187 
188  //- Add to phase change mass
189  void addToPhaseChangeMass(const scalar dMass);
190 
191 
192  // I-O
193 
194  //- Write injection info to stream
195  virtual void info(Ostream& os);
196 };
197 
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 } // End namespace Foam
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 #define makePhaseChangeModel(CloudType) \
206  \
207  typedef Foam::CloudType::reactingCloudType reactingCloudType; \
208  defineNamedTemplateTypeNameAndDebug \
209  ( \
210  Foam::PhaseChangeModel<reactingCloudType>, \
211  0 \
212  ); \
213  namespace Foam \
214  { \
215  defineTemplateRunTimeSelectionTable \
216  ( \
217  PhaseChangeModel<reactingCloudType>, \
218  dictionary \
219  ); \
220  }
221 
222 
223 #define makePhaseChangeModelType(SS, CloudType) \
224  \
225  typedef Foam::CloudType::reactingCloudType reactingCloudType; \
226  defineNamedTemplateTypeNameAndDebug(Foam::SS<reactingCloudType>, 0); \
227  \
228  Foam::PhaseChangeModel<reactingCloudType>:: \
229  adddictionaryConstructorToTable<Foam::SS<reactingCloudType>> \
230  add##SS##CloudType##reactingCloudType##ConstructorToTable_;
231 
232 
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 
235 #ifdef NoRepository
236  #include "PhaseChangeModel.C"
237 #endif
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 #endif
242 
243 // ************************************************************************* //
Base class for cloud sub-models.
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:225
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Templated phase change model class.
scalar dMass_
Mass of lagrangian phase converted.
static autoPtr< PhaseChangeModel< CloudType > > New(const dictionary &dict, CloudType &owner)
Selector.
virtual ~PhaseChangeModel()
Destructor.
virtual autoPtr< PhaseChangeModel< CloudType > > clone() const =0
Construct and return a clone.
TypeName("phaseChangeModel")
Runtime type information.
virtual scalar dh(const label idc, const label idl, const scalar p, const scalar T) const
Return the enthalpy per unit mass.
const enthalpyTransferType & enthalpyTransfer() const
Return the enthalpy transfer type enumeration.
virtual void info(Ostream &os)
Write injection info to stream.
virtual void calculate(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar dt, const scalar Re, const scalar Pr, const scalar d, const scalar nu, const scalar T, const scalar Ts, const scalar pc, const scalar Tc, const scalarField &X, scalarField &dMassPC) const =0
Update model.
enthalpyTransferType enthalpyTransfer_
Enthalpy transfer type enumeration.
void addToPhaseChangeMass(const scalar dMass)
Add to phase change mass.
enthalpyTransferType
Enthalpy transfer type.
virtual scalar Tvap(const scalarField &X) const
Return vapourisation temperature.
scalar Sh() const
Sherwood number.
enthalpyTransferType wordToEnthalpyTransfer(const word &etName) const
Convert word to enthalpy transfer type.
PhaseChangeModel(CloudType &owner)
Construct null from owner.
virtual scalar TMax(const scalar p, const scalarField &X) const
Return maximum/limiting temperature.
static const Foam::wordList enthalpyTransferTypeNames
Name representations of enthalpy transfer types.
declareRunTimeSelectionTable(autoPtr, PhaseChangeModel, dictionary,(const dictionary &dict, CloudType &owner),(dict, owner))
Declare runtime constructor selection table.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:110
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
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
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Macros to ease declaration of run-time selection tables.
volScalarField & p