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-2018 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
118  PhaseChangeModel(CloudType& 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,
143  CloudType& owner
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 scalar dt,
159  const label celli,
160  const scalar Re,
161  const scalar Pr,
162  const scalar d,
163  const scalar nu,
164  const scalar T,
165  const scalar Ts,
166  const scalar pc,
167  const scalar Tc,
168  const scalarField& X,
169  scalarField& dMassPC
170  ) const = 0;
171 
172  //- Return the enthalpy per unit mass
173  virtual scalar dh
174  (
175  const label idc,
176  const label idl,
177  const scalar p,
178  const scalar T
179  ) const;
180 
181  //- Return vapourisation temperature
182  virtual scalar Tvap(const scalarField& X) const;
183 
184  //- Return maximum/limiting temperature
185  virtual scalar TMax(const scalar p, const scalarField& X) const;
186 
187  //- Add to phase change mass
188  void addToPhaseChangeMass(const scalar dMass);
189 
190 
191  // I-O
192 
193  //- Write injection info to stream
194  virtual void info(Ostream& os);
195 };
196 
197 
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 
200 } // End namespace Foam
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 #define makePhaseChangeModel(CloudType) \
205  \
206  typedef Foam::CloudType::reactingCloudType reactingCloudType; \
207  defineNamedTemplateTypeNameAndDebug \
208  ( \
209  Foam::PhaseChangeModel<reactingCloudType>, \
210  0 \
211  ); \
212  namespace Foam \
213  { \
214  defineTemplateRunTimeSelectionTable \
215  ( \
216  PhaseChangeModel<reactingCloudType>, \
217  dictionary \
218  ); \
219  }
220 
222 #define makePhaseChangeModelType(SS, CloudType) \
223  \
224  typedef Foam::CloudType::reactingCloudType reactingCloudType; \
225  defineNamedTemplateTypeNameAndDebug(Foam::SS<reactingCloudType>, 0); \
226  \
227  Foam::PhaseChangeModel<reactingCloudType>:: \
228  adddictionaryConstructorToTable<Foam::SS<reactingCloudType>> \
229  add##SS##CloudType##reactingCloudType##ConstructorToTable_;
230 
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 #ifdef NoRepository
235  #include "PhaseChangeModel.C"
236 #endif
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 #endif
241 
242 // ************************************************************************* //
virtual scalar TMax(const scalar p, const scalarField &X) const
Return maximum/limiting temperature.
dimensionedScalar Pr("Pr", dimless, laminarTransport)
static const Foam::wordList enthalpyTransferTypeNames
Name representations of enthalpy transfer types.
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
Templated phase change model class.
Definition: ReactingCloud.H:55
static autoPtr< PhaseChangeModel< CloudType > > New(const dictionary &dict, CloudType &owner)
Selector.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
void addToPhaseChangeMass(const scalar dMass)
Add to phase change mass.
const enthalpyTransferType & enthalpyTransfer() const
Return the enthalpy transfer type enumeration.
virtual ~PhaseChangeModel()
Destructor.
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:110
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.
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 string.
Definition: word.H:59
virtual void info(Ostream &os)
Write injection info to stream.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
enthalpyTransferType enthalpyTransfer_
Enthalpy transfer type enumeration.
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.
virtual void calculate(const scalar dt, const label celli, 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.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
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:69
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
volScalarField & nu
Namespace for OpenFOAM.
scalar Sh() const
Sherwood number.
enthalpyTransferType wordToEnthalpyTransfer(const word &etName) const
Convert word to enthalpy transfer type.