phaseModel.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 | Copyright (C) 2011-2017 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 \*---------------------------------------------------------------------------*/
25 
26 #include "phaseModel.H"
27 #include "diameterModel.H"
29 #include "slipFvPatchFields.H"
31 #include "surfaceInterpolate.H"
32 #include "fvcFlux.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
37 (
38  const word& phaseName,
39  const dictionary& phaseDict,
40  const fvMesh& mesh
41 )
42 :
44  (
45  IOobject
46  (
47  IOobject::groupName("alpha", phaseName),
48  mesh.time().timeName(),
49  mesh,
50  IOobject::MUST_READ,
51  IOobject::AUTO_WRITE
52  ),
53  mesh
54  ),
55  name_(phaseName),
56  phaseDict_(phaseDict),
57  nu_
58  (
59  "nu",
60  dimensionSet(0, 2, -1, 0, 0),
61  phaseDict_
62  ),
63  kappa_
64  (
65  "kappa",
66  dimensionSet(1, 1, -3, -1, 0),
67  phaseDict_
68  ),
69  Cp_
70  (
71  "Cp",
72  dimensionSet(0, 2, -2, -1, 0),
73  phaseDict_
74  ),
75  rho_
76  (
77  "rho",
78  dimDensity,
79  phaseDict_
80  ),
81  U_
82  (
83  IOobject
84  (
85  IOobject::groupName("U", phaseName),
86  mesh.time().timeName(),
87  mesh,
88  IOobject::MUST_READ,
89  IOobject::AUTO_WRITE
90  ),
91  mesh
92  ),
93  DDtU_
94  (
95  IOobject
96  (
97  IOobject::groupName("DDtU", phaseName),
98  mesh.time().timeName(),
99  mesh
100  ),
101  mesh,
103  ),
104  alphaPhi_
105  (
106  IOobject
107  (
108  IOobject::groupName("alphaPhi", phaseName),
109  mesh.time().timeName(),
110  mesh
111  ),
112  mesh,
113  dimensionedScalar("0", dimensionSet(0, 3, -1, 0, 0), 0)
114  )
115 {
116  const word phiName = IOobject::groupName("phi", name_);
117 
118  IOobject phiHeader
119  (
120  phiName,
121  mesh.time().timeName(),
122  mesh,
124  );
125 
126  if (phiHeader.typeHeaderOk<surfaceScalarField>(true))
127  {
128  Info<< "Reading face flux field " << phiName << endl;
129 
130  phiPtr_.reset
131  (
133  (
134  IOobject
135  (
136  phiName,
137  mesh.time().timeName(),
138  mesh,
141  ),
142  mesh
143  )
144  );
145  }
146  else
147  {
148  Info<< "Calculating face flux field " << phiName << endl;
149 
150  wordList phiTypes
151  (
152  U_.boundaryField().size(),
153  calculatedFvPatchScalarField::typeName
154  );
155 
156  forAll(U_.boundaryField(), i)
157  {
158  if
159  (
160  isA<fixedValueFvPatchVectorField>(U_.boundaryField()[i])
161  || isA<slipFvPatchVectorField>(U_.boundaryField()[i])
162  || isA<partialSlipFvPatchVectorField>(U_.boundaryField()[i])
163  )
164  {
165  phiTypes[i] = fixedValueFvPatchScalarField::typeName;
166  }
167  }
168 
169  phiPtr_.reset
170  (
172  (
173  IOobject
174  (
175  phiName,
176  mesh.time().timeName(),
177  mesh,
180  ),
181  fvc::flux(U_),
182  phiTypes
183  )
184  );
185  }
186 
187  dPtr_ = diameterModel::New
188  (
189  phaseDict_,
190  *this
191  );
192 }
193 
194 
195 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
196 
198 {}
199 
200 
201 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
202 
204 {
206  return autoPtr<phaseModel>(nullptr);
207 }
208 
209 
210 
212 {
213  //nuModel_->correct();
214 }
215 
216 
217 bool Foam::phaseModel::read(const dictionary& phaseDict)
218 {
219  phaseDict_ = phaseDict;
220 
221  //if (nuModel_->read(phaseDict_))
222  {
223  phaseDict_.lookup("nu") >> nu_.value();
224  phaseDict_.lookup("kappa") >> kappa_.value();
225  phaseDict_.lookup("Cp") >> Cp_.value();
226  phaseDict_.lookup("rho") >> rho_.value();
227 
228  return true;
229  }
230  // else
231  // {
232  // return false;
233  // }
234 
235  return true;
236 }
237 
238 
240 {
241  surfaceScalarField::Boundary& alphaPhiBf = alphaPhi.boundaryFieldRef();
242  const volScalarField::Boundary& alphaBf = boundaryField();
243  const surfaceScalarField::Boundary& phiBf = phi().boundaryField();
244 
245  forAll(alphaPhiBf, patchi)
246  {
247  fvsPatchScalarField& alphaPhip = alphaPhiBf[patchi];
248 
249  if (!alphaPhip.coupled())
250  {
251  alphaPhip = phiBf[patchi]*alphaBf[patchi];
252  }
253  }
254 }
255 
256 
258 {
259  return dPtr_().d();
260 }
261 
262 
263 // ************************************************************************* //
fvsPatchField< scalar > fvsPatchScalarField
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
const surfaceScalarField & phi() const
Return the volumetric flux.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:179
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
autoPtr< phaseModel > clone() const
Return clone.
tmp< volScalarField > d() const
phaseModel(const word &phaseName, const volScalarField &p, const volScalarField &T)
Construct from components.
static autoPtr< diameterModel > New(const dictionary &dict, const phaseModel &phase)
dictionary()
Construct top-level dictionary null.
Definition: dictionary.C:238
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual ~phaseModel()
Destructor.
dynamicFvMesh & mesh
virtual bool read()
Read phase properties dictionary.
Calculate the face-flux of the given field.
static word groupName(Name name, const word &group)
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:91
void correct()
Correct the laminar viscosity.
const dimensionSet dimDensity
List< word > wordList
A List of words.
Definition: fileName.H:54
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
const dimensionSet dimVelocity