Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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"
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 (
49 (
50 IOobject::groupName("alpha", phaseName),
51 mesh.time().timeName(),
52 mesh,
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",
82 phaseDict_
83 ),
84 U_
85 (
87 (
88 IOobject::groupName("U", phaseName),
89 mesh.time().timeName(),
90 mesh,
93 ),
94 mesh
95 ),
96 DDtU_
97 (
99 (
100 IOobject::groupName("DDtU", phaseName),
101 mesh.time().timeName(),
102 mesh
103 ),
104 mesh,
106 ),
107 alphaPhi_
108 (
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
122 (
123 IOobject::groupName("phi", name_),
124 mesh.time().timeName(),
125 mesh,
129 );
130
131 if (io.typeHeaderOk<surfaceScalarField>(true))
132 {
133 Info<< "Reading face flux field " << io.name() << endl;
134
135 phiPtr_.reset(new surfaceScalarField(io, mesh));
136 }
137 else
138 {
139 Info<< "Calculating face flux field " << io.name() << endl;
140
141 io.readOpt(IOobject::NO_READ);
142
143 wordList phiTypes
144 (
145 U_.boundaryField().size(),
147 );
148
149 forAll(U_.boundaryField(), i)
150 {
151 if
152 (
153 isA<fixedValueFvPatchVectorField>(U_.boundaryField()[i])
154 || isA<slipFvPatchVectorField>(U_.boundaryField()[i])
155 || isA<partialSlipFvPatchVectorField>(U_.boundaryField()[i])
156 )
157 {
158 phiTypes[i] = fixedValueFvPatchScalarField::typeName;
159 }
160 }
161
162 phiPtr_.reset
163 (
165 (
166 io,
167 fvc::flux(U_),
168 phiTypes
169 )
170 );
171 }
172
174 (
175 phaseDict_,
176 *this
177 );
178}
179
180
181// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
187// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
188
192 return nullptr;
193}
194
195
198{
199 //nuModel_->correct();
200}
201
202
203bool Foam::phaseModel::read(const dictionary& phaseDict)
204{
205 phaseDict_ = phaseDict;
206
207 //if (nuModel_->read(phaseDict_))
208 {
209 phaseDict_.readEntry("nu", nu_.value());
210 phaseDict_.readEntry("kappa", kappa_.value());
211 phaseDict_.readEntry("Cp", Cp_.value());
212 phaseDict_.readEntry("rho", rho_.value());
213
214 return true;
215 }
216
217 // return false;
218}
219
220
222{
223 surfaceScalarField::Boundary& alphaPhiBf = alphaPhi.boundaryFieldRef();
224 const volScalarField::Boundary& alphaBf = boundaryField();
225 const tmp<surfaceScalarField> tphi(phi());
226 const auto& phiBf = tphi().boundaryField();
227
228 forAll(alphaPhiBf, patchi)
229 {
230 fvsPatchScalarField& alphaPhip = alphaPhiBf[patchi];
231
232 if (!alphaPhip.coupled())
234 alphaPhip = phiBf[patchi]*alphaBf[patchi];
235 }
236 }
237}
238
239
240Foam::tmp<Foam::volScalarField> Foam::phaseModel::d() const
241{
242 return dPtr_().d();
243}
244
245
246// ************************************************************************* //
const Mesh & mesh() const noexcept
Return const reference to mesh.
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
const Boundary & boundaryField() const noexcept
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
virtual bool coupled() const
True if the patch field is coupled.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
static autoPtr< diameterModel > New(const dictionary &dict, const phaseModel &phase)
void correct()
Correct the phase properties.
Definition phaseModel.C:190
virtual ~phaseModel()
Destructor.
Definition phaseModel.C:176
phaseModel(const word &phaseName, const dictionary &phaseDict, const fvMesh &mesh)
Definition phaseModel.C:33
tmp< volScalarField > d() const
Definition phaseModel.C:233
const surfaceScalarField & alphaPhi() const
Definition phaseModel.H:228
autoPtr< phaseModel > clone() const
Return clone.
Definition phaseModel.C:182
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Definition phaseModel.C:214
virtual bool read()
Read phase properties dictionary.
Definition phaseModel.C:193
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
const auto & io
Calculate the face-flux of the given field.
word timeName
Definition getTimeIndex.H:3
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition fvcFlux.C:27
const dimensionSet dimViscosity
List< word > wordList
List of word.
Definition fileName.H:60
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const dimensionSet dimSpecificHeatCapacity(dimGasConstant)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
fvsPatchField< scalar > fvsPatchScalarField
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))