inclinedFilmNusseltInletVelocityFvPatchVectorField.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) 2012-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
30 #include "volFields.H"
31 #include "kinematicSingleLayer.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
43  fixedValueFvPatchVectorField(p, iF),
44  filmRegionName_("surfaceFilmProperties"),
45  GammaMean_(),
46  a_(),
47  omega_()
48 {}
49 
50 
53 (
55  const fvPatch& p,
57  const fvPatchFieldMapper& mapper
58 )
59 :
60  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
61  filmRegionName_(ptf.filmRegionName_),
62  GammaMean_(ptf.GammaMean_.clone()),
63  a_(ptf.a_.clone()),
64  omega_(ptf.omega_.clone())
65 {}
66 
67 
70 (
71  const fvPatch& p,
73  const dictionary& dict
74 )
75 :
76  fixedValueFvPatchVectorField(p, iF, dict),
77  filmRegionName_
78  (
79  dict.getOrDefault<word>("filmRegion", "surfaceFilmProperties")
80  ),
81  GammaMean_(Function1<scalar>::New("GammaMean", dict, &db())),
82  a_(Function1<scalar>::New("a", dict, &db())),
83  omega_(Function1<scalar>::New("omega", dict, &db()))
84 {}
85 
86 
89 (
91 )
92 :
93  fixedValueFvPatchVectorField(fmfrpvf),
94  filmRegionName_(fmfrpvf.filmRegionName_),
95  GammaMean_(fmfrpvf.GammaMean_.clone()),
96  a_(fmfrpvf.a_.clone()),
97  omega_(fmfrpvf.omega_.clone())
98 {}
99 
100 
103 (
106 )
107 :
108  fixedValueFvPatchVectorField(fmfrpvf, iF),
109  filmRegionName_(fmfrpvf.filmRegionName_),
110  GammaMean_(fmfrpvf.GammaMean_.clone()),
111  a_(fmfrpvf.a_.clone()),
112  omega_(fmfrpvf.omega_.clone())
113 {}
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 {
120  if (updated())
121  {
122  return;
123  }
124 
125  const label patchi = patch().index();
126 
127  // Retrieve the film region from the database
128 
129  const regionModels::regionModel& region =
130  db().time().lookupObject<regionModels::regionModel>(filmRegionName_);
131 
133  dynamic_cast
134  <
136  >(region);
137 
138 
139  // Calculate the vector tangential to the patch
140  // note: normal pointing into the domain
141  const vectorField n(-patch().nf());
142 
143  const scalarField gTan(film.gTan(patchi) & n);
144 
145  if (patch().size() && (max(mag(gTan)) < SMALL))
146  {
148  << "is designed to operate on patches inclined with respect to "
149  << "gravity"
150  << endl;
151  }
152 
153  const volVectorField& nHat = film.nHat();
154 
155  const vectorField nHatp(nHat.boundaryField()[patchi].patchInternalField());
156 
157  vectorField nTan(nHatp ^ n);
158  nTan /= mag(nTan) + ROOTVSMALL;
159 
160 
161  // Calculate distance in patch tangential direction
162 
163  const vectorField& Cf = patch().Cf();
164  scalarField d(nTan & Cf);
165 
166 
167  // Calculate the wavy film height
168 
169  const scalar t = db().time().timeOutputValue();
170 
171  const scalar GMean = GammaMean_->value(t);
172  const scalar a = a_->value(t);
173  const scalar omega = omega_->value(t);
174 
175  const scalarField G(GMean + a*sin(omega*constant::mathematical::twoPi*d));
176 
177  const volScalarField& mu = film.mu();
178  const scalarField mup(mu.boundaryField()[patchi].patchInternalField());
179 
180  const volScalarField& rho = film.rho();
181  const scalarField rhop(rho.boundaryField()[patchi].patchInternalField());
182 
183  const scalarField Re(max(G, scalar(0))/mup);
184 
185  operator==(n*cbrt(gTan*mup/(3*rhop))*pow(Re, 2.0/3.0));
186 
187  fixedValueFvPatchVectorField::updateCoeffs();
188 }
189 
190 
192 (
193  Ostream& os
194 ) const
195 {
198  (
199  "filmRegion",
200  "surfaceFilmProperties",
201  filmRegionName_
202  );
203  GammaMean_->writeData(os);
204  a_->writeData(os);
205  omega_->writeData(os);
206  writeEntry("value", os);
207 }
208 
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
212 namespace Foam
213 {
215  (
218  );
219 }
220 
221 
222 // ************************************************************************* //
Foam::fvPatchField< vector >
volFields.H
Foam::fvPatchField< vector >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:248
Foam::inclinedFilmNusseltInletVelocityFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: inclinedFilmNusseltInletVelocityFvPatchVectorField.C:192
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::inclinedFilmNusseltInletVelocityFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: inclinedFilmNusseltInletVelocityFvPatchVectorField.C:118
Foam::regionModels::regionModel::time
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:40
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer
Kinematic form of single-cell layer surface film model.
Definition: kinematicSingleLayer.H:67
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
rho
rho
Definition: readInitialConditions.H:88
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: propellerInfo.H:291
Foam::regionModels::regionModel
Base class for region models.
Definition: regionModel.H:60
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::Field< vector >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::inclinedFilmNusseltInletVelocityFvPatchVectorField::inclinedFilmNusseltInletVelocityFvPatchVectorField
inclinedFilmNusseltInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: inclinedFilmNusseltInletVelocityFvPatchVectorField.C:38
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
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::constant::mathematical::twoPi
constexpr scalar twoPi(2 *M_PI)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::inclinedFilmNusseltInletVelocityFvPatchVectorField
Film velocity boundary condition for inclined films that imposes a sinusoidal perturbation on top of ...
Definition: inclinedFilmNusseltInletVelocityFvPatchVectorField.H:54
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
kinematicSingleLayer.H
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::Re
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
Foam::GeometricField< vector, fvPatchField, volMesh >
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
inclinedFilmNusseltInletVelocityFvPatchVectorField.H
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54