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) 2016-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"
29#include "phaseCompressibleTurbulenceModel.H"
31#include "fvmSup.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace diameterModels
39{
40namespace IATEsources
41{
44}
45}
46}
47
48
49// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50
52(
53 const IATE& iate,
54 const dictionary& dict
55)
56:
57 IATEsource(iate)
58{}
59
60
61// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62
65(
66 const volScalarField& alphai,
67 volScalarField& kappai
68) const
69{
71 (
73 (
74 "wallBoiling:R",
75 phase().time().timeName(),
76 phase().mesh()
77 ),
78 phase().mesh(),
80 );
81
83 (
85 (
86 "wallBoiling:Rdk",
87 phase().time().timeName(),
88 phase().mesh()
89 ),
90 phase().mesh(),
92 );
93
96 (
98 (
100 otherPhase().name()
101 )
102 );
103
104 const tmp<volScalarField> talphat(turbulence.alphat());
105 const volScalarField::Boundary& alphatBf = talphat().boundaryField();
106
107 const scalarField& rho = phase().rho();
108
110 alphatWallBoilingWallFunction;
111
112 forAll(alphatBf, patchi)
113 {
114 if
115 (
116 isA<alphatWallBoilingWallFunction>(alphatBf[patchi])
117 )
118 {
119 const alphatWallBoilingWallFunction& alphatw =
120 refCast<const alphatWallBoilingWallFunction>(alphatBf[patchi]);
121
122 const scalarField& dmdt = alphatw.dmdt();
123 const scalarField& dDep = alphatw.dDeparture();
124
125 const labelList& faceCells = alphatw.patch().faceCells();
126
127 forAll(alphatw, facei)
128 {
129 if (dmdt[facei] > SMALL)
130 {
131 const label faceCelli = faceCells[facei];
132 R[faceCelli] =
133 dmdt[facei]/(alphai[faceCelli]*rho[faceCelli]);
134 Rdk[faceCelli] = R[faceCelli]*(6.0/dDep[facei]);
135 }
136 }
137 }
138 }
139
140 return Rdk - fvm::Sp(R, kappai);
141}
142
143
144// ************************************************************************* //
#define R(A, B, C, D, E, F, K, M)
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const dimensionSet & dimensions() const
Return dimensions.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:500
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
IATE (Interfacial Area Transport Equation) bubble diameter model.
Definition: IATE.H:71
IATE (Interfacial Area Transport Equation) bubble diameter model run-time selectable sources.
Definition: IATEsource.H:57
Wall-boiling model which requires a velocityGroup (i.e. phase) to be specified in which the nucleatio...
Definition: wallBoiling.H:63
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
Type & lookupObjectRef(const word &name, const bool recursive=false) const
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:57
const dimensionedScalar & rho() const
Return const-access to phase1 density.
Definition: phase.H:140
A class for managing temporary objects.
Definition: tmp.H:65
static const word propertiesName
Default name of the turbulence properties dictionary.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
compressible::turbulenceModel & turbulence
Calculate the finiteVolume matrix for implicit and explicit sources.
word timeName
Definition: getTimeIndex.H:3
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333