atmAmbientTurbSourceTemplates.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 
29 #include "atmAmbientTurbSource.H"
30 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
33 
34 template<class AlphaFieldType, class RhoFieldType>
35 void Foam::fv::atmAmbientTurbSource::atmAmbientTurbSourceEpsilon
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  const volScalarField& epsilon = turbPtr->epsilon();
49 
50  // (Heuristically derived from RS:Eq. 4, rhs-term:5)
51  eqn +=
52  fvm::Sp(alpha()*rho()*C2_*sqr(epsilonAmb_)/(kAmb_*epsilon()), epsilon);
53 }
54 
55 
56 template<class AlphaFieldType, class RhoFieldType>
57 void Foam::fv::atmAmbientTurbSource::atmAmbientTurbSourceOmega
58 (
59  const AlphaFieldType& alpha,
60  const RhoFieldType& rho,
61  fvMatrix<scalar>& eqn,
62  const label fieldi
63 ) const
64 {
65  const auto* turbPtr =
66  mesh_.findObject<turbulenceModel>
67  (
69  );
70  const volScalarField& omega = turbPtr->omega();
72  mesh_.lookupObjectRef<volScalarField::Internal>
73  (
74  word(turbPtr->type() + ":beta")
75  );
76 
77  // (RS:Eq. 4, rhs-term:5)
78  eqn += fvm::Sp(alpha()*rho()*Cmu_*beta*sqr(omegaAmb_)/omega(), omega);
79 }
80 
81 
82 template<class AlphaFieldType, class RhoFieldType>
83 void Foam::fv::atmAmbientTurbSource::atmAmbientTurbSourceK
84 (
85  const AlphaFieldType& alpha,
86  const RhoFieldType& rho,
87  fvMatrix<scalar>& eqn,
88  const label fieldi
89 ) const
90 {
91  const auto* turbPtr =
92  mesh_.findObject<turbulenceModel>
93  (
95  );
96  const volScalarField& k = turbPtr->k();
97 
98  if (isEpsilon_)
99  {
100  // (Heuristically derived from RS:Eq. 3, rhs-term:4)
101  eqn += fvm::Sp(alpha()*rho()*epsilonAmb_/k(), k);
102  }
103  else
104  {
105  // (RS:Eq. 3, rhs-term:4)
106  eqn += fvm::Sp(alpha()*rho()*Cmu_*omegaAmb_*kAmb_/k(), k);
107  }
108 }
109 
110 
111 // ************************************************************************* //
volFields.H
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
atmAmbientTurbSource.H
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::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
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)
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
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
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
epsilon
scalar epsilon
Definition: evaluateNearWall.H:7