atmBuoyancyTurbSourceTemplates.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 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::atmBuoyancyTurbSource::atmBuoyancyTurbSourceEpsilon
36(
37 const AlphaFieldType& alpha,
38 const RhoFieldType& rho,
39 fvMatrix<scalar>& eqn,
40 const label fieldi
41) const
42{
43 const auto* turbPtr =
45 (
47 );
48
49 // Fetch required fields from the epsilon-based model
50 const volScalarField& k = turbPtr->k();
51 const volScalarField& epsilon = turbPtr->epsilon();
52 const volScalarField::Internal& GbyNu =
54 (
55 word(turbPtr->type() + ":GbyNu")
56 );
57 const volScalarField::Internal G(GbyNu*Cmu_*sqr(k())/epsilon());
58
59 // (ARAL:Eq. 5, rhs-term:3)
60 eqn += fvm::Sp(alpha()*rho()*calcC3(k(), epsilon(), G)*B_/k(), epsilon);
61}
62
63
64template<class AlphaFieldType, class RhoFieldType>
65void Foam::fv::atmBuoyancyTurbSource::atmBuoyancyTurbSourceOmega
66(
67 const AlphaFieldType& alpha,
68 const RhoFieldType& rho,
69 fvMatrix<scalar>& eqn,
70 const label fieldi
71) const
72{
73 const auto* turbPtr =
74 mesh_.findObject<turbulenceModel>
75 (
77 );
78
79 // Fetch required fields from the omega-based model
80 const volScalarField& k = turbPtr->k();
81 const volScalarField& omega = turbPtr->omega();
82 const volScalarField::Internal& GbyNu =
83 mesh_.lookupObjectRef<volScalarField::Internal>
84 (
85 word(turbPtr->type() + ":GbyNu")
86 );
87 const volScalarField::Internal G(GbyNu*Cmu_*k()/omega());
89 mesh_.lookupObjectRef<volScalarField::Internal>
90 (
91 word(turbPtr->type() + ":gamma")
92 );
94 mesh_.lookupObjectRef<volScalarField::Internal>
95 (
96 word(turbPtr->type() + ":beta")
97 );
98
99 // (ARAL:Eq. 5, rhs-term:3)
100 eqn +=
101 fvm::Sp
102 (
103 alpha()*rho()*calcC3(k(), omega(), G, gamma, beta)*B_/k(),
104 omega
105 );
106}
107
108
109template<class AlphaFieldType, class RhoFieldType>
110void Foam::fv::atmBuoyancyTurbSource::atmBuoyancyTurbSourceK
111(
112 const AlphaFieldType& alpha,
113 const RhoFieldType& rho,
114 fvMatrix<scalar>& eqn,
115 const label fieldi
116) const
117{
118 const auto* turbPtr =
119 mesh_.findObject<turbulenceModel>
120 (
122 );
123
124 const volScalarField& k = turbPtr->k();
125
126 eqn += fvm::Sp(alpha()*rho()*B_/k(), k);
127}
128
129
130// ************************************************************************* //
label k
DimensionedField< scalar, volMesh > Internal
The internal field type from which this GeometricField is derived.
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
Type & lookupObjectRef(const word &name, const bool recursive=false) const
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
static const word propertiesName
Default name of the turbulence properties dictionary.
const scalar gamma
Definition: EEqn.H:9
scalar epsilon
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
const dimensionedScalar G
Newtonian constant of gravitation.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
volScalarField & alpha
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)