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-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "atmBuoyancyTurbSource.H"
30 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
33 
34 template<class AlphaFieldType, class RhoFieldType>
35 void 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 tmp<volScalarField> tk(turbPtr->k());
51  const volScalarField& k = tk();
52  const tmp<volScalarField> tepsilon(turbPtr->epsilon());
53  const volScalarField& epsilon = tepsilon();
54  const volScalarField::Internal& GbyNu =
56  (
57  IOobject::scopedName(turbPtr->type(), "GbyNu")
58  );
59  const volScalarField::Internal G(GbyNu*Cmu_*sqr(k())/epsilon());
60 
61  // (ARAL:Eq. 5, rhs-term:3)
62  eqn += fvm::Sp(alpha()*rho()*calcC3(k(), epsilon(), G)*B_/k(), epsilon);
63 }
64 
65 
66 template<class AlphaFieldType, class RhoFieldType>
67 void Foam::fv::atmBuoyancyTurbSource::atmBuoyancyTurbSourceOmega
68 (
69  const AlphaFieldType& alpha,
70  const RhoFieldType& rho,
71  fvMatrix<scalar>& eqn,
72  const label fieldi
73 ) const
74 {
75  const auto* turbPtr =
76  mesh_.findObject<turbulenceModel>
77  (
79  );
80 
81  // Fetch required fields from the omega-based model
82  const tmp<volScalarField> tk(turbPtr->k());
83  const volScalarField& k = tk();
84  const tmp<volScalarField> tomega(turbPtr->omega());
85  const volScalarField& omega = tomega();
86  const volScalarField::Internal& GbyNu =
87  mesh_.lookupObjectRef<volScalarField::Internal>
88  (
89  IOobject::scopedName(turbPtr->type(), "GbyNu")
90  );
91  const volScalarField::Internal G(GbyNu*Cmu_*k()/omega());
93  mesh_.lookupObjectRef<volScalarField::Internal>
94  (
95  IOobject::scopedName(turbPtr->type(), "gamma")
96  );
98  mesh_.lookupObjectRef<volScalarField::Internal>
99  (
100  IOobject::scopedName(turbPtr->type(), "beta")
101  );
102 
103  // (ARAL:Eq. 5, rhs-term:3)
104  eqn +=
105  fvm::Sp
106  (
107  alpha()*rho()*calcC3(k(), omega(), G, gamma, beta)*B_/k(),
108  omega
109  );
110 }
111 
112 
113 template<class AlphaFieldType, class RhoFieldType>
114 void Foam::fv::atmBuoyancyTurbSource::atmBuoyancyTurbSourceK
115 (
116  const AlphaFieldType& alpha,
117  const RhoFieldType& rho,
118  fvMatrix<scalar>& eqn,
119  const label fieldi
120 ) const
121 {
122  const auto* turbPtr =
123  mesh_.findObject<turbulenceModel>
124  (
126  );
127 
128  const tmp<volScalarField> tk(turbPtr->k());
129  const volScalarField& k = tk();
130 
131  eqn += fvm::Sp(alpha()*rho()*B_/k(), k);
132 }
133 
134 
135 // ************************************************************************* //
Type & lookupObjectRef(const word &name, const bool recursive=false) const
Lookup and return non-const reference to the object of the given Type. Fatal if not found or the wron...
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
const dimensionedScalar G
Newtonian constant of gravitation.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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.
DimensionedField< scalar, volMesh > Internal
The internal field type from which this GeometricField is derived.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
static const word propertiesName
Default name of the turbulence properties dictionary.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
scalar epsilon
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
const scalar gamma
Definition: EEqn.H:9
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].