kinematicThinFilm.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 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 "kinematicThinFilm.H"
31 #include "volFields.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace regionModels
38 {
39 namespace areaSurfaceFilmModels
40 {
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const word& modelType,
52  const fvPatch& patch,
53  const dictionary& dict
54 )
55 :
56  liquidFilmModel(modelType, patch, dict)
57 {}
58 
59 
60 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
61 
63 {
67 
68  // Update mass exchange sources
70 
71  // gas pressure map from primary region
72  ppf_ = pg();
73 }
74 
75 
77 {
78  if (debug)
79  {
81  }
82 
84 
85  const areaVectorField gs(g_ - ns*(ns & g_));
86 
88 
89  for (int oCorr=1; oCorr<=nOuterCorr_; oCorr++)
90  {
92 
93  faVectorMatrix UsEqn
94  (
95  fam::ddt(h_, Uf_)
96  + fam::div(phi2s_, Uf_)
97  ==
98  gs*h_
99  + turbulence_->Su(Uf_)
100  + faOptions()(h_, Uf_, sqr(dimVelocity))
101  + forces_.correct(Uf_)
102  + USp_
103  );
104 
105  UsEqn.relax();
106 
107  faOptions().constrain(UsEqn);
108 
109  if (momentumPredictor_)
110  {
111  solve(UsEqn == - fac::grad(pf_*h_)/rho_ + pf_*fac::grad(h_)/rho_);
112  }
113 
114  for (int corr=1; corr<=nCorr_; corr++)
115  {
116  areaScalarField UsA(UsEqn.A());
117 
118  Uf_ = UsEqn.H()/UsA;
120  faOptions().correct(Uf_);
121 
122  phif_ =
124  - fac::interpolate(1.0/(rho_*UsA))
126  + fac::interpolate(pf_/(rho_*UsA))
128 
129  for (int nFilm=1; nFilm<=nFilmCorr_; nFilm++)
130  {
131  faScalarMatrix hEqn
132  (
133  fam::ddt(h_)
134  + fam::div(phif_, h_)
135  ==
137  + rhoSp_
138  );
139 
140  hEqn.relax();
141  faOptions().constrain(hEqn);
142  hEqn.solve();
143  faOptions().correct(h_);
144 
145  if (nFilm == nFilmCorr_)
146  {
147  phi2s_ = hEqn.flux();
148  }
149  }
150 
151  // Bound h_
152  h_ = max(h_, h0_);
153 
156  pf_.relax();
157 
158  Uf_ -= (1.0/(rho_*UsA))*fac::grad(pf_*h_)
159  - (pf_/(rho_*UsA))*fac::grad(h_);
161  faOptions().correct(Uf_);
162  }
163  }
164 
165  Info<< "Film h min/max = " << min(h_).value() << ", "
166  << max(h_).value() << endl;
167 
168  Info<< "Film U min/max = " << max(mag(Uf_)).value() << endl;
169 }
170 
171 
173 {
174  // Reset sources
176 
177  // Correct thermo
179 
180  // Correct turbulence
181  turbulence_->correct();
182 }
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 } // End namespace areaSurfaceFilmModels
187 } // End namespace regionModels
188 } // End namespace Foam
189 
190 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
volFields.H
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase
Definition: liquidFilmBase.H:56
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::ppf_
areaScalarField ppf_
Definition: liquidFilmBase.H:110
InfoInFunction
#define InfoInFunction
Definition: messageStream.H:387
Foam::faMesh::Le
const edgeVectorField & Le() const
Definition: faMesh.C:580
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::phif_
edgeScalarField phif_
Definition: liquidFilmBase.H:113
Foam::faMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Foam::regionModels::areaSurfaceFilmModels::liquidFilmModel::forces_
forceList forces_
Definition: liquidFilmModel.H:134
Foam::faMatrix
A special matrix type and solver, designed for finite area solutions of scalar equations....
Definition: faMatricesFwd.H:37
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::regionModels::areaSurfaceFilmModels::liquidFilmModel::preEvolveRegion
virtual void preEvolveRegion()
Definition: liquidFilmModel.C:307
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::nFilmCorr_
label nFilmCorr_
Definition: liquidFilmBase.H:76
Foam::regionModels::areaSurfaceFilmModels::liquidFilmModel::postEvolveRegion
virtual void postEvolveRegion()
Definition: liquidFilmModel.C:341
Foam::regionModels::areaSurfaceFilmModels::liquidFilmModel::pnSp_
areaScalarField pnSp_
Definition: liquidFilmModel.H:104
Foam::fac::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:49
Foam::faMatrix::flux
tmp< GeometricField< Type, faePatchField, edgeMesh > > flux() const
Definition: faMatrix.C:698
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::gn_
areaScalarField gn_
Definition: liquidFilmBase.H:119
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::regionModels::areaSurfaceFilmModels::kinematicThinFilm::evolveRegion
virtual void evolveRegion()
Definition: kinematicThinFilm.C:69
Foam::regionModels::areaSurfaceFilmModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(liquidFilmBase, kinematicThinFilm, dictionary)
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::momentumPredictor_
Switch momentumPredictor_
Definition: liquidFilmBase.H:67
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::h0_
dimensionedScalar h0_
Definition: liquidFilmBase.H:82
Foam::regionModels::areaSurfaceFilmModels::kinematicThinFilm
Definition: kinematicThinFilm.H:54
Foam::regionModels::areaSurfaceFilmModels::kinematicThinFilm::postEvolveRegion
virtual void postEvolveRegion()
Definition: kinematicThinFilm.C:165
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Definition: hashSets.C:26
Foam::fam::ddt
tmp< faMatrix< Type > > ddt(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famDdt.C:41
Foam::regionModels::areaSurfaceFilmModels::liquidFilmModel::correctThermoFields
void correctThermoFields()
Definition: liquidFilmModel.C:43
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::pf_
areaScalarField pf_
Definition: liquidFilmBase.H:107
Foam::regionModels::areaSurfaceFilmModels::forceList::correct
tmp< faVectorMatrix > correct(areaVectorField &U)
Definition: forceList.C:75
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::pg
tmp< areaScalarField > pg() const
Definition: liquidFilmBase.C:390
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::g_
uniformDimensionedVectorField g_
Definition: liquidFilmBase.H:122
Foam::regionModels::areaSurfaceFilmModels::liquidFilmModel::rhoSp_
areaScalarField rhoSp_
Definition: liquidFilmModel.H:98
Foam::Info
messageStream Info
Foam::faMesh::faceAreaNormals
const areaVectorField & faceAreaNormals() const
Definition: faMesh.C:675
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::fa::optionList::constrain
void constrain(faMatrix< Type > &eqn)
Definition: faOptionListTemplates.C:242
Foam::fac::lnGrad
tmp< GeometricField< Type, faePatchField, edgeMesh > > lnGrad(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLnGrad.C:40
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Foam::regionModels::areaSurfaceFilmModels::liquidFilmModel
Definition: liquidFilmModel.H:59
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::Uf_
areaVectorField Uf_
Definition: liquidFilmBase.H:104
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
Foam::faMesh::magLe
const edgeScalarField & magLe() const
Definition: faMesh.C:591
kinematicThinFilm
Thin film model.
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::nOuterCorr_
label nOuterCorr_
Definition: liquidFilmBase.H:70
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
Foam::regionModels::areaSurfaceFilmModels::kinematicThinFilm::preEvolveRegion
virtual void preEvolveRegion()
Definition: kinematicThinFilm.C:55
Foam::faMatrix::H
tmp< GeometricField< Type, faPatchField, areaMesh > > H() const
Definition: faMatrix.C:654
Foam::regionModels::areaSurfaceFilmModels::kinematicThinFilm::kinematicThinFilm
kinematicThinFilm(const word &modelType, const fvPatch &patch, const dictionary &dict)
Definition: kinematicThinFilm.C:43
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::regionModels::regionFaModel::regionMesh
const faMesh & regionMesh() const
Definition: regionFaModelI.H:56
Foam
Definition: atmBoundaryLayer.C:26
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Definition: GeometricField.C:933
Foam::regionModels::areaSurfaceFilmModels::liquidFilmModel::turbulence_
autoPtr< filmTurbulenceModel > turbulence_
Definition: liquidFilmModel.H:119
uniformDimensionedFields.H
Foam::GeometricField::storePrevIter
void storePrevIter() const
Definition: GeometricField.C:893
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:44
Foam::regionModels::areaSurfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicThinFilm, 0)
Foam::foamVersion::patch
const std::string patch
Foam::regionModels::areaSurfaceFilmModels::liquidFilmModel::USp_
areaVectorField USp_
Definition: liquidFilmModel.H:101
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::GeometricField::relax
void relax(const scalar alpha)
Definition: GeometricField.C:965
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::faOptions
Foam::fa::options & faOptions()
Definition: liquidFilmBase.C:488
kinematicThinFilm.H
Foam::faMatrix::A
tmp< areaScalarField > A() const
Definition: faMatrix.C:627
Foam::regionModels::areaSurfaceFilmModels::liquidFilmModel::sigma_
areaScalarField sigma_
Definition: liquidFilmModel.H:89
Foam::regionModels::areaSurfaceFilmModels::liquidFilmModel::rho_
areaScalarField rho_
Definition: liquidFilmModel.H:77
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::h_
areaScalarField h_
Definition: liquidFilmBase.H:101
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:40
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:49
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::phi2s_
edgeScalarField phi2s_
Definition: liquidFilmBase.H:116
Foam::fa::optionList::correct
void correct(GeometricField< Type, faPatchField, areaMesh > &field)
Definition: faOptionListTemplates.C:283
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Foam::faMatrix::relax
void relax(const scalar alpha)
Definition: faMatrix.C:512
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::nCorr_
label nCorr_
Definition: liquidFilmBase.H:73
Foam::fam::div
tmp< faMatrix< Type > > div(const edgeScalarField &flux, const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: famDiv.C:41
liquidFilmBase
Base class for thermal 2D shells.