wallBoilingPhaseChangeRateFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 
34 const hashedWordList
35 wallBoilingPhaseChangeRateFvPatchScalarField::propertyNames_
36 ({
37  "wetFraction",
38  "dDeparture",
39  "fDeparture",
40  "nucleationSiteDensity",
41  "qQuenching",
42  "qEvaporative"
43 });
44 
45 const List<const scalarField wallBoilingPhaseChangeRateFvPatchScalarField::*>
46 wallBoilingPhaseChangeRateFvPatchScalarField::propertyPtrs_
47 ({
48  &wallBoilingPhaseChangeRateFvPatchScalarField::wetFraction_,
49  &wallBoilingPhaseChangeRateFvPatchScalarField::dDeparture_,
50  &wallBoilingPhaseChangeRateFvPatchScalarField::fDeparture_,
51  &wallBoilingPhaseChangeRateFvPatchScalarField::nucleationSiteDensity_,
52  &wallBoilingPhaseChangeRateFvPatchScalarField::qQuenching_,
53  &wallBoilingPhaseChangeRateFvPatchScalarField::qEvaporative_
54 });
55 
58 
60 wallBoilingPhaseChangeRateFvPatchScalarField::propertyDimensions_
61 ({
62  &dimless,
63  &dimLength,
64  &dimRate,
65  &dimInvArea,
66  &dimHeatFlux,
68 });
69 
70 }
71 
72 
73 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * /
74 
77 (
78  const fvPatch& p,
80 )
81 :
82  calculatedFvPatchScalarField(p, iF),
83  wetFraction_(p.size(), scalar(0)),
84  dDeparture_(p.size(), vGreat),
85  fDeparture_(p.size(), scalar(0)),
86  nucleationSiteDensity_(p.size(), scalar(0)),
87  qQuenching_(p.size(), scalar(0)),
88  qEvaporative_(p.size(), scalar(0)),
89  alphatLiquid_(p.size(), scalar(0)),
90  alphatVapour_(p.size(), scalar(0))
91 {}
92 
93 
96 (
97  const fvPatch& p,
99  const dictionary& dict
100 )
101 :
102  calculatedFvPatchScalarField(p, iF, dict),
103  wetFraction_("wetFraction", dimless, dict, p.size()),
104  dDeparture_("dDeparture", dimLength, dict, p.size()),
105  fDeparture_("fDeparture", dimRate, dict, p.size()),
106  nucleationSiteDensity_("nucleationSiteDensity", dimInvArea, dict, p.size()),
107  qQuenching_("qQuenching", dimHeatFlux, dict, p.size()),
108  qEvaporative_("qEvaporative", dimHeatFlux, dict, p.size()),
109  alphatLiquid_("alphatLiquid", dimMass/dimTime/dimLength, dict, p.size()),
110  alphatVapour_("alphatVapour", dimMass/dimTime/dimLength, dict, p.size())
111 {}
112 
113 
116 (
118  const fvPatch& p,
120  const fieldMapper& mapper
121 )
122 :
123  calculatedFvPatchScalarField(psf, p, iF, mapper),
124  wetFraction_(mapper(psf.wetFraction_)),
125  dDeparture_(mapper(psf.dDeparture_)),
126  fDeparture_(mapper(psf.fDeparture_)),
127  nucleationSiteDensity_(mapper(psf.nucleationSiteDensity_)),
128  qQuenching_(mapper(psf.qQuenching_)),
129  qEvaporative_(mapper(psf.qEvaporative_)),
130  alphatLiquid_(mapper(psf.alphatLiquid_)),
131  alphatVapour_(mapper(psf.alphatVapour_))
132 {}
133 
134 
137 (
140 )
141 :
142  calculatedFvPatchScalarField(psf, iF),
143  wetFraction_(psf.wetFraction_),
144  dDeparture_(psf.dDeparture_),
145  fDeparture_(psf.fDeparture_),
146  nucleationSiteDensity_(psf.nucleationSiteDensity_),
147  qQuenching_(psf.qQuenching_),
148  qEvaporative_(psf.qEvaporative_),
149  alphatLiquid_(psf.alphatLiquid_),
150  alphatVapour_(psf.alphatVapour_)
151 {}
152 
153 
154 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
155 
156 const Foam::scalarField&
158 (
159  const word& name
160 ) const
161 {
162  return this->*propertyPtrs_[propertyNames_[name]];
163 }
164 
165 
166 const Foam::dimensionSet&
168 (
169  const word& name
170 )
171 {
172  return *propertyDimensions_[propertyNames_[name]];
173 }
174 
175 
176 const Foam::scalarField&
178 {
179  return alphatLiquid_;
180 }
181 
182 
183 const Foam::scalarField&
185 {
186  return alphatVapour_;
187 }
188 
189 
191 (
192  const fvPatchScalarField& ptf,
193  const fieldMapper& mapper
194 )
195 {
196  calculatedFvPatchScalarField::map(ptf, mapper);
197 
199  refCast<const wallBoilingPhaseChangeRateFvPatchScalarField>(ptf);
200 
201  mapper(wetFraction_, tiptf.wetFraction_);
202  mapper(dDeparture_, tiptf.dDeparture_);
203  mapper(fDeparture_, tiptf.fDeparture_);
204  mapper(nucleationSiteDensity_, tiptf.nucleationSiteDensity_);
205  mapper(qQuenching_, tiptf.qQuenching_);
206  mapper(qEvaporative_, tiptf.qEvaporative_);
207  mapper(alphatLiquid_, tiptf.alphatLiquid_);
208  mapper(alphatVapour_, tiptf.alphatVapour_);
209 }
210 
211 
213 (
214  const fvPatchScalarField& ptf
215 )
216 {
217  calculatedFvPatchScalarField::reset(ptf);
218 
220  refCast<const wallBoilingPhaseChangeRateFvPatchScalarField>(ptf);
221 
222  wetFraction_.reset(tiptf.wetFraction_);
223  dDeparture_.reset(tiptf.dDeparture_);
224  fDeparture_.reset(tiptf.fDeparture_);
225  nucleationSiteDensity_.reset(tiptf.nucleationSiteDensity_);
226  qQuenching_.reset(tiptf.qQuenching_);
227  qEvaporative_.reset(tiptf.qEvaporative_);
228  alphatLiquid_.reset(tiptf.alphatLiquid_);
229  alphatVapour_.reset(tiptf.alphatVapour_);
230 }
231 
232 
234 {
236 }
237 
238 
240 (
241  Ostream& os
242 ) const
243 {
245 
246  writeEntry(os, "wetFraction", wetFraction_);
247  writeEntry(os, "dDeparture", dDeparture_);
248  writeEntry(os, "fDeparture", fDeparture_);
249  writeEntry(os, "nucleationSiteDensity", nucleationSiteDensity_);
250  writeEntry(os, "qQuenching", qQuenching_);
251  writeEntry(os, "qEvaporative", qEvaporative_);
252  writeEntry(os, "alphatLiquid", alphatLiquid_);
253  writeEntry(os, "alphatVapour", alphatVapour_);
254 }
255 
256 
257 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258 
259 namespace Foam
260 {
262  (
265  );
266 }
267 
268 
269 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:91
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
This boundary condition is used for the phase change rate field of the wall boiling fvModel....
const scalarField & alphatVapour() const
Access the vapour turbulent thermal diffusivity.
static const dimensionSet & propertyDimensions(const word &name)
Access one of the property fields' dimensions by name.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
const scalarField & alphatLiquid() const
Access the liquid turbulent thermal diffusivity.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
const scalarField & property(const word &name) const
Access one of the property fields by name.
wallBoilingPhaseChangeRateFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
const dimensionSet dimEnergy
const dimensionSet dimRate
const dimensionSet dimless
const dimensionSet dimLength
makeNullConstructablePatchTypeField(fvPatchVectorField, noSlipFvPatchVectorField)
const dimensionSet dimTime
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
const dimensionSet dimHeatFlux(dimEnergy *inv(dimTime *dimArea))
const dimensionSet dimArea
void inv(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
const dimensionSet dimInvArea(inv(dimArea))
dictionary dict
volScalarField & p