buoyancyTurbSourceTemplates.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 OpenCFD Ltd
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "buoyancyTurbSource.H"
29
30// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
31
32template<class AlphaFieldType, class RhoFieldType>
33void Foam::fv::buoyancyTurbSource::buoyancyTurbSourceK
34(
35 const AlphaFieldType& alpha,
36 const RhoFieldType& rho,
37 fvMatrix<scalar>& eqn,
38 const label fieldi
39) const
40{
41 const volScalarField& k = eqn.psi();
42 const dimensionedScalar k0(k.dimensions(), SMALL);
43
44 const auto* turbPtr =
46 (
48 );
49 const volScalarField& nut = turbPtr->nut();
50
51 const dictionary& turbDict = turbPtr->coeffDict();
52 const scalar Prt
53 (
54 turbDict.getCheckOrDefault<scalar>("Prt", 0.85, scalarMinMax::ge(SMALL))
55 );
56
57 // (DTR:Eq. 21, buoyancy correction term)
58 const tmp<volScalarField> GbByK((nut/Prt)*(fvc::grad(rho) & g_)/max(k, k0));
59
60 eqn -= fvm::Sp(GbByK, k);
61}
62
63
64// ************************************************************************* //
label k
static MinMax< scalar > ge(const scalar &minVal)
A semi-infinite range from minVal to the type max.
Definition: MinMaxI.H:31
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
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.
scalar nut
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
volScalarField & alpha
dimensionedScalar Prt("Prt", dimless, laminarTransport)