Loading...
Searching...
No Matches
actuationDiskSource.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2020 ENERCON GmbH
10 Copyright (C) 2018-2022 OpenCFD Ltd
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "actuationDiskSource.H"
31#include "geometricOneField.H"
32#include "cellSet.H"
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace fv
40{
43}
44}
45
46
47const Foam::Enum
48<
50>
52({
53 { forceMethodType::FROUDE, "Froude" },
54 { forceMethodType::VARIABLE_SCALING, "variableScaling" },
55});
56
57
58const Foam::Enum
59<
61>
63({
65 { monitorMethodType::CELLSET, "cellSet" },
66});
67
68
69// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
70
72{
73 writeFile::writeHeader(os, "Actuation disk source");
74 writeFile::writeCommented(os, "Time");
75 writeFile::writeCommented(os, "Uref");
76 writeFile::writeCommented(os, "Cp");
77 writeFile::writeCommented(os, "Ct");
78
80 {
81 writeFile::writeCommented(os, "a");
82 writeFile::writeCommented(os, "T");
83 }
85 {
86 writeFile::writeCommented(os, "Udisk");
87 writeFile::writeCommented(os, "CpStar");
88 writeFile::writeCommented(os, "CtStar");
89 writeFile::writeCommented(os, "T");
90 writeFile::writeCommented(os, "P");
91 }
92
93 writeFile::writeCommented(os, "diskDir");
94
95 os << endl;
96}
97
98
99// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
100
101void Foam::fv::actuationDiskSource::setMonitorCells(const dictionary& dict)
102{
103 switch (monitorMethod_)
104 {
105 case monitorMethodType::POINTS:
106 {
107 Info<< " - selecting cells using points" << endl;
108
109 labelHashSet selectedCells;
110
111 List<point> monitorPoints;
112
113 const dictionary* coeffsDictPtr = dict.findDict("monitorCoeffs");
114 if (coeffsDictPtr)
115 {
116 coeffsDictPtr->readIfPresent("points", monitorPoints);
117 }
118 else
119 {
120 monitorPoints.resize(1);
121 dict.readEntry("upstreamPoint", monitorPoints.first());
122 }
123
124 for (const point& p : monitorPoints)
125 {
126 const label celli = mesh_.findCell(p);
127
128 const bool found = (celli >= 0);
129
130 if (found)
131 {
132 selectedCells.insert(celli);
133 }
134
135 if (!returnReduceOr(found))
136 {
138 << "No owner cell found for point "
139 << p << endl;
140 }
141 }
142
143 monitorCells_ = selectedCells.sortedToc();
144 break;
145 }
146 case monitorMethodType::CELLSET:
147 {
148 Info<< " - selecting cells using cellSet "
149 << zoneName() << endl;
150
151 monitorCells_ = cellSet(mesh_, zoneName()).sortedToc();
152 break;
153 }
154 default:
155 {
157 << "Unknown type for monitoring of incoming velocity"
158 << monitorMethodTypeNames[monitorMethod_]
159 << ". Valid monitor method types : "
160 << monitorMethodTypeNames
161 << exit(FatalError);
163 }
164}
165
166
167// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
168
170(
171 const word& name,
172 const word& modelType,
173 const dictionary& dict,
174 const fvMesh& mesh
175)
176:
177 fv::cellSetOption(name, modelType, dict, mesh),
178 writeFile(mesh, name, modelType, coeffs_),
179 forceMethod_
180 (
181 forceMethodTypeNames.getOrDefault
182 (
183 "variant",
184 coeffs_,
185 forceMethodType::FROUDE
186 )
187 ),
188 monitorMethod_
189 (
190 monitorMethodTypeNames.getOrDefault
191 (
192 "monitorMethod",
193 coeffs_,
194 monitorMethodType::POINTS
195 )
196 ),
197 sink_
198 (
199 coeffs_.getOrDefault<bool>("sink", true)
200 ? 1
201 : -1
202 ),
203 writeFileStart_(coeffs_.getOrDefault<scalar>("writeFileStart", 0)),
204 writeFileEnd_(coeffs_.getOrDefault<scalar>("writeFileEnd", VGREAT)),
205 diskArea_
206 (
207 coeffs_.getCheck<scalar>
208 (
209 "diskArea",
210 scalarMinMax::ge(VSMALL)
211 )
212 ),
213 diskDir_(Function1<vector>::New("diskDir", coeffs_, &mesh)),
214 UvsCpPtr_(Function1<scalar>::New("Cp", coeffs_, &mesh)),
215 UvsCtPtr_(Function1<scalar>::New("Ct", coeffs_, &mesh)),
216 monitorCells_()
217{
218 setMonitorCells(coeffs_);
219
220 fieldNames_.resize(1, "U");
221
223
224 Info<< " - creating actuation disk zone: " << this->name() << endl;
225
226 Info<< " - force computation method: "
230}
231
232
233// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
234
236{
237 const scalar t = mesh_.time().timeOutputValue();
238 const vector dir(diskDir_->value(t));
239 const scalar magDir = Foam::mag(dir);
240
241 if (magDir < SMALL)
242 {
244 << "Actuator disk surface-normal vector is zero: " << dir
245 << " at time=" << t
246 << exit(FatalError);
248
249 // normalised:
250 return dir/magDir;
251}
252
253
255(
256 fvMatrix<vector>& eqn,
257 const label fieldi
258)
259{
260 if (V() > VSMALL)
261 {
262 calc(geometricOneField(), geometricOneField(), eqn);
263 }
264}
265
266
268(
269 const volScalarField& rho,
270 fvMatrix<vector>& eqn,
271 const label fieldi
272)
273{
274 if (V() > VSMALL)
275 {
276 calc(geometricOneField(), rho, eqn);
277 }
278}
279
280
282(
283 const volScalarField& alpha,
284 const volScalarField& rho,
285 fvMatrix<vector>& eqn,
286 const label fieldi
287)
288{
289 if (V() > VSMALL)
290 {
291 calc(alpha, rho, eqn);
292 }
293}
294
295
297{
298 if (fv::cellSetOption::read(dict) && writeFile::read(dict))
299 {
300 dict.readIfPresent("sink", sink_);
301 dict.readIfPresent("writeFileStart", writeFileStart_);
302 dict.readIfPresent("writeFileEnd", writeFileEnd_);
303 dict.readIfPresent("diskArea", diskArea_);
304 if (diskArea_ < VSMALL)
305 {
307 << "Actuator disk has zero area: "
308 << "diskArea = " << diskArea_
309 << exit(FatalIOError);
310 }
311
312 // TBD: runTime re-reading of "diskDir" ?
313
314 return true;
315 }
316 return false;
317}
318
319
320// ************************************************************************* //
bool found
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition HashTable.C:156
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
T & first()
Access first element of the list, position [0].
Definition UList.H:957
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and it is a dictionary) otherwise return nullptr...
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
wordList sortedToc() const
Return the sorted table of contents.
Definition dictionary.C:601
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
Definition writeFile.C:200
virtual OFstream & file()
Return access to the file (if only 1).
Definition writeFile.C:270
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
Applies sources on velocity (i.e. U) within a specified region to enable actuator disk models for aer...
actuationDiskSource()=delete
No default construct.
static const Enum< monitorMethodType > monitorMethodTypeNames
Names for monitorMethodType.
autoPtr< Function1< scalar > > UvsCtPtr_
Velocity vs thrust coefficients.
vector diskDir() const
Normal disk direction, evaluated at timeOutputValue.
static const Enum< forceMethodType > forceMethodTypeNames
Names for forceMethodType.
monitorMethodType
Options for the incoming velocity monitoring method types.
scalar writeFileEnd_
End time for file output.
virtual bool read(const dictionary &dict)
Read dictionary.
autoPtr< Function1< vector > > diskDir_
Surface-normal vector of the actuator disk pointing downstream as a function of time.
scalar writeFileStart_
Start time for file output.
enum monitorMethodType monitorMethod_
The type of incoming velocity monitoring method.
enum forceMethodType forceMethod_
The type of the force computation method.
forceMethodType
Options for the force computation method types.
@ FROUDE
"Froude's ideal actuator disk method"
@ VARIABLE_SCALING
"Variable-scaling actuator disk method"
labelList monitorCells_
Set of cells whereat the incoming velocity is monitored.
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Source term to momentum equation.
virtual void writeFileHeader(Ostream &os)
Output file header information.
autoPtr< Function1< scalar > > UvsCpPtr_
Velocity vs power coefficients.
label sink_
Flag for body forces to act as a source (false) or a sink (true).
scalar diskArea_
Actuator disk planar surface area [m2].
scalar V() const noexcept
Return const access to the total cell volume.
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 word & name() const noexcept
Return const access to the source name.
Definition fvOptionI.H:24
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
static autoPtr< option > New(const word &name, const dictionary &dict, const fvMesh &mesh)
Return a reference to the selected fvOption model.
Definition fvOption.C:78
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition fvOption.H:157
dictionary coeffs_
Dictionary containing source coefficients.
Definition fvOption.H:152
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...
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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
volScalarField & p
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
auto & name
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for finite-volume.
Namespace for OpenFOAM.
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
vector point
Point is a vector.
Definition point.H:37
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
labelList fv(nPoints)
volScalarField & alpha
dictionary dict