Loading...
Searching...
No Matches
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-2024 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
32// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace fv
37{
40}
41}
42
43
44// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
45
46Foam::volScalarField& Foam::fv::atmPlantCanopyTurbSource::getOrReadField
47(
48 const word& fieldName
49) const
50{
51 auto* ptr = mesh_.getObjectPtr<volScalarField>(fieldName);
52
53 if (!ptr)
54 {
55 ptr = new volScalarField
56 (
58 (
59 fieldName,
61 mesh_.thisDb(),
65 ),
66 mesh_
67 );
69 }
70
71 return *ptr;
72}
73
74
75Foam::tmp<Foam::volScalarField::Internal>
76Foam::fv::atmPlantCanopyTurbSource::calcPlantCanopyTerm
77(
79) const
80{
81 const volScalarField& Cd = getOrReadField(CdName_);
82 const volScalarField& LAD = getOrReadField(LADname_);
83
84 // (SP:Eq. 42)
85 return 12.0*Foam::sqrt(Cmu_)*Cd()*LAD()*mag(U);
86}
87
88
89// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90
92(
93 const word& sourceName,
94 const word& modelType,
95 const dictionary& dict,
96 const fvMesh& mesh
97)
98:
99 fv::cellSetOption(sourceName, modelType, dict, mesh),
100 isEpsilon_(true),
101 Cmu_(Zero),
102 C1_(Zero),
103 C2_(Zero),
104 CdName_(),
105 LADname_()
106{
107 read(dict);
108
109 const auto* turbPtr =
111 (
113 );
114
115 if (!turbPtr)
116 {
118 << "Unable to find a turbulence model."
119 << abort(FatalError);
120 }
121
123
124 tmp<volScalarField> tepsilon = turbPtr->epsilon();
125 tmp<volScalarField> tomega = turbPtr->omega();
126
127 if (tepsilon.is_reference())
128 {
129 isEpsilon_ = true;
130 fieldNames_[0] = tepsilon().name();
131
132 const dictionary& turbDict = turbPtr->coeffDict();
133 Cmu_.read("Cmu", turbDict);
134 C1_.read("C1", turbDict);
135 C2_.read("C2", turbDict);
136 }
137 else if (tomega.is_reference())
138 {
139 isEpsilon_ = false;
140 fieldNames_[0] = tomega().name();
141
142 const dictionary& turbDict = turbPtr->coeffDict();
143 Cmu_.read("betaStar", turbDict);
144 }
145 else
146 {
148 << "Needs either epsilon or omega field."
149 << abort(FatalError);
150 }
151
153
154 Log << " Applying atmPlantCanopyTurbSource to: " << fieldNames_[0]
155 << endl;
156}
157
158
159// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
160
162(
163 fvMatrix<scalar>& eqn,
164 const label fieldi
165)
166{
167 if (isEpsilon_)
168 {
169 atmPlantCanopyTurbSourceEpsilon
170 (
173 eqn,
174 fieldi
175 );
176 }
177 else
178 {
179 atmPlantCanopyTurbSourceOmega
180 (
183 eqn,
184 fieldi
185 );
186 }
187}
188
189
191(
192 const volScalarField& rho,
193 fvMatrix<scalar>& eqn,
194 const label fieldi
195)
196{
197 if (isEpsilon_)
198 {
199 atmPlantCanopyTurbSourceEpsilon(geometricOneField(), rho, eqn, fieldi);
200 }
201 else
202 {
203 atmPlantCanopyTurbSourceOmega(geometricOneField(), rho, eqn, fieldi);
204 }
205}
206
207
209(
210 const volScalarField& alpha,
211 const volScalarField& rho,
212 fvMatrix<scalar>& eqn,
213 const label fieldi
214)
215{
216 if (isEpsilon_)
217 {
218 atmPlantCanopyTurbSourceEpsilon(alpha, rho, eqn, fieldi);
219 }
220 else
221 {
222 atmPlantCanopyTurbSourceOmega(alpha, rho, eqn, fieldi);
223 }
224}
225
226
228{
230 {
231 return false;
232 }
233
234 CdName_ = dict.getOrDefault<word>("Cd", "Cd");
235 LADname_ = dict.getOrDefault<word>("LAD", "LAD");
236
237 (void) getOrReadField(CdName_);
238 (void) getOrReadField(LADname_);
239
240 return true;
241}
242
243
244// ************************************************************************* //
#define Log
Definition PDRblock.C:28
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
DimensionedField< vector, volMesh > Internal
@ REGISTER
Request registration (bool: true).
@ MUST_READ
Reading required.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool read(Istream &is)
Read dictionary from Istream (discards the header). Reads entries until EOF or when the first token i...
bool read(const dictionary &dict)
Update the value of dimensioned<Type>, lookup in dictionary with the name().
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition fvMesh.H:376
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
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)
Add explicit contribution to epsilon or omega equation for incompressible flow computations.
virtual bool read(const dictionary &dict)
Read source dictionary.
atmPlantCanopyTurbSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
Intermediate abstract class for handling cell-set options for the derived fvOptions.
virtual bool read(const dictionary &dict)
Read source dictionary.
cellSetOption(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition fvOption.H:124
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition fvOption.H:157
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition fvOption.C:41
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition fvOptionI.H:30
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.
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
bool store()
Register object with its registry and transfer ownership to the registry.
A class for managing temporary objects.
Definition tmp.H:75
bool is_reference() const noexcept
True if this is a reference (not a pointer).
Definition tmpI.H:206
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:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for finite-volume.
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
volScalarField & alpha
dictionary dict