atmLengthScaleTurbSourceTemplates.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 -------------------------------------------------------------------------------
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 
30 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
33 
34 template<class AlphaFieldType, class RhoFieldType>
35 void Foam::fv::atmLengthScaleTurbSource::atmLengthScaleTurbSourceEpsilon
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::Internal& k = turbPtr->k()();
51  const volScalarField::Internal& epsilon = turbPtr->epsilon()();
52  const volScalarField::Internal& GbyNu =
54  (
55  word(turbPtr->type() + ":GbyNu")
56  );
57 
58  eqn += alpha()*rho()*calcC1Star(k, epsilon)*GbyNu*Cmu_*k;
59 }
60 
61 
62 template<class AlphaFieldType, class RhoFieldType>
63 void Foam::fv::atmLengthScaleTurbSource::atmLengthScaleTurbSourceOmega
64 (
65  const AlphaFieldType& alpha,
66  const RhoFieldType& rho,
67  fvMatrix<scalar>& eqn,
68  const label fieldi
69 ) const
70 {
71  const auto* turbPtr =
72  mesh_.findObject<turbulenceModel>
73  (
75  );
76 
77  // Fetch required fields from the omega-based model
78  const volScalarField::Internal& k = turbPtr->k()();
79  const volScalarField::Internal& omega = turbPtr->omega()();
80  const volScalarField::Internal& GbyNu =
81  mesh_.lookupObjectRef<volScalarField::Internal>
82  (
83  word(turbPtr->type() + ":GbyNu")
84  );
86  mesh_.lookupObjectRef<volScalarField::Internal>
87  (
88  word(turbPtr->type() + ":gamma")
89  );
91  mesh_.lookupObjectRef<volScalarField::Internal>
92  (
93  word(turbPtr->type() + ":beta")
94  );
95 
96  eqn += alpha()*rho()*calcGammaStar(k, omega, gamma, beta)*GbyNu;
97 }
98 
99 
100 // ************************************************************************* //
volFields.H
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
rho
rho
Definition: readInitialConditions.H:88
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
Foam::compressible::turbulenceModel
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Definition: turbulentFluidThermoModel.H:63
Foam::GeometricField< scalar, fvPatchField, volMesh >::Internal
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
Definition: GeometricField.H:107
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
atmLengthScaleTurbSource.H
Foam::objectRegistry::lookupObjectRef
Type & lookupObjectRef(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:478
Foam::objectRegistry::findObject
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
Definition: objectRegistryTemplates.C:401
gamma
const scalar gamma
Definition: EEqn.H:9
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
epsilon
scalar epsilon
Definition: evaluateNearWall.H:7