Loading...
Searching...
No Matches
atmPlantCanopyTurbSourceTemplates.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-------------------------------------------------------------------------------
11License
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
34template<class AlphaFieldType, class RhoFieldType>
35void Foam::fv::atmPlantCanopyTurbSource::atmPlantCanopyTurbSourceEpsilon
36(
37 const AlphaFieldType& alpha,
38 const RhoFieldType& rho,
39 fvMatrix<scalar>& eqn,
40 const label fieldi
41) const
42{
43 const auto* turbPtr =
44 mesh_.findObject<turbulenceModel>
45 (
47 );
48 const tmp<volScalarField> tepsilon(turbPtr->epsilon());
49 const volScalarField& epsilon = tepsilon();
50 const volVectorField::Internal& U = turbPtr->U()();
51
52 eqn -= fvm::Sp(alpha()*rho()*(C1_ - C2_)*calcPlantCanopyTerm(U), epsilon);
53}
54
55
56template<class AlphaFieldType, class RhoFieldType>
57void Foam::fv::atmPlantCanopyTurbSource::atmPlantCanopyTurbSourceOmega
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 tmp<volScalarField> tomega(turbPtr->omega());
71 const volScalarField& omega = tomega();
72 const volVectorField::Internal& U = turbPtr->U()();
74 mesh_.lookupObjectRef<volScalarField::Internal>
75 (
76 IOobject::scopedName(turbPtr->type(), "gamma")
77 );
79 mesh_.lookupObjectRef<volScalarField::Internal>
80 (
81 IOobject::scopedName(turbPtr->type(), "beta")
82 );
83
84 eqn -= fvm::Sp(alpha()*rho()*(gamma - beta)*calcPlantCanopyTerm(U), omega);
85}
86
87
88// ************************************************************************* //
DimensionedField< vector, volMesh > Internal
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
static const word propertiesName
Default name of the turbulence properties dictionary.
U
Definition pEqn.H:72
const scalar gamma
Definition EEqn.H:9
scalar epsilon
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
volScalarField & alpha
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)