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-2018 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 (
51  const fvPatch& p,
53  const fvPatchFieldMapper& mapper
54 )
55 :
56  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
57  GammaMean_(ptf.GammaMean_().clone().ptr()),
58  a_(ptf.a_().clone().ptr()),
59  omega_(ptf.omega_().clone().ptr())
60 {}
61 
62 
65 (
66  const fvPatch& p,
68  const dictionary& dict
69 )
70 :
71  fixedValueFvPatchScalarField(p, iF, dict),
72  GammaMean_(Function1<scalar>::New("GammaMean", dict)),
73  a_(Function1<scalar>::New("a", dict)),
74  omega_(Function1<scalar>::New("omega", dict))
75 {}
76 
77 
80 (
82 )
83 :
84  fixedValueFvPatchScalarField(wmfrhpsf),
85  GammaMean_(wmfrhpsf.GammaMean_().clone().ptr()),
86  a_(wmfrhpsf.a_().clone().ptr()),
87  omega_(wmfrhpsf.omega_().clone().ptr())
88 {}
89 
90 
93 (
96 )
97 :
98  fixedValueFvPatchScalarField(wmfrhpsf, iF),
99  GammaMean_(wmfrhpsf.GammaMean_().clone().ptr()),
100  a_(wmfrhpsf.a_().clone().ptr()),
101  omega_(wmfrhpsf.omega_().clone().ptr())
102 {}
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
108 {
109  if (updated())
110  {
111  return;
112  }
113 
114  const label patchi = patch().index();
115 
116  // retrieve the film region from the database
117 
118  const regionModels::regionModel& region =
120  (
121  "surfaceFilmProperties"
122  );
123 
125  dynamic_cast
126  <
128  >(region);
129 
130  // calculate the vector tangential to the patch
131 
132  // note: normal pointing into the domain
133  const vectorField n(-patch().nf());
134 
135  // TODO: currently re-evaluating the entire gTan field to return this patch
136  const scalarField gTan(film.gTan()().boundaryField()[patchi] & n);
137 
138  if (patch().size() && (max(mag(gTan)) < small))
139  {
141  << "is designed to operate on patches inclined with respect to "
142  << "gravity"
143  << endl;
144  }
145 
146  const volVectorField& nHat = film.nHat();
147 
148  const vectorField nHatp(nHat.boundaryField()[patchi].patchInternalField());
149 
150  vectorField nTan(nHatp ^ n);
151  nTan /= mag(nTan) + rootVSmall;
152 
153  // calculate distance in patch tangential direction
154 
155  const vectorField& Cf = patch().Cf();
156  scalarField d(nTan & Cf);
157 
158  // calculate the wavy film height
159 
160  const scalar t = db().time().timeOutputValue();
161 
162  const scalar GMean = GammaMean_->value(t);
163  const scalar a = a_->value(t);
164  const scalar omega = omega_->value(t);
165 
166  const scalarField G(GMean + a*sin(omega*constant::mathematical::twoPi*d));
167 
168  const volScalarField& mu = film.mu();
169  const scalarField mup(mu.boundaryField()[patchi].patchInternalField());
170 
171  const volScalarField& rho = film.rho();
172  const scalarField rhop(rho.boundaryField()[patchi].patchInternalField());
173 
174  const scalarField Re(max(G, scalar(0))/mup);
175 
176  operator==
177  (
178  pow(3*sqr(mup/rhop)/(gTan + rootVSmall), 1.0/3.0)*pow(Re, 1.0/3.0)
179  );
180 
181  fixedValueFvPatchScalarField::updateCoeffs();
182 }
183 
184 
186 (
187  Ostream& os
188 ) const
189 {
191  GammaMean_->writeData(os);
192  a_->writeData(os);
193  omega_->writeData(os);
194  writeEntry("value", os);
195 }
196 
197 
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 
200 namespace Foam
201 {
203  (
206  );
207 }
208 
209 
210 // ************************************************************************* //
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
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:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const dimensionedScalar G
Newtonian constant of gravitation.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
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 ...
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
Foam::fvPatchFieldMapper.
const scalar twoPi(2 *pi)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
dimensionedScalar sin(const dimensionedScalar &ds)
const dimensionedScalar mu
Atomic mass unit.
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:57
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.