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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2022 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 \*---------------------------------------------------------------------------*/
28 
29 #include "phaseModel.H"
30 #include "diameterModel.H"
32 #include "slipFvPatchFields.H"
34 #include "surfaceInterpolate.H"
35 #include "fvcFlux.H"
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
40 (
41  const word& phaseName,
42  const dictionary& phaseDict,
43  const fvMesh& mesh
44 )
45 :
47  (
48  IOobject
49  (
50  IOobject::groupName("alpha", phaseName),
51  mesh.time().timeName(),
52  mesh,
53  IOobject::MUST_READ,
54  IOobject::AUTO_WRITE
55  ),
56  mesh
57  ),
58  name_(phaseName),
59  phaseDict_(phaseDict),
60  nu_
61  (
62  "nu",
64  phaseDict_
65  ),
66  kappa_
67  (
68  "kappa",
69  dimensionSet(1, 1, -3, -1, 0),
70  phaseDict_
71  ),
72  Cp_
73  (
74  "Cp",
76  phaseDict_
77  ),
78  rho_
79  (
80  "rho",
81  dimDensity,
82  phaseDict_
83  ),
84  U_
85  (
86  IOobject
87  (
88  IOobject::groupName("U", phaseName),
89  mesh.time().timeName(),
90  mesh,
91  IOobject::MUST_READ,
92  IOobject::AUTO_WRITE
93  ),
94  mesh
95  ),
96  DDtU_
97  (
98  IOobject
99  (
100  IOobject::groupName("DDtU", phaseName),
101  mesh.time().timeName(),
102  mesh
103  ),
104  mesh,
106  ),
107  alphaPhi_
108  (
109  IOobject
110  (
111  IOobject::groupName("alphaPhi", phaseName),
112  mesh.time().timeName(),
113  mesh
114  ),
115  mesh,
116  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero)
117  )
118 {
119  alphaPhi_.setOriented();
120 
121  const word phiName = IOobject::groupName("phi", name_);
122 
123  IOobject phiHeader
124  (
125  phiName,
126  mesh.time().timeName(),
127  mesh,
129  );
130 
131  if (phiHeader.typeHeaderOk<surfaceScalarField>(true))
132  {
133  Info<< "Reading face flux field " << phiName << endl;
134 
135  phiPtr_.reset
136  (
138  (
139  IOobject
140  (
141  phiName,
142  mesh.time().timeName(),
143  mesh,
146  ),
147  mesh
148  )
149  );
150  }
151  else
152  {
153  Info<< "Calculating face flux field " << phiName << endl;
154 
155  wordList phiTypes
156  (
157  U_.boundaryField().size(),
159  );
160 
161  forAll(U_.boundaryField(), i)
162  {
163  if
164  (
165  isA<fixedValueFvPatchVectorField>(U_.boundaryField()[i])
166  || isA<slipFvPatchVectorField>(U_.boundaryField()[i])
167  || isA<partialSlipFvPatchVectorField>(U_.boundaryField()[i])
168  )
169  {
170  phiTypes[i] = fixedValueFvPatchScalarField::typeName;
171  }
172  }
173 
174  phiPtr_.reset
175  (
177  (
178  IOobject
179  (
180  phiName,
181  mesh.time().timeName(),
182  mesh,
185  ),
186  fvc::flux(U_),
187  phiTypes
188  )
189  );
190  }
191 
193  (
194  phaseDict_,
195  *this
196  );
197 }
198 
199 
200 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
203 {}
204 
205 
206 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
207 
209 {
211  return nullptr;
212 }
213 
214 
217 {
218  //nuModel_->correct();
219 }
220 
221 
222 bool Foam::phaseModel::read(const dictionary& phaseDict)
223 {
224  phaseDict_ = phaseDict;
225 
226  //if (nuModel_->read(phaseDict_))
227  {
228  phaseDict_.readEntry("nu", nu_.value());
229  phaseDict_.readEntry("kappa", kappa_.value());
230  phaseDict_.readEntry("Cp", Cp_.value());
231  phaseDict_.readEntry("rho", rho_.value());
232 
233  return true;
234  }
235 
236  // return false;
237 }
238 
239 
241 {
243  const volScalarField::Boundary& alphaBf = boundaryField();
244  const surfaceScalarField::Boundary& phiBf = phi().boundaryField();
245 
246  forAll(alphaPhiBf, patchi)
247  {
248  fvsPatchScalarField& alphaPhip = alphaPhiBf[patchi];
249 
250  if (!alphaPhip.coupled())
251  {
252  alphaPhip = phiBf[patchi]*alphaBf[patchi];
253  }
254  }
255 }
256 
257 
259 {
260  return dPtr_().d();
261 }
262 
263 
264 // ************************************************************************* //
fvsPatchField< scalar > fvsPatchScalarField
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Definition: phaseModel.C:233
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
autoPtr< phaseModel > clone() const
Return clone.
Definition: phaseModel.C:201
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
const dimensionSet dimViscosity
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const dimensionSet dimSpecificHeatCapacity(dimGasConstant)
Definition: dimensionSets.H:77
static const word & calculatedType() noexcept
The type name for calculated patch fields.
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
Type of boundary fields.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
word timeName
Definition: getTimeIndex.H:3
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
Calculate the face-flux of the given field.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
void setOriented(bool on=true) noexcept
Set the oriented flag: on/off.
void correct()
Correct the phase properties.
Definition: phaseModel.C:209
const Mesh & mesh() const noexcept
Return mesh.
static autoPtr< diameterModel > New(const dictionary &dict, const phaseModel &phase)
Definition: diameterModel.C:53
const dimensionSet dimDensity
phaseModel(const word &phaseName, const dictionary &phaseDict, const fvMesh &mesh)
Definition: phaseModel.C:33
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
Automatically write from objectRegistry::writeObject()
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
virtual bool read()
Read phase properties dictionary.
Definition: phaseModel.C:193
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
tmp< volScalarField > d() const
Definition: phaseModel.C:251
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
virtual ~phaseModel()
Destructor.
Definition: phaseModel.C:195
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity