Loading...
Searching...
No Matches
EddyDiffusivity.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) 2015-2017 OpenFOAM Foundation
9 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "EddyDiffusivity.H"
30
31// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
32
33template<class BasicTurbulenceModel>
35{
36 // Read Prt if provided
37 Prt_ = dimensionedScalar("Prt", dimless, 1.0, this->coeffDict());
38 alphat_ = this->rho_*this->nut()/Prt_;
39 alphat_.correctBoundaryConditions();
40}
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
45template<class BasicTurbulenceModel>
47(
48 const word& type,
49 const alphaField& alpha,
50 const volScalarField& rho,
51 const volVectorField& U,
52 const surfaceScalarField& alphaRhoPhi,
54 const transportModel& transport,
55 const word& propertiesName
56)
57:
58 BasicTurbulenceModel
59 (
60 type,
61 alpha,
62 rho,
63 U,
64 alphaRhoPhi,
65 phi,
66 transport,
68 ),
69
70 // Cannot read Prt yet
71 Prt_("Prt", dimless, 1.0),
72
74 (
76 (
77 IOobject::groupName("alphat", alphaRhoPhi.group()),
78 this->runTime_.timeName(),
79 this->mesh_,
82 ),
83 this->mesh_
84 )
85{}
86
87
88// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89
90template<class BasicTurbulenceModel>
92{
93 if (BasicTurbulenceModel::read())
94 {
95 Prt_.readIfPresent(this->coeffDict());
96
97 return true;
98 }
99
100 return false;
101}
102
103
104template<class BasicTurbulenceModel>
108}
109
110
111// ************************************************************************* //
BasicTurbulenceModel::alphaField alphaField
virtual void correctEnergyTransport()
Correct the turbulence thermal diffusivity for energy transport.
dimensionedScalar Prt_
virtual void correctNut()
BasicTurbulenceModel::transportModel transportModel
EddyDiffusivity(const word &type, const alphaField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
Construct.
virtual bool read()
Re-read model coefficients if they have changed.
@ MUST_READ
Reading required.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word group(const word &name)
Return group (extension part of name).
Definition IOobject.C:300
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Base-class for all transport models used by the incompressible turbulence models.
static const word propertiesName
Default name of the turbulence properties dictionary.
const surfaceScalarField & alphaRhoPhi() const
Access function to phase flux field.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
scalar nut
word timeName
Definition getTimeIndex.H:3
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & alpha