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 -------------------------------------------------------------------------------
12 License
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"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace fv
40 {
41  defineTypeNameAndDebug(actuationDiskSource, 0);
42  addToRunTimeSelectionTable(option, actuationDiskSource, dictionary);
43 }
44 }
45 
46 
47 const Foam::Enum
48 <
50 >
52 ({
53  { forceMethodType::FROUDE, "Froude" },
54  { forceMethodType::VARIABLE_SCALING, "variableScaling" },
55 });
56 
57 
58 const Foam::Enum
59 <
61 >
63 ({
64  { monitorMethodType::POINTS, "points" },
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 
79  if (forceMethod_ == forceMethodType::FROUDE)
80  {
81  writeFile::writeCommented(os, "a");
82  writeFile::writeCommented(os, "T");
83  }
84  else if (forceMethod_ == forceMethodType::VARIABLE_SCALING)
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  os << endl;
94 }
95 
96 
97 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
98 
99 void Foam::fv::actuationDiskSource::setMonitorCells(const dictionary& dict)
100 {
101  switch (monitorMethod_)
102  {
103  case monitorMethodType::POINTS:
104  {
105  Info<< " - selecting cells using points" << endl;
106 
107  labelHashSet selectedCells;
108 
109  List<point> monitorPoints;
110 
111  const dictionary* coeffsDictPtr = dict.findDict("monitorCoeffs");
112  if (coeffsDictPtr)
113  {
114  coeffsDictPtr->readIfPresent("points", monitorPoints);
115  }
116  else
117  {
118  monitorPoints.resize(1);
119  dict.readEntry("upstreamPoint", monitorPoints.first());
120  }
121 
122  for (const point& p : monitorPoints)
123  {
124  const label celli = mesh_.findCell(p);
125 
126  const bool found = (celli >= 0);
127 
128  if (found)
129  {
130  selectedCells.insert(celli);
131  }
132 
133  if (!returnReduceOr(found))
134  {
136  << "No owner cell found for point "
137  << p << endl;
138  }
139  }
140 
141  monitorCells_ = selectedCells.sortedToc();
142  break;
143  }
144  case monitorMethodType::CELLSET:
145  {
146  Info<< " - selecting cells using cellSet "
147  << zoneName() << endl;
148 
149  monitorCells_ = cellSet(mesh_, zoneName()).sortedToc();
150  break;
151  }
152  default:
153  {
155  << "Unknown type for monitoring of incoming velocity"
156  << monitorMethodTypeNames[monitorMethod_]
157  << ". Valid monitor method types : "
158  << monitorMethodTypeNames
159  << exit(FatalError);
160  }
161  }
162 }
163 
164 
165 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
166 
168 (
169  const word& name,
170  const word& modelType,
171  const dictionary& dict,
172  const fvMesh& mesh
173 )
174 :
175  fv::cellSetOption(name, modelType, dict, mesh),
176  writeFile(mesh, name, modelType, coeffs_),
177  forceMethod_
178  (
179  forceMethodTypeNames.getOrDefault
180  (
181  "variant",
182  coeffs_,
183  forceMethodType::FROUDE
184  )
185  ),
186  monitorMethod_
187  (
188  monitorMethodTypeNames.getOrDefault
189  (
190  "monitorMethod",
191  coeffs_,
192  monitorMethodType::POINTS
193  )
194  ),
195  sink_
196  (
197  coeffs_.getOrDefault<bool>("sink", true)
198  ? 1
199  : -1
200  ),
201  writeFileStart_(coeffs_.getOrDefault<scalar>("writeFileStart", 0)),
202  writeFileEnd_(coeffs_.getOrDefault<scalar>("writeFileEnd", VGREAT)),
203  diskArea_
204  (
205  coeffs_.getCheck<scalar>
206  (
207  "diskArea",
208  scalarMinMax::ge(VSMALL)
209  )
210  ),
211  diskDir_
212  (
213  coeffs_.getCheck<vector>
214  (
215  "diskDir",
216  [&](const vector& vec){ return mag(vec) > VSMALL; }
217  ).normalise()
218  ),
219  UvsCpPtr_(Function1<scalar>::New("Cp", coeffs_, &mesh)),
220  UvsCtPtr_(Function1<scalar>::New("Ct", coeffs_, &mesh)),
221  monitorCells_()
222 {
223  setMonitorCells(coeffs_);
224 
225  fieldNames_.resize(1, "U");
226 
228 
229  Info<< " - creating actuation disk zone: " << this->name() << endl;
230 
231  Info<< " - force computation method: "
232  << forceMethodTypeNames[forceMethod_] << endl;
233 
234  writeFileHeader(file());
235 }
236 
237 
238 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
239 
241 (
242  fvMatrix<vector>& eqn,
243  const label fieldi
244 )
245 {
246  if (V() > VSMALL)
247  {
248  calc(geometricOneField(), geometricOneField(), eqn);
249  }
250 }
251 
252 
254 (
255  const volScalarField& rho,
256  fvMatrix<vector>& eqn,
257  const label fieldi
258 )
259 {
260  if (V() > VSMALL)
261  {
262  calc(geometricOneField(), rho, eqn);
263  }
264 }
265 
266 
268 (
269  const volScalarField& alpha,
270  const volScalarField& rho,
271  fvMatrix<vector>& eqn,
272  const label fieldi
273 )
274 {
275  if (V() > VSMALL)
276  {
277  calc(alpha, rho, eqn);
278  }
279 }
280 
281 
283 {
285  {
286  dict.readIfPresent("sink", sink_);
287  dict.readIfPresent("writeFileStart", writeFileStart_);
288  dict.readIfPresent("writeFileEnd", writeFileEnd_);
289  dict.readIfPresent("diskArea", diskArea_);
290  if (diskArea_ < VSMALL)
291  {
293  << "Actuator disk has zero area: "
294  << "diskArea = " << diskArea_
295  << exit(FatalIOError);
296  }
297 
298  dict.readIfPresent("diskDir", diskDir_);
299  diskDir_.normalise();
300  if (mag(diskDir_) < VSMALL)
301  {
303  << "Actuator disk surface-normal vector is zero: "
304  << "diskDir = " << diskDir_
305  << exit(FatalIOError);
306  }
307 
308  return true;
309  }
310 
311  return false;
312 }
313 
314 
315 // ************************************************************************* //
dictionary dict
static void writeHeader(Ostream &os, const word &fieldName)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
forceMethodType
Options for the force computation method types.
static const Enum< forceMethodType > forceMethodTypeNames
Names for forceMethodType.
virtual void writeFileHeader(Ostream &os)
Output file header information.
virtual bool read(const dictionary &dict)
Read source dictionary.
Macros for easy insertion into run-time selection tables.
enum forceMethodType forceMethod_
The type of the force computation method.
actuationDiskSource()=delete
No default construct.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
dynamicFvMesh & mesh
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Source term to momentum equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
static const Enum< monitorMethodType > monitorMethodTypeNames
Names for monitorMethodType.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
OBJstream os(runTime.globalPath()/outputName)
virtual bool read(const dictionary &dict)
Read dictionary.
vector point
Point is a vector.
Definition: point.H:37
#define WarningInFunction
Report a warning using Foam::Warning.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
messageStream Info
Information stream (stdout output on master, null elsewhere)
monitorMethodType
Options for the incoming velocity monitoring method types.
volScalarField & p
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
bool found
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...