Loading...
Searching...
No Matches
atmAmbientTurbSourceTemplates.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 ENERCON GmbH
9 Copyright (C) 2020-2023 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
30#include "volFields.H"
31
32// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
33
34template<class AlphaFieldType, class RhoFieldType>
35void Foam::fv::atmAmbientTurbSource::atmAmbientTurbSourceEpsilon
36(
37 const AlphaFieldType& alpha,
38 const RhoFieldType& rho,
39 fvMatrix<scalar>& eqn,
40 const label fieldi
41) const
42{
43 const auto* turbPtr =
44 mesh_.findObject<turbulenceModel>
45 (
47 );
48 const tmp<volScalarField> tepsilon(turbPtr->epsilon());
49 const volScalarField& epsilon = tepsilon();
50
51 // (Heuristically derived from RS:Eq. 4, rhs-term:5)
52 eqn +=
53 fvm::Sp(alpha()*rho()*C2_*sqr(epsilonAmb_)/(kAmb_*epsilon()), epsilon);
54}
55
56
57template<class AlphaFieldType, class RhoFieldType>
58void Foam::fv::atmAmbientTurbSource::atmAmbientTurbSourceOmega
59(
60 const AlphaFieldType& alpha,
61 const RhoFieldType& rho,
62 fvMatrix<scalar>& eqn,
63 const label fieldi
64) const
65{
66 const auto* turbPtr =
67 mesh_.findObject<turbulenceModel>
68 (
70 );
71 const tmp<volScalarField> tomega(turbPtr->omega());
72 const volScalarField& omega = tomega();
73
75 mesh_.lookupObjectRef<volScalarField::Internal>
76 (
77 IOobject::scopedName(turbPtr->type(), "beta")
78 );
79
80 // (RS:Eq. 4, rhs-term:5)
81 eqn += fvm::Sp(alpha()*rho()*Cmu_*beta*sqr(omegaAmb_)/omega(), omega);
82}
83
84
85template<class AlphaFieldType, class RhoFieldType>
86void Foam::fv::atmAmbientTurbSource::atmAmbientTurbSourceK
87(
88 const AlphaFieldType& alpha,
89 const RhoFieldType& rho,
90 fvMatrix<scalar>& eqn,
91 const label fieldi
92) const
93{
94 const auto* turbPtr =
95 mesh_.findObject<turbulenceModel>
96 (
98 );
99 const tmp<volScalarField> tk(turbPtr->k());
100 const volScalarField& k = tk();
101
102 if (isEpsilon_)
103 {
104 // (Heuristically derived from RS:Eq. 3, rhs-term:4)
105 eqn += fvm::Sp(alpha()*rho()*epsilonAmb_/k(), k);
106 }
107 else
108 {
109 // (RS:Eq. 3, rhs-term:4)
110 eqn += fvm::Sp(alpha()*rho()*Cmu_*omegaAmb_*kAmb_/k(), k);
111 }
112}
113
114
115// ************************************************************************* //
label k
DimensionedField< scalar, volMesh > Internal
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
static const word propertiesName
Default name of the turbulence properties dictionary.
scalar epsilon
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
volScalarField & alpha
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)