Loading...
Searching...
No Matches
wallBoiling.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) 2018 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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 "wallBoiling.H"
30#include "phaseSystem.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35namespace Foam
36{
37namespace diameterModels
38{
39namespace nucleationModels
40{
43 (
47 );
49}
50}
51
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
57(
58 const populationBalanceModel& popBal,
59 const dictionary& dict
60)
61:
62 nucleationModel(popBal, dict),
63 velGroup_
64 (
65 refCast<const velocityGroup>
66 (
67 popBal.mesh().lookupObject<phaseModel>
68 (
69 IOobject::groupName
70 (
71 "alpha",
72 dict.get<word>("velocityGroup")
73 )
74 ).dPtr()()
75 )
76 ),
77 turbulence_
78 (
79 popBal_.mesh().lookupObjectRef<phaseCompressibleTurbulenceModel>
80 (
81 IOobject::groupName
82 (
83 turbulenceModel::propertiesName,
84 popBal_.continuousPhase().name()
85 )
86 )
87 )
88{}
89
90
91// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92
94{
95 const tmp<volScalarField> talphat(turbulence_.alphat());
96 const volScalarField::Boundary& alphatBf = talphat().boundaryField();
97
99 alphatWallBoilingWallFunction;
100
101 forAll(alphatBf, patchi)
102 {
103 if
104 (
105 isA<alphatWallBoilingWallFunction>(alphatBf[patchi])
106 )
107 {
108 const alphatWallBoilingWallFunction& alphatw =
110
111 const scalarField& dDep = alphatw.dDeparture();
112
113 if (min(dDep) < velGroup_.sizeGroups().first().d().value())
114 {
115 Warning
116 << "Minimum departure diameter " << min(dDep)
117 << " m outside of range ["
118 << velGroup_.sizeGroups().first().d().value() << ", "
119 << velGroup_.sizeGroups().last().d().value() << "] m"
120 << " at patch " << alphatw.patch().name()
121 << endl
122 << " The nucleation rate in populationBalance "
123 << popBal_.name() << " is set to zero." << endl
124 << " Adjust discretization over property space to"
125 << " suppress this warning."
126 << endl;
127 }
128 else if (max(dDep) > velGroup_.sizeGroups().last().d().value())
129 {
130 Warning
131 << "Maximum departure diameter " << max(dDep)
132 << " m outside of range ["
133 << velGroup_.sizeGroups().first().d().value() << ", "
134 << velGroup_.sizeGroups().last().d().value() << "] m"
135 << " at patch " << alphatw.patch().name()
136 << endl
137 << " The nucleation rate in populationBalance "
138 << popBal_.name() << " is set to zero." << endl
139 << " Adjust discretization over property space to"
140 << " suppress this warning."
141 << endl;
143 }
144 }
145}
146
147
148void
150(
151 volScalarField& nucleationRate,
152 const label i
153)
154{
155 const sizeGroup& fi = popBal_.sizeGroups()[i];
156 const phaseModel& phase = fi.phase();
157 const tmp<volScalarField> trho(phase.rho());
158 const volScalarField& rho = trho();
159 const tmp<volScalarField> talphat(turbulence_.alphat());
160 const volScalarField::Boundary& alphatBf = talphat().boundaryField();
161
162 typedef compressible::alphatWallBoilingWallFunctionFvPatchScalarField
163 alphatWallBoilingWallFunction;
164
165 forAll(alphatBf, patchi)
166 {
167 if
168 (
169 isA<alphatWallBoilingWallFunction>(alphatBf[patchi])
170 )
171 {
172 const alphatWallBoilingWallFunction& alphatw =
174
175 const scalarField& dmdt = alphatw.dmdt();
176 const scalarField& dDep = alphatw.dDeparture();
177
178 const labelUList& faceCells = alphatw.patch().faceCells();
179
180 dimensionedScalar unitLength("unitLength", dimLength, 1.0);
181
182 forAll(alphatw, facei)
183 {
184 if (dmdt[facei] > SMALL)
185 {
186 const label faceCelli = faceCells[facei];
187
188 nucleationRate[faceCelli] +=
189 popBal_.gamma
190 (
191 i,
192 velGroup_.formFactor()*pow3(dDep[facei]*unitLength)
193 ).value()
194 *dmdt[facei]/rho[faceCelli]/fi.x().value();
195 }
196 }
197 }
198 }
199}
200
201
202// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Base class for nucleation models.
const populationBalanceModel & popBal() const
Return reference to the populationBalanceModel.
nucleationModel(const populationBalanceModel &popBal, const dictionary &dict)
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Wall-boiling model which requires a velocityGroup (i.e. phase) to be specified in which the nucleatio...
Definition wallBoiling.H:59
virtual void correct()
Correct diameter independent expressions.
Definition wallBoiling.C:86
virtual void addToNucleationRate(volScalarField &nucleationRate, const label i)
Add to nucleationRate.
wallBoiling(const populationBalanceModel &popBal, const dictionary &dict)
Definition wallBoiling.C:50
Class that solves the univariate population balance equation by means of a class method (also called ...
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
Definition sizeGroup.H:95
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition sizeGroupI.H:52
const phaseModel & phase() const
Return const-reference to the phase.
Definition sizeGroupI.H:31
This diameterModel is intended for use with a populationBalanceModel in order to simulate polydispers...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const Type & value() const noexcept
Return const reference to value.
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
const dimensionedScalar & rho() const
Return const-access to phase1 density.
Definition phase.H:151
A class for managing temporary objects.
Definition tmp.H:75
Abstract base class for turbulence models (RAS, LES and laminar).
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
messageStream Warning
Warning stream (stdout output on master, null elsewhere), with additional 'FOAM Warning' header text.
dictionary dict
tmp< volScalarField > trho
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299