atmPlantCanopyTurbSource.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-2021 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
31
32// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace fv
37{
40}
41}
42
43
44// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
45
47Foam::fv::atmPlantCanopyTurbSource::calcPlantCanopyTerm
48(
49 const volVectorField::Internal& U
50) const
51{
52 // (SP:Eq. 42)
53 return 12.0*Foam::sqrt(Cmu_)*plantCd_()*leafAreaDensity_()*mag(U);
54}
55
56
57// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58
60(
61 const word& sourceName,
62 const word& modelType,
63 const dictionary& dict,
64 const fvMesh& mesh
65)
66:
67 fv::cellSetOption(sourceName, modelType, dict, mesh),
68 isEpsilon_(true),
69 rhoName_(coeffs_.getOrDefault<word>("rho", "rho")),
70 Cmu_(Zero),
71 C1_(Zero),
72 C2_(Zero),
73 plantCd_
74 (
76 (
77 "plantCd",
78 mesh.time().timeName(),
79 mesh,
80 IOobject::MUST_READ,
81 IOobject::AUTO_WRITE
82 ),
83 mesh
84 ),
85 leafAreaDensity_
86 (
88 (
89 "leafAreaDensity",
90 mesh.time().timeName(),
91 mesh,
92 IOobject::MUST_READ,
93 IOobject::AUTO_WRITE
94 ),
95 mesh
96 )
97{
98 const auto* turbPtr =
100 (
102 );
103
104 if (!turbPtr)
105 {
107 << "Unable to find a turbulence model."
108 << abort(FatalError);
109 }
110
112
113 tmp<volScalarField> tepsilon = turbPtr->epsilon();
114 tmp<volScalarField> tomega = turbPtr->omega();
115
116 if (!tepsilon.isTmp())
117 {
118 fieldNames_[0] = tepsilon().name();
119
120 const dictionary& turbDict = turbPtr->coeffDict();
121 Cmu_.read("Cmu", turbDict);
122 C1_.read("C1", turbDict);
123 C2_.read("C2", turbDict);
124 }
125 else if (!tomega.isTmp())
126 {
127 isEpsilon_ = false;
128 fieldNames_[0] = tomega().name();
129
130 const dictionary& turbDict = turbPtr->coeffDict();
131 Cmu_.read("betaStar", turbDict);
132 }
133 else
134 {
136 << "Unable to find neither epsilon nor omega field." << nl
137 << "atmPlantCanopyTurbSource needs either epsilon or omega field."
138 << abort(FatalError);
139 }
140
142
143 Log << " Applying atmPlantCanopyTurbSource to: " << fieldNames_[0]
144 << endl;
145}
146
147
148// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
149
151(
152 fvMatrix<scalar>& eqn,
153 const label fieldi
154)
155{
156 if (isEpsilon_)
157 {
158 atmPlantCanopyTurbSourceEpsilon
159 (
162 eqn,
163 fieldi
164 );
165 }
166 else
167 {
168 atmPlantCanopyTurbSourceOmega
169 (
172 eqn,
173 fieldi
174 );
175 }
176}
177
178
180(
181 const volScalarField& rho,
182 fvMatrix<scalar>& eqn,
183 const label fieldi
184)
185{
186 if (isEpsilon_)
187 {
188 atmPlantCanopyTurbSourceEpsilon(geometricOneField(), rho, eqn, fieldi);
189 }
190 else
191 {
192 atmPlantCanopyTurbSourceOmega(geometricOneField(), rho, eqn, fieldi);
193 }
194}
195
196
198(
199 const volScalarField& alpha,
200 const volScalarField& rho,
201 fvMatrix<scalar>& eqn,
202 const label fieldi
203)
204{
205 if (isEpsilon_)
206 {
207 atmPlantCanopyTurbSourceEpsilon(alpha, rho, eqn, fieldi);
208 }
209 else
210 {
211 atmPlantCanopyTurbSourceOmega(alpha, rho, eqn, fieldi);
212 }
213}
214
215
216// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool read(const dictionary &dict)
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Applies sources on either epsilon or omega to incorporate effects of plant canopy for atmospheric bou...
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Intermediate abstract class for handling cell-set options for the derived fvOptions.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:127
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: fvOption.H:148
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:48
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
A class for managing temporary objects.
Definition: tmp.H:65
bool isTmp() const noexcept
Identical to is_pointer()
Definition: tmp.H:295
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList fv(nPoints)
volScalarField & alpha
dictionary dict