laminar.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) 2011-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 
26 #include "laminar.H"
28 #include "fvMesh.H"
29 #include "fvMatrices.H"
30 #include "Time.H"
31 #include "volFields.H"
32 #include "fvmSup.H"
33 #include "kinematicSingleLayer.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace regionModels
41 {
42 namespace surfaceFilmModels
43 {
44 
45 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
46 
47 defineTypeNameAndDebug(laminar, 0);
48 addToRunTimeSelectionTable(filmTurbulenceModel, laminar, dictionary);
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
55  const dictionary& dict
56 )
57 :
58  filmTurbulenceModel(type(), film, dict),
59  Cf_(readScalar(coeffDict_.lookup("Cf")))
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
64 
66 {}
67 
68 
69 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
70 
72 {
74  (
76  (
77  typeName + ":Us",
80  extrapolatedCalculatedFvPatchVectorField::typeName
81  )
82  );
83 
84  // apply quadratic profile
85  tUs.ref() = Foam::sqrt(2.0)*filmModel_.U();
86  tUs.ref().correctBoundaryConditions();
87 
88  return tUs;
89 }
90 
91 
93 {
94  return tmp<volScalarField>
95  (
97  (
98  typeName + ":mut",
101  )
102  );
103 }
104 
105 
107 {}
108 
109 
111 {
112  // local reference to film model
113  const kinematicSingleLayer& film =
114  static_cast<const kinematicSingleLayer&>(filmModel_);
115 
116  // local references to film fields
117  const volScalarField& mu = film.mu();
118  const volVectorField& Uw = film.Uw();
119  const volScalarField& delta = film.delta();
120  const volVectorField& Up = film.UPrimary();
121  const volScalarField& rhop = film.rhoPrimary();
122 
123  // employ simple coeff-based model
124  volScalarField Cs("Cs", Cf_*rhop*mag(Up - U));
125  volScalarField Cw("Cw", mu/((1.0/3.0)*(delta + film.deltaSmall())));
126  Cw.min(5000.0);
127 
128  return
129  (
130  - fvm::Sp(Cs, U) + Cs*Up // surface contribution
131  - fvm::Sp(Cw, U) + Cw*Uw // wall contribution
132  );
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 
138 } // End namespace surfaceFilmModels
139 } // End namespace regionModels
140 } // End namespace Foam
141 
142 // ************************************************************************* //
scalar delta
dictionary dict
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:158
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const volScalarField & rhoPrimary() const
Density [kg/m^3].
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
laminar(surfaceFilmRegionModel &film, const dictionary &dict)
Construct from surface film model.
Definition: laminar.C:53
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Macros for easy insertion into run-time selection tables.
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< vector >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
const volVectorField & UPrimary() const
Velocity [m/s].
const volScalarField & mu() const
Return const access to the dynamic viscosity [Pa.s].
const dimensionedScalar & deltaSmall() const
Return small delta.
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:61
void min(const dimensioned< Type > &)
virtual tmp< fvVectorMatrix > Su(volVectorField &U) const
Return the source for the film momentum equation.
Definition: laminar.C:110
static const zero Zero
Definition: zero.H:97
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
virtual const volVectorField & U() const =0
Return the film velocity [m/s].
virtual tmp< volVectorField > Us() const
Return the film surface velocity.
Definition: laminar.C:71
const volScalarField & delta() const
Return const access to the film thickness [m].
const dimensionedScalar mu
Atomic mass unit.
U
Definition: pEqn.H:72
virtual const volVectorField & Uw() const
Return the film wall velocity [m/s].
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volScalarField > mut() const
Return the film turbulence viscosity.
Definition: laminar.C:92
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual void correct()
Correct/update the model.
Definition: laminar.C:106
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A special matrix type and solver, designed for finite volume solutions of scalar equations.
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:53
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.
const dimensionSet dimVelocity