inclinedFilmNusseltHeightFvPatchScalarField.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) 2012-2021 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 
27 #include "volFields.H"
28 #include "kinematicSingleLayer.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
35 (
36  const fvPatch& p,
38 )
39 :
40  fixedValueFvPatchScalarField(p, iF),
41  GammaMean_(),
42  a_(),
43  omega_()
44 {}
45 
46 
49 (
50  const fvPatch& p,
52  const dictionary& dict
53 )
54 :
55  fixedValueFvPatchScalarField(p, iF, dict),
56  GammaMean_(Function1<scalar>::New("GammaMean", dict)),
57  a_(Function1<scalar>::New("a", dict)),
58  omega_(Function1<scalar>::New("omega", dict))
59 {}
60 
61 
64 (
66  const fvPatch& p,
68  const fvPatchFieldMapper& mapper
69 )
70 :
71  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
72  GammaMean_(ptf.GammaMean_().clone().ptr()),
73  a_(ptf.a_().clone().ptr()),
74  omega_(ptf.omega_().clone().ptr())
75 {}
76 
77 
80 (
83 )
84 :
85  fixedValueFvPatchScalarField(wmfrhpsf, iF),
86  GammaMean_(wmfrhpsf.GammaMean_().clone().ptr()),
87  a_(wmfrhpsf.a_().clone().ptr()),
88  omega_(wmfrhpsf.omega_().clone().ptr())
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
95 {
96  if (updated())
97  {
98  return;
99  }
100 
101  const label patchi = patch().index();
102 
103  // retrieve the film region from the database
104 
105  const regionModels::regionModel& region =
107  (
108  "surfaceFilmProperties"
109  );
110 
112  dynamic_cast
113  <
115  >(region);
116 
117  // calculate the vector tangential to the patch
118 
119  // note: normal pointing into the domain
120  const vectorField n(-patch().nf());
121 
122  const scalarField gTan
123  (
124  (
125  film.g().value()
126  - film.nHat().boundaryField()[patchi]
127  *(film.g().value() & film.nHat().boundaryField()[patchi])
128  ) & n
129  );
130 
131  if (patch().size() && (max(mag(gTan)) < small))
132  {
134  << "is designed to operate on patches inclined with respect to "
135  << "gravity"
136  << endl;
137  }
138 
139  const volVectorField& nHat = film.nHat();
140 
141  const vectorField nHatp(nHat.boundaryField()[patchi].patchInternalField());
142 
143  vectorField nTan(nHatp ^ n);
144  nTan /= mag(nTan) + rootVSmall;
145 
146  // calculate distance in patch tangential direction
147 
148  const vectorField& Cf = patch().Cf();
149  scalarField d(nTan & Cf);
150 
151  // calculate the wavy film height
152 
153  const scalar t = db().time().timeOutputValue();
154 
155  const scalar GMean = GammaMean_->value(t);
156  const scalar a = a_->value(t);
157  const scalar omega = omega_->value(t);
158 
159  const scalarField G(GMean + a*sin(omega*constant::mathematical::twoPi*d));
160 
161  const volScalarField& mu = film.mu();
162  const scalarField mup(mu.boundaryField()[patchi].patchInternalField());
163 
164  const volScalarField& rho = film.rho();
165  const scalarField rhop(rho.boundaryField()[patchi].patchInternalField());
166 
167  const scalarField Re(max(G, scalar(0))/mup);
168 
169  operator==
170  (
171  pow(3*sqr(mup/rhop)/(gTan + rootVSmall), 1.0/3.0)*pow(Re, 1.0/3.0)
172  );
173 
174  fixedValueFvPatchScalarField::updateCoeffs();
175 }
176 
177 
179 (
180  Ostream& os
181 ) const
182 {
184  writeEntry(os, GammaMean_());
185  writeEntry(os, a_());
186  writeEntry(os, omega_());
187  writeEntry(os, "value", *this);
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 namespace Foam
194 {
196  (
199  );
200 }
201 
202 
203 // ************************************************************************* //
Run-time selectable general function of one variable.
Definition: Function1.H:52
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Kinematic form of single-cell layer surface film model.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
T * ptr()
Return object pointer for reuse.
Definition: autoPtrI.H:90
const dimensionedScalar G
Newtonian constant of gravitation.
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
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 ...
virtual Type value(const scalar x) const =0
Return value as a function of scalar x.
Foam::fvPatchFieldMapper.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
const scalar twoPi(2 *pi)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
dimensionedScalar sin(const dimensionedScalar &ds)
const dimensionedScalar mu
Atomic mass unit.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
inclinedFilmNusseltHeightFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label patchi
#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...
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
Base class for region models.
Definition: regionModel.H:55
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Namespace for OpenFOAM.