Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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"
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
132 const regionModels::surfaceFilmModels::kinematicSingleLayer& film =
133 dynamic_cast
134 <
135 const regionModels::surfaceFilmModels::kinematicSingleLayer&
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);
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{
197
198 os.writeEntryIfDifferent<word>
199 (
200 "filmRegion",
201 "surfaceFilmProperties",
202 filmRegionName_
203 );
204 GammaMean_->writeData(os);
205 a_->writeData(os);
206 omega_->writeData(os);
207
210
211
212// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213
214namespace Foam
215{
217 (
219 inclinedFilmNusseltInletVelocityFvPatchVectorField
220 );
221}
222
223
224// ************************************************************************* //
label n
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...
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
Film velocity boundary condition for inclined films that imposes a sinusoidal perturbation on top of ...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
inclinedFilmNusseltInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Base class for region models.
Definition regionModel.H:59
const Time & time() const noexcept
Return the reference to the time database.
Kinematic form of single-cell layer surface film model.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
#define WarningInFunction
Report a warning using Foam::Warning.
constexpr scalar twoPi(2 *M_PI)
const dimensionedScalar G
Newtonian constant of gravitation.
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar sin(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
scalarField Re(const UList< complex > &cmplx)
Extract real component.
dimensionedScalar cbrt(const dimensionedScalar &ds)
fvPatchField< vector > fvPatchVectorField
dictionary dict