inclinedFilmNusseltHeightFvPatchScalarField.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"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
43  fixedValueFvPatchScalarField(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  fixedValueFvPatchScalarField(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  fixedValueFvPatchScalarField(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  fixedValueFvPatchScalarField(wmfrhpsf),
94  filmRegionName_(wmfrhpsf.filmRegionName_),
95  GammaMean_(wmfrhpsf.GammaMean_.clone()),
96  a_(wmfrhpsf.a_.clone()),
97  omega_(wmfrhpsf.omega_.clone())
98 {}
99 
100 
103 (
106 )
107 :
108  fixedValueFvPatchScalarField(wmfrhpsf, iF),
109  filmRegionName_(wmfrhpsf.filmRegionName_),
110  GammaMean_(wmfrhpsf.GammaMean_.clone()),
111  a_(wmfrhpsf.a_.clone()),
112  omega_(wmfrhpsf.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  // Calculate the vector tangential to the patch
139 
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  // Calculate distance in patch tangential direction
161 
162  const vectorField& Cf = patch().Cf();
163  scalarField d(nTan & Cf);
164 
165  // Calculate the wavy film height
166 
167  const scalar t = db().time().timeOutputValue();
168 
169  const scalar GMean = GammaMean_->value(t);
170  const scalar a = a_->value(t);
171  const scalar omega = omega_->value(t);
172 
173  const scalarField G(GMean + a*sin(omega*constant::mathematical::twoPi*d));
174 
175  const volScalarField& mu = film.mu();
176  const scalarField mup(mu.boundaryField()[patchi].patchInternalField());
177 
178  const volScalarField& rho = film.rho();
179  const scalarField rhop(rho.boundaryField()[patchi].patchInternalField());
180 
181  const scalarField Re(max(G, scalar(0))/mup);
182 
183  operator==
184  (
185  cbrt(3*sqr(mup/rhop)/(gTan + ROOTVSMALL))*cbrt(Re)
186  );
187 
188  fixedValueFvPatchScalarField::updateCoeffs();
189 }
190 
191 
193 (
194  Ostream& os
195 ) const
196 {
198 
200  (
201  "filmRegion",
202  "surfaceFilmProperties",
203  filmRegionName_
204  );
205  GammaMean_->writeData(os);
206  a_->writeData(os);
207  omega_->writeData(os);
208 
210 }
211 
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 namespace Foam
216 {
218  (
220  inclinedFilmNusseltHeightFvPatchScalarField
221  );
222 }
223 
224 
225 // ************************************************************************* //
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const dimensionedScalar G
Newtonian constant of gravitation.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
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.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Macros for easy insertion into run-time selection tables.
Film height boundary condition for inclined films that imposes a sinusoidal perturbation on top of a ...
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
fvPatchField< scalar > fvPatchScalarField
constexpr scalar twoPi(2 *M_PI)
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
dimensionedScalar cbrt(const dimensionedScalar &ds)
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:336
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
dimensionedScalar sin(const dimensionedScalar &ds)
OBJstream os(runTime.globalPath()/outputName)
const dimensionedScalar mu
Atomic mass unit.
inclinedFilmNusseltHeightFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
scalarField Re(const UList< complex > &cmplx)
Extract real component.
Definition: complexField.C:207
#define WarningInFunction
Report a warning using Foam::Warning.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.