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,2023 OpenCFD Ltd
9 -------------------------------------------------------------------------------
10 License
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 
32 template<class AlphaFieldType, class RhoFieldType>
33 void 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 tmp<volScalarField> tnut(turbPtr->nut());
50  const volScalarField& nut = tnut();
51 
52  const dictionary& turbDict = turbPtr->coeffDict();
53  const scalar Prt
54  (
55  turbDict.getCheckOrDefault<scalar>("Prt", 0.85, scalarMinMax::ge(SMALL))
56  );
57 
58  // (DTR:Eq. 21, buoyancy correction term)
59  const tmp<volScalarField> GbByK((nut/Prt)*(fvc::grad(rho) & g_)/max(k, k0));
60 
61  eqn -= fvm::Sp(GbByK, k);
62 }
63 
64 
65 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
label k
Boltzmann constant.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
static const word propertiesName
Default name of the turbulence properties dictionary.
static MinMax< scalar > ge(const scalar &minVal)
A semi-infinite range from minVal to the type max.
Definition: MinMaxI.H:24
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
dimensionedScalar Prt("Prt", dimless, laminarTransport)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
scalar nut