laminar.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) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "laminar.H"
31 #include "fvMesh.H"
32 #include "fvMatrices.H"
33 #include "Time.H"
34 #include "volFields.H"
35 #include "fvmSup.H"
36 #include "kinematicSingleLayer.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 namespace regionModels
44 {
45 namespace surfaceFilmModels
46 {
47 
48 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
49 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
55 laminar::laminar
56 (
58  const dictionary& dict
59 )
60 :
62  Cf_(coeffDict_.get<scalar>("Cf"))
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
69 {}
70 
71 
72 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
73 
75 {
76  auto tUs = volVectorField::New
77  (
78  IOobject::scopedName(typeName, "Us"),
83  );
84 
85  // apply quadratic profile
86  tUs.ref() = Foam::sqrt(2.0)*filmModel_.U();
87  tUs.ref().correctBoundaryConditions();
88 
89  return tUs;
90 }
91 
92 
94 {
95  return volScalarField::New
96  (
97  IOobject::scopedName(typeName, "mut"),
101  );
102 }
103 
104 
105 void laminar::correct()
106 {}
107 
108 
110 {
111  // local reference to film model
112  const kinematicSingleLayer& film =
113  static_cast<const kinematicSingleLayer&>(filmModel_);
114 
115  // local references to film fields
116  const volScalarField& mu = film.mu();
117  const volVectorField& Uw = film.Uw();
118  const volScalarField& delta = film.delta();
119  const volVectorField& Up = film.UPrimary();
120  const volScalarField& rhop = film.rhoPrimary();
121 
122  // employ simple coeff-based model
123  volScalarField Cs("Cs", Cf_*rhop*mag(Up - U));
124  volScalarField Cw("Cw", mu/((1.0/3.0)*(delta + film.deltaSmall())));
125  Cw.clamp_max(5000.0);
126 
127  return
128  (
129  - fvm::Sp(Cs, U) + Cs*Up // surface contribution
130  - fvm::Sp(Cw, U) + Cw*Uw // wall contribution
131  );
132 }
133 
134 
135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
136 
137 } // End namespace surfaceFilmModels
138 } // End namespace regionModels
139 } // End namespace Foam
140 
141 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
scalar delta
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
dimensionedScalar sqrt(const dimensionedScalar &ds)
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Macros for easy insertion into run-time selection tables.
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
virtual const volVectorField & Uw() const =0
Return the film wall velocity [m/s].
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
Film laminar turbulence model.
Definition: laminar.H:50
virtual const volScalarField & delta() const =0
Return the film thickness [m].
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:26
virtual tmp< fvVectorMatrix > Su(volVectorField &U) const
Return the source for the film momentum equation.
Definition: laminar.C:102
virtual const volVectorField & U() const =0
Return the film velocity [m/s].
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< vector >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
virtual tmp< volVectorField > Us() const
Return the film surface velocity.
Definition: laminar.C:67
void clamp_max(const Type &upper)
Impose upper (ceiling) clamp on the field values (in-place)
const dimensionedScalar mu
Atomic mass unit.
U
Definition: pEqn.H:72
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:86
virtual void correct()
Correct/update the model.
Definition: laminar.C:98
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.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Definition: fvPatchField.H:213
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Do not request registration (bool: false)
Calculate the finiteVolume matrix for implicit and explicit sources.
const scalarField & Cs
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity