liquidFilmBase.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) 2020-2021 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "liquidFilmBase.H"
29 #include "faMesh.H"
30 #include "gravityMeshObject.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace regionModels
41 {
42 namespace areaSurfaceFilmModels
43 {
44 
45 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
46 
48 
50 
51 const Foam::word liquidFilmName("liquidFilm");
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
56 (
57  const word& modelType,
58  const fvPatch& p,
59  const dictionary& dict
60 )
61 :
62  regionFaModel(p, liquidFilmName, modelType, dict, true),
63 
64  momentumPredictor_
65  (
66  this->solution().subDict("PIMPLE").get<bool>("momentumPredictor")
67  ),
68  nOuterCorr_
69  (
70  this->solution().subDict("PIMPLE").get<label>("nOuterCorr")
71  ),
72  nCorr_(this->solution().subDict("PIMPLE").get<label>("nCorr")),
73  nFilmCorr_
74  (
75  this->solution().subDict("PIMPLE").get<label>("nFilmCorr")
76  ),
77 
78  h0_("h0", dimLength, 1e-7, dict),
79 
80  deltaWet_("deltaWet", dimLength, 1e-4, dict),
81 
82  UName_(dict.get<word>("U")),
83 
84  pName_(dict.lookupOrDefault<word>("p", word::null)),
85 
86  pRef_(dict.get<scalar>("pRef")),
87 
88  h_
89  (
90  IOobject
91  (
92  "hf_" + regionName_,
93  primaryMesh().time().timeName(),
94  primaryMesh(),
97  ),
98  regionMesh()
99  ),
100 
101  Uf_
102  (
103  IOobject
104  (
105  "Uf_" + regionName_,
106  primaryMesh().time().timeName(),
107  primaryMesh(),
110  ),
111  regionMesh()
112  ),
113  pf_
114  (
115  IOobject
116  (
117  "pf_" + regionName_,
118  primaryMesh().time().timeName(),
119  primaryMesh(),
122  ),
123  regionMesh(),
125  ),
126  ppf_
127  (
128  IOobject
129  (
130  "ppf_" + regionName_,
131  primaryMesh().time().timeName(),
132  primaryMesh(),
135  ),
136  regionMesh(),
138  ),
139  phif_
140  (
141  IOobject
142  (
143  "phif_" + regionName_,
144  primaryMesh().time().timeName(),
145  primaryMesh(),
148  ),
149  fac::interpolate(Uf_) & regionMesh().Le()
150  ),
151 
152  phi2s_
153  (
154  IOobject
155  (
156  "phi2s_" + regionName_,
157  primaryMesh().time().timeName(),
158  primaryMesh(),
161  ),
162  fac::interpolate(h_*Uf_) & regionMesh().Le()
163  ),
164 
165 
166  gn_
167  (
168  IOobject
169  (
170  "gn",
171  primaryMesh().time().timeName(),
172  primaryMesh(),
175  ),
176  regionMesh(),
178  ),
179 
180  g_(meshObjects::gravity::New(primaryMesh().time())),
181 
182  massSource_
183  (
184  IOobject
185  (
186  "massSource",
187  primaryMesh().time().timeName(),
188  primaryMesh()
189  ),
190  primaryMesh(),
193  ),
194 
195  momentumSource_
196  (
197  IOobject
198  (
199  "momentumSource",
200  primaryMesh().time().timeName(),
201  primaryMesh()
202  ),
203  primaryMesh(),
206  ),
207 
208  pnSource_
209  (
210  IOobject
211  (
212  "pnSource",
213  primaryMesh().time().timeName(),
214  primaryMesh()
215  ),
216  primaryMesh(),
219  ),
220 
221  energySource_
222  (
223  IOobject
224  (
225  "energySource",
226  primaryMesh().time().timeName(),
227  primaryMesh()
228  ),
229  primaryMesh(),
232  ),
233 
234  addedMassTotal_(0),
235 
236  faOptions_(Foam::fa::options::New(p))
237 {
238  const areaVectorField& ns = regionMesh().faceAreaNormals();
239 
240  gn_ = g_ & ns;
241 
242  if (!faOptions_.optionList::size())
243  {
244  Info << "No finite area options present" << endl;
245  }
246 }
247 
248 
249 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
250 
252 {}
253 
254 
255 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
256 
258 {
259  scalar CoNum = 0.0;
260  scalar velMag = 0.0;
261 
262  edgeScalarField SfUfbyDelta
263  (
265  );
266 
267  CoNum =
268  max(SfUfbyDelta/regionMesh().magLe()).value()*time().deltaT().value();
269 
270  velMag = max(mag(phif_)/regionMesh().magLe()).value();
271 
274 
275  Info<< "Max film Courant Number: " << CoNum
276  << " Film velocity magnitude: " << velMag << endl;
277 
278  return CoNum;
279 }
280 
281 
283 {
285  (
286  new areaVectorField
287  (
288  IOobject
289  (
290  "tUw",
291  primaryMesh().time().timeName(),
292  primaryMesh()
293  ),
294  regionMesh(),
296  )
297  );
298 
299  areaVectorField& Uw = tUw.ref();
300 
301  const polyPatch& pp = primaryMesh().boundaryMesh()[patch_.index()];
302 
303  if
304  (
305  primaryMesh().moving()
306  && isA<movingWallVelocityFvPatchVectorField>(pp)
307  )
308  {
310  refCast<const movingWallVelocityFvPatchVectorField>(pp);
311 
312  tmp<vectorField> tUwall = wpp.Uwall();
313 
314  // Map Up to area
315  tmp<vectorField> UsWall = vsmPtr_->mapToSurface(tUwall());
316 
317  const vectorField& nHat =
319 
320  Uw.primitiveFieldRef() = UsWall() - nHat*(UsWall() & nHat);
321  }
322 
323  return tUw;
324 }
325 
326 
328 {
330  (
331  new areaVectorField
332  (
333  IOobject
334  (
335  "tUs",
336  primaryMesh().time().timeName(),
337  primaryMesh()
338  ),
339  regionMesh(),
341  )
342  );
343 
344  // Uf quadratic profile
345  tUs.ref() = Foam::sqrt(2.0)*Uf_;
346 
347  return tUs;
348 }
349 
350 
352 {
353  const label patchi = patch_.index();
354 
355  const volVectorField& Uprimary =
357 
358  const fvPatchVectorField& Uw = Uprimary.boundaryField()[patchi];
359 
361  (
362  new areaVectorField
363  (
364  IOobject
365  (
366  "tUp",
367  primaryMesh().time().timeName(),
368  primaryMesh()
369  ),
370  regionMesh(),
372  )
373  );
374 
375  areaVectorField& Up = tUp.ref();
376 
377  scalarField hp(patch_.size(), Zero);
378 
379  // map areas h to hp on primary
380  vsmPtr_->mapToField(h_, hp);
381 
383 
384  // U tangential on primary
385  const vectorField Ust(-Uw.snGrad()*hp);
386 
387  // Map U tang to surface
388  Up.primitiveFieldRef() = vsmPtr_->mapToSurface(Ust);
389 
390  // U tangent on surface
391  Up.primitiveFieldRef() -= nHat*(Up.primitiveField() & nHat);
392 
393  return tUp;
394 }
395 
396 
398 {
400  (
401  new areaScalarField
402  (
403  IOobject
404  (
405  "tpg",
406  primaryMesh().time().timeName(),
407  primaryMesh()
408  ),
409  regionMesh(),
411  )
412  );
413 
414  areaScalarField& pfg = tpg.ref();
415 
416  if (pName_ != word::null)
417  {
418  const volScalarField& pp =
420 
421  volScalarField::Boundary& pw =
422  const_cast<volScalarField::Boundary&>(pp.boundaryField());
423 
424  //pw -= pRef_;
425 
426  pfg.primitiveFieldRef() = vsmPtr_->mapInternalToSurface<scalar>(pw)();
427  }
428  return tpg;
429 }
430 
431 
433 {
434  tmp<areaScalarField> talpha
435  (
436  new areaScalarField
437  (
438  IOobject
439  (
440  "talpha",
441  primaryMesh().time().timeName(),
442  primaryMesh()
443  ),
444  regionMesh(),
446  )
447  );
448  areaScalarField& alpha = talpha.ref();
449 
450  alpha = pos0(h_ - deltaWet_);
451 
452  return talpha;
453 }
454 
455 
457 (
458  const label patchi,
459  const label facei,
460  const scalar massSource,
461  const vector& momentumSource,
462  const scalar pressureSource,
463  const scalar energySource
464 )
465 {
466  massSource_.boundaryFieldRef()[patchi][facei] += massSource;
467  pnSource_.boundaryFieldRef()[patchi][facei] += pressureSource;
468  momentumSource_.boundaryFieldRef()[patchi][facei] += momentumSource;
469 }
470 
471 
473 {
475 }
476 
477 
479 {
480  if (debug && primaryMesh().time().writeTime())
481  {
482  massSource_.write();
483  pnSource_.write();
484  momentumSource_.write();
485  }
486 
490 
492 }
493 
494 
496 {
497  return faOptions_;
498 }
499 
500 
502 {
503  return Uf_;
504 }
505 
506 
508 {
509  return gn_;
510 }
511 
512 
514 {
515  return g_;
516 }
517 
518 
520 {
521  return h_;
522 }
523 
524 
526 {
527  return phif_;
528 }
529 
530 
532 {
533  return phi2s_;
534 }
535 
536 
538 {
539  return h0_;
540 }
541 
542 
544 {
545  return *this;
546 }
547 
548 
550 {
551  return pRef_;
552 }
553 
554 
556 {
557  return UName_;
558 }
559 
560 
561 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
562 
563 } // End namespace areaSurfaceFilmModels
564 } // End namespace regionModels
565 } // End namespace Foam
566 
567 // ************************************************************************* //
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::UName_
word UName_
Name of the velocity field.
Definition: liquidFilmBase.H:93
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< vector >
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::dimPressure
const dimensionSet dimPressure
Foam::maxOp
Definition: ops.H:223
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::phif_
edgeScalarField phif_
Film momentum flux.
Definition: liquidFilmBase.H:117
Foam::solution
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:55
Foam::fa::options::New
static options & New(const fvPatch &p)
Definition: faOptions.C:102
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::TimeState::deltaT
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:55
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::regionModels::regionFaModel::primaryMesh
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionFaModelI.H:33
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimEnergy
const dimensionSet dimEnergy
calculatedFvPatchFields.H
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::postEvolveRegion
virtual void postEvolveRegion()
Post-evolve film.
Definition: liquidFilmBase.C:478
turbulentTransportModel.H
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::alpha
tmp< areaScalarField > alpha() const
Wet indicator using h0.
Definition: liquidFilmBase.C:432
gravityMeshObject.H
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::pnSource_
volScalarField pnSource_
Normal pressure by particles.
Definition: liquidFilmBase.H:138
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::phi2s
const edgeScalarField & phi2s() const
Access continuity flux.
Definition: liquidFilmBase.C:531
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::Up
tmp< areaVectorField > Up() const
Primary region velocity at film hight. Assume the film to be.
Definition: liquidFilmBase.C:351
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::UName
word UName() const
Name of the U field.
Definition: liquidFilmBase.C:555
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::gn_
areaScalarField gn_
Normal gravity field.
Definition: liquidFilmBase.H:123
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
faMesh.H
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::Us
tmp< areaVectorField > Us() const
Film surface film velocity.
Definition: liquidFilmBase.C:327
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1258
Foam::GeometricField::internalField
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
Definition: GeometricFieldI.H:43
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve film.
Definition: liquidFilmBase.C:472
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::h0_
dimensionedScalar h0_
Smallest numerical thickness.
Definition: liquidFilmBase.H:86
movingWallVelocityFvPatchVectorField.H
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::liquidFilmBase
liquidFilmBase(const word &modelType, const fvPatch &patch, const dictionary &dict)
Construct from type name and mesh and dict.
Definition: liquidFilmBase.C:56
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::Uf
const areaVectorField & Uf() const
Access const reference Uf.
Definition: liquidFilmBase.C:501
Foam::fvPatch::size
virtual label size() const
Return size.
Definition: fvPatch.H:179
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::pg
tmp< areaScalarField > pg() const
Map primary static pressure.
Definition: liquidFilmBase.C:397
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::g
const uniformDimensionedVectorField & g() const
Gravity.
Definition: liquidFilmBase.C:513
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::g_
uniformDimensionedVectorField g_
Gravity.
Definition: liquidFilmBase.H:126
Foam::UniformDimensionedField< vector >
Foam::Field< vector >
Foam::edgeInterpolation::deltaCoeffs
const edgeScalarField & deltaCoeffs() const
Return reference to difference factors array.
Definition: edgeInterpolation.C:101
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::gn
const areaScalarField & gn() const
Access const reference gn.
Definition: liquidFilmBase.C:507
Foam::faMesh::faceAreaNormals
const areaVectorField & faceAreaNormals() const
Return face area normals.
Definition: faMesh.C:682
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::massSource_
volScalarField massSource_
Mass.
Definition: liquidFilmBase.H:132
Foam::regionModels::regionFaModel::vsmPtr_
autoPtr< volSurfaceMapping > vsmPtr_
Volume-to surface mapping.
Definition: regionFaModel.H:161
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:187
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::momentumSource_
volVectorField momentumSource_
Momentum.
Definition: liquidFilmBase.H:135
Foam::regionModels::areaSurfaceFilmModels::defineRunTimeSelectionTable
defineRunTimeSelectionTable(liquidFilmBase, dictionary)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::regionModels::regionFaModel::postEvolveRegion
virtual void postEvolveRegion()
Post-evolve region.
Definition: regionFaModel.C:195
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::Uf_
areaVectorField Uf_
Film velocity.
Definition: liquidFilmBase.H:108
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::regionModels::regionFaModel::regionMesh
const faMesh & regionMesh() const
Return the region mesh database.
Definition: regionFaModelI.H:63
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::fvPatch::index
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:197
liquidFilmBase.H
Foam::dimensioned< scalar >
Foam::calculatedFvPatchField
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
Definition: calculatedFvPatchField.H:66
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::region
const regionFaModel & region() const
Access to this region.
Definition: liquidFilmBase.C:543
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::regionModels::regionFaModel::time
const Time & time() const
Return the reference to the time database.
Definition: regionFaModelI.H:39
CoNum
scalar CoNum
Definition: compressibleCourantNo.H:31
Foam::regionModels::regionFaModel::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve region.
Definition: regionFaModel.C:187
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::pName_
word pName_
Name of the pressure field.
Definition: liquidFilmBase.H:96
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::dimAcceleration
const dimensionSet dimAcceleration
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::addSources
virtual void addSources(const label patchi, const label facei, const scalar massSource, const vector &momentumSource, const scalar pressureSource, const scalar energySource=0)
Add sources.
Definition: liquidFilmBase.C:457
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::deltaWet_
dimensionedScalar deltaWet_
Delta wet for sub-models.
Definition: liquidFilmBase.H:89
Foam::fa::options
Finite-area options.
Definition: faOptions.H:54
Foam::movingWallVelocityFvPatchVectorField::Uwall
tmp< vectorField > Uwall() const
Return wall velocity field.
Definition: movingWallVelocityFvPatchVectorField.C:97
Foam::regionModels::areaSurfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicThinFilm, 0)
Foam::movingWallVelocityFvPatchVectorField
This boundary condition provides a velocity condition for cases with moving walls.
Definition: movingWallVelocityFvPatchVectorField.H:69
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::h0
const dimensionedScalar & h0() const
Return h0.
Definition: liquidFilmBase.C:537
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::Vector< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::regionModels::regionFaModel
Base class for area region models.
Definition: regionFaModel.H:113
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::faOptions_
Foam::fa::options & faOptions_
faOptions
Definition: liquidFilmBase.H:148
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::h
const areaScalarField & h() const
Access const reference h.
Definition: liquidFilmBase.C:519
Foam::word::null
static const word null
An empty word.
Definition: word.H:80
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::pRef
scalar pRef()
Access to pRef.
Definition: liquidFilmBase.C:549
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::faOptions
Foam::fa::options & faOptions()
Return faOptions.
Definition: liquidFilmBase.C:495
Foam::regionModels::regionFaModel::patch_
const fvPatch & patch_
Reference to fvPatch.
Definition: regionFaModel.H:137
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::pRef_
scalar pRef_
Reference absolute pressure.
Definition: liquidFilmBase.H:99
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::h_
areaScalarField h_
Film hight.
Definition: liquidFilmBase.H:105
Foam::GeometricField< vector, faPatchField, areaMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::~liquidFilmBase
virtual ~liquidFilmBase()
Destructor.
Definition: liquidFilmBase.C:251
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::phi2s_
edgeScalarField phi2s_
Film height flux.
Definition: liquidFilmBase.H:120
Foam::regionModels::areaSurfaceFilmModels::liquidFilmName
const Foam::word liquidFilmName("liquidFilm")
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
turbulentFluidThermoModel.H
Foam::meshObjects::gravity::New
static const gravity & New(const Time &runTime)
Construct on Time.
Definition: gravityMeshObject.H:93
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
velMag
scalar velMag
Definition: surfaceCourantNo.H:37
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::Uw
tmp< areaVectorField > Uw() const
Wall velocity.
Definition: liquidFilmBase.C:282
liquidFilmBase
Base class for thermal 2D shells.
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::CourantNumber
virtual scalar CourantNumber() const
Courant number evaluation.
Definition: liquidFilmBase.C:257
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::phif
const edgeScalarField & phif() const
Access to momentum flux.
Definition: liquidFilmBase.C:525