PhaseChangeModel.C
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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "PhaseChangeModel.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 template<class CloudType>
34 enthalpyTransferTypeNames
35 {
36  "latentHeat", "enthalpyDifference"
37 };
38 
39 
40 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
41 
42 template<class CloudType>
45 const
46 {
47  forAll(enthalpyTransferTypeNames, i)
48  {
49  if (etName == enthalpyTransferTypeNames[i])
50  {
51  return enthalpyTransferType(i);
52  }
53  }
54 
56  << "Unknown enthalpyType " << etName << ". Valid selections are:" << nl
57  << enthalpyTransferTypeNames << exit(FatalError);
58 
59  return enthalpyTransferType(0);
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
65 template<class CloudType>
67 (
68  CloudType& owner
69 )
70 :
72  enthalpyTransfer_(etLatentHeat),
73  dMass_(0.0)
74 {}
75 
76 
77 template<class CloudType>
79 (
81 )
82 :
84  enthalpyTransfer_(pcm.enthalpyTransfer_),
85  dMass_(pcm.dMass_)
86 {}
87 
88 
89 template<class CloudType>
91 (
92  const dictionary& dict,
93  CloudType& owner,
94  const word& type
95 )
96 :
97  CloudSubModelBase<CloudType>(owner, dict, typeName, type),
98  enthalpyTransfer_
99  (
100  wordToEnthalpyTransfer(this->coeffDict().getWord("enthalpyTransfer"))
101  ),
102  dMass_(0.0)
103 {}
104 
105 
106 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 
108 template<class CloudType>
111 {
112  return enthalpyTransfer_;
113 }
114 
115 
116 template<class CloudType>
118 (
119  const label idc,
120  const label idl,
121  const scalar p,
122  const scalar T
123 ) const
124 {
125  return 0.0;
126 }
127 
128 
129 template<class CloudType>
131 (
132  const scalar p,
133  const scalarField& X
134 ) const
135 {
136  return GREAT;
137 }
138 
139 
140 template<class CloudType>
142 {
143  return -GREAT;
144 }
145 
146 
147 template<class CloudType>
149 {
150  dMass_ += dMass;
151 }
152 
153 
154 template<class CloudType>
156 {
157  const scalar mass0 = this->template getBaseProperty<scalar>("mass");
158  const scalar massTotal = mass0 + returnReduce(dMass_, sumOp<scalar>());
159 
160  Info<< " Mass transfer phase change = " << massTotal << nl;
161 
162  if (this->writeTime())
163  {
164  this->setBaseProperty("mass", massTotal);
165  dMass_ = 0.0;
166  }
167 }
168 
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 #include "PhaseChangeModelNew.C"
173 
174 // ************************************************************************* //
Foam::PhaseChangeModel
Templated phase change model class.
Definition: ReactingCloud.H:61
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::PhaseChangeModel::enthalpyTransfer_
enthalpyTransferType enthalpyTransfer_
Enthalpy transfer type enumeration.
Definition: PhaseChangeModel.H:83
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::PhaseChangeModel::dMass_
scalar dMass_
Mass of lagrangian phase converted.
Definition: PhaseChangeModel.H:89
Foam::CloudSubModelBase
Base class for cloud sub-models.
Definition: CloudSubModelBase.H:51
Foam::PhaseChangeModel::PhaseChangeModel
PhaseChangeModel(CloudType &owner)
Construct null from owner.
Definition: PhaseChangeModel.C:67
PhaseChangeModel.H
Foam::PhaseChangeModel::info
virtual void info(Ostream &os)
Write injection info to stream.
Definition: PhaseChangeModel.C:155
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::PhaseChangeModel::dh
virtual scalar dh(const label idc, const label idl, const scalar p, const scalar T) const
Return the enthalpy per unit mass.
Definition: PhaseChangeModel.C:118
Foam::PhaseChangeModel::enthalpyTransferType
enthalpyTransferType
Enthalpy transfer type.
Definition: PhaseChangeModel.H:68
PhaseChangeModelNew.C
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::PhaseChangeModel::wordToEnthalpyTransfer
enthalpyTransferType wordToEnthalpyTransfer(const word &etName) const
Convert word to enthalpy transfer type.
Definition: PhaseChangeModel.C:44
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::PhaseChangeModel::TMax
virtual scalar TMax(const scalar p, const scalarField &X) const
Return maximum/limiting temperature.
Definition: PhaseChangeModel.C:131
Foam::List< word >
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::PhaseChangeModel::Tvap
virtual scalar Tvap(const scalarField &X) const
Return vapourisation temperature.
Definition: PhaseChangeModel.C:141
Foam::PhaseChangeModel::enthalpyTransfer
const enthalpyTransferType & enthalpyTransfer() const
Return the enthalpy transfer type enumeration.
Definition: PhaseChangeModel.C:110
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::PhaseChangeModel::addToPhaseChangeMass
void addToPhaseChangeMass(const scalar dMass)
Add to phase change mass.
Definition: PhaseChangeModel.C:148