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 -------------------------------------------------------------------------------
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 "wallBoiling.H"
30 #include "phaseSystem.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace diameterModels
38 {
39 namespace nucleationModels
40 {
41  defineTypeNameAndDebug(wallBoiling, 0);
43  (
44  nucleationModel,
45  wallBoiling,
46  dictionary
47  );
48 }
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 =
109  refCast<const alphatWallBoilingWallFunction>(alphatBf[patchi]);
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;
142  }
143  }
144  }
145 }
146 
147 
148 void
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 volScalarField& rho = phase.rho();
158  const tmp<volScalarField> talphat(turbulence_.alphat());
159  const volScalarField::Boundary& alphatBf = talphat().boundaryField();
160 
162  alphatWallBoilingWallFunction;
163 
164  forAll(alphatBf, patchi)
165  {
166  if
167  (
168  isA<alphatWallBoilingWallFunction>(alphatBf[patchi])
169  )
170  {
171  const alphatWallBoilingWallFunction& alphatw =
172  refCast<const alphatWallBoilingWallFunction>(alphatBf[patchi]);
173 
174  const scalarField& dmdt = alphatw.dmdt();
175  const scalarField& dDep = alphatw.dDeparture();
176 
177  const labelList& faceCells = alphatw.patch().faceCells();
178 
179  dimensionedScalar unitLength("unitLength", dimLength, 1.0);
180 
181  forAll(alphatw, facei)
182  {
183  if (dmdt[facei] > SMALL)
184  {
185  const label faceCelli = faceCells[facei];
186 
187  nucleationRate[faceCelli] +=
188  popBal_.gamma
189  (
190  i,
191  velGroup_.formFactor()*pow3(dDep[facei]*unitLength)
192  ).value()
193  *dmdt[facei]/rho[faceCelli]/fi.x().value();
194  }
195  }
196  }
197  }
198 }
199 
200 
201 // ************************************************************************* //
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
Foam::diameterModels::sizeGroup::phase
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:38
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::phase
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:54
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::diameterModels::nucleationModels::wallBoiling::correct
virtual void correct()
Correct diameter independent expressions.
Definition: wallBoiling.C:93
alphatWallBoilingWallFunctionFvPatchScalarField.H
Foam::Warning
messageStream Warning
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::diameterModels::nucleationModels::wallBoiling::addToNucleationRate
virtual void addToNucleationRate(volScalarField &nucleationRate, const label i)
Add to nucleationRate.
Definition: wallBoiling.C:150
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
wallBoiling.H
rho
rho
Definition: readInitialConditions.H:88
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField
A thermal wall function for simulation of boiling wall.
Definition: alphatWallBoilingWallFunctionFvPatchScalarField.H:438
Foam::diameterModels::nucleationModel
Base class for nucleation models.
Definition: nucleationModel.H:52
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::diameterModels::populationBalanceModel::mesh
const fvMesh & mesh() const
Return reference to the mesh.
Definition: populationBalanceModelI.H:56
Foam::diameterModels::nucleationModels::wallBoiling::wallBoiling
wallBoiling(const populationBalanceModel &popBal, const dictionary &dict)
Definition: wallBoiling.C:57
Foam::Field< scalar >
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::diameterModels::sizeGroup
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
Definition: sizeGroup.H:96
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::diameterModels::populationBalanceModel
Class that solves the univariate population balance equation by means of a class method (also called ...
Definition: populationBalanceModel.H:179
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::diameterModels::nucleationModel::popBal_
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Definition: nucleationModel.H:59
Foam::IOobject::name
const word & name() const noexcept
Return name.
Definition: IOobjectI.H:65
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: phaseCompressibleTurbulenceModelFwd.H:47
Foam::ThermalDiffusivity::alphat
virtual tmp< volScalarField > alphat() const
Return the turbulent thermal diffusivity for enthalpy [kg/m/s].
Definition: ThermalDiffusivity.C:121
Foam::List< label >
Foam::diameterModels::nucleationModels::defineTypeNameAndDebug
defineTypeNameAndDebug(constantNucleation, 0)
Foam::diameterModels::nucleationModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(nucleationModel, constantNucleation, dictionary)
Foam::diameterModels::velocityGroup::sizeGroups
const PtrList< sizeGroup > & sizeGroups() const
Return sizeGroups belonging to this velocityGroup.
Definition: velocityGroupI.H:52
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::diameterModels::sizeGroup::x
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:59
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
Foam::phase::rho
const dimensionedScalar & rho() const
Return const-access to phase1 density.
Definition: phase.H:140