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 | Copyright (C) 2012-2016 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  fixedValueFvPatchVectorField(p, iF),
41  GammaMean_(),
42  a_(),
43  omega_()
44 {}
45 
46 
49 (
51  const fvPatch& p,
53  const fvPatchFieldMapper& mapper
54 )
55 :
56  fixedValueFvPatchVectorField(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  fixedValueFvPatchVectorField(p, iF),
72  GammaMean_(Function1<scalar>::New("GammaMean", dict)),
73  a_(Function1<scalar>::New("a", dict)),
74  omega_(Function1<scalar>::New("omega", dict))
75 {
76  fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
77 }
78 
79 
82 (
84 )
85 :
86  fixedValueFvPatchVectorField(fmfrpvf),
87  GammaMean_(fmfrpvf.GammaMean_().clone().ptr()),
88  a_(fmfrpvf.a_().clone().ptr()),
89  omega_(fmfrpvf.omega_().clone().ptr())
90 {}
91 
92 
95 (
98 )
99 :
100  fixedValueFvPatchVectorField(fmfrpvf, iF),
101  GammaMean_(fmfrpvf.GammaMean_().clone().ptr()),
102  a_(fmfrpvf.a_().clone().ptr()),
103  omega_(fmfrpvf.omega_().clone().ptr())
104 {}
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
110 {
111  if (updated())
112  {
113  return;
114  }
115 
116  const label patchi = patch().index();
117 
118  // retrieve the film region from the database
119 
120  const regionModels::regionModel& region =
122  (
123  "surfaceFilmProperties"
124  );
125 
127  dynamic_cast
128  <
130  >(region);
131 
132  // calculate the vector tangential to the patch
133  // note: normal pointing into the domain
134  const vectorField n(-patch().nf());
135 
136  // TODO: currently re-evaluating the entire gTan field to return this patch
137  const scalarField gTan(film.gTan()().boundaryField()[patchi] & n);
138 
139  if (patch().size() && (max(mag(gTan)) < SMALL))
140  {
142  << "is designed to operate on patches inclined with respect to "
143  << "gravity"
144  << endl;
145  }
146 
147  const volVectorField& nHat = film.nHat();
148 
149  const vectorField nHatp(nHat.boundaryField()[patchi].patchInternalField());
150 
151  vectorField nTan(nHatp ^ n);
152  nTan /= mag(nTan) + ROOTVSMALL;
153 
154  // calculate distance in patch tangential direction
155 
156  const vectorField& Cf = patch().Cf();
157  scalarField d(nTan & Cf);
158 
159  // calculate the wavy film height
160 
161  const scalar t = db().time().timeOutputValue();
162 
163  const scalar GMean = GammaMean_->value(t);
164  const scalar a = a_->value(t);
165  const scalar omega = omega_->value(t);
166 
167  const scalarField G(GMean + a*sin(omega*constant::mathematical::twoPi*d));
168 
169  const volScalarField& mu = film.mu();
170  const scalarField mup(mu.boundaryField()[patchi].patchInternalField());
171 
172  const volScalarField& rho = film.rho();
173  const scalarField rhop(rho.boundaryField()[patchi].patchInternalField());
174 
175  const scalarField Re(max(G, scalar(0.0))/mup);
176 
177  operator==(n*pow(gTan*mup/(3.0*rhop), 0.333)*pow(Re, 0.666));
178 
179  fixedValueFvPatchVectorField::updateCoeffs();
180 }
181 
182 
184 (
185  Ostream& os
186 ) const
187 {
189  GammaMean_->writeData(os);
190  a_->writeData(os);
191  omega_->writeData(os);
192  writeEntry("value", os);
193 }
194 
195 
196 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197 
198 namespace Foam
199 {
201  (
204  );
205 }
206 
207 
208 // ************************************************************************* //
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:53
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
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
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 dimensionedScalar G
Newtonian constant of gravitation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
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:65
volVectorField vectorField(fieldObject, mesh)
inclinedFilmNusseltInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Macros for easy insertion into run-time selection tables.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:161
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.
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
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Film velocity boundary condition for inclined films that imposes a sinusoidal perturbation on top of ...
Namespace for OpenFOAM.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363