Loading...
Searching...
No Matches
actuationDiskSourceTemplates.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-2020 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 "fvMesh.H"
32#include "fvMatrix.H"
33#include "volFields.H"
34
35// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36
37template<class AlphaFieldType, class RhoFieldType>
38void Foam::fv::actuationDiskSource::calc
39(
40 const AlphaFieldType& alpha,
41 const RhoFieldType& rho,
42 fvMatrix<vector>& eqn
43)
44{
45 switch (forceMethod_)
46 {
48 {
49 calcFroudeMethod(alpha, rho, eqn);
50 break;
51 }
52
54 {
55 calcVariableScalingMethod(alpha, rho, eqn);
56 break;
57 }
58
59 default:
60 break;
61 }
62}
63
64
65template<class AlphaFieldType, class RhoFieldType>
66void Foam::fv::actuationDiskSource::calcFroudeMethod
67(
68 const AlphaFieldType& alpha,
69 const RhoFieldType& rho,
70 fvMatrix<vector>& eqn
71)
72{
73 const vectorField& U = eqn.psi();
74 vectorField& Usource = eqn.source();
75 const scalarField& cellsV = mesh_.V();
76
77 // Compute upstream U and rho, spatial-averaged over monitor-region
78 vector Uref(Zero);
79 scalar rhoRef = 0.0;
80 label szMonitorCells = monitorCells_.size();
81
82 for (const label celli : monitorCells_)
83 {
84 Uref += U[celli];
85 rhoRef = rhoRef + rho[celli];
86 }
87 reduce(Uref, sumOp<vector>());
88 reduce(rhoRef, sumOp<scalar>());
89 reduce(szMonitorCells, sumOp<label>());
90
91 if (szMonitorCells == 0)
92 {
94 << "No cell is available for incoming velocity monitoring."
95 << exit(FatalError);
96 }
97
98 Uref /= szMonitorCells;
99 rhoRef /= szMonitorCells;
100
101 const scalar Ct = sink_*UvsCtPtr_->value(mag(Uref));
102 const scalar Cp = sink_*UvsCpPtr_->value(mag(Uref));
103
104 if (Cp <= VSMALL || Ct <= VSMALL)
105 {
107 << "Cp and Ct must be greater than zero." << nl
108 << "Cp = " << Cp << ", Ct = " << Ct
109 << exit(FatalError);
110 }
111
112 // (BJSB:Eq. 3.9)
113 const vector diskDir = this->diskDir();
114 const scalar a = 1.0 - Cp/Ct;
115 const scalar T = 2.0*rhoRef*diskArea_*magSqr(Uref & diskDir)*a*(1 - a);
116
117 for (const label celli : cells_)
118 {
119 Usource[celli] += ((cellsV[celli]/V())*T)*diskDir;
120 }
121
122 if
123 (
124 mesh_.time().timeOutputValue() >= writeFileStart_
125 && mesh_.time().timeOutputValue() <= writeFileEnd_
126 )
127 {
128 Ostream& os = file();
129 writeCurrentTime(os);
130
131 os << Uref << tab << Cp << tab << Ct << tab << a << tab << T
132 << tab << diskDir << endl;
133 }
134}
135
136
137template<class AlphaFieldType, class RhoFieldType>
138void Foam::fv::actuationDiskSource::calcVariableScalingMethod
139(
140 const AlphaFieldType& alpha,
141 const RhoFieldType& rho,
142 fvMatrix<vector>& eqn
143)
144{
145 const vectorField& U = eqn.psi();
146 vectorField& Usource = eqn.source();
147 const scalarField& cellsV = mesh_.V();
148
149 // Monitor and average monitor-region U and rho
150 vector Uref(Zero);
151 scalar rhoRef = 0.0;
152 label szMonitorCells = monitorCells_.size();
153
154 for (const label celli : monitorCells_)
155 {
156 Uref += U[celli];
157 rhoRef = rhoRef + rho[celli];
158 }
159 reduce(Uref, sumOp<vector>());
160 reduce(rhoRef, sumOp<scalar>());
161 reduce(szMonitorCells, sumOp<label>());
162
163 if (szMonitorCells == 0)
164 {
166 << "No cell is available for incoming velocity monitoring."
167 << exit(FatalError);
168 }
169
170 Uref /= szMonitorCells;
171 const scalar magUref = mag(Uref);
172 rhoRef /= szMonitorCells;
173
174 // Monitor and average U and rho on actuator disk
175 vector Udisk(Zero);
176 scalar rhoDisk = 0.0;
177 scalar totalV = 0.0;
178
179 for (const label celli : cells_)
180 {
181 Udisk += U[celli]*cellsV[celli];
182 rhoDisk += rho[celli]*cellsV[celli];
183 totalV += cellsV[celli];
184 }
185 reduce(Udisk, sumOp<vector>());
186 reduce(rhoDisk, sumOp<scalar>());
187 reduce(totalV, sumOp<scalar>());
188
189 if (totalV < SMALL)
190 {
192 << "No cell in the actuator disk."
193 << exit(FatalError);
194 }
195
196 Udisk /= totalV;
197 const scalar magUdisk = mag(Udisk);
198 rhoDisk /= totalV;
199
200 if (mag(Udisk) < SMALL)
201 {
203 << "Velocity spatial-averaged on actuator disk is zero." << nl
204 << "Please check if the initial U field is zero."
205 << exit(FatalError);
206 }
207
208 // Interpolated thrust/power coeffs from power/thrust curves
209 const scalar Ct = sink_*UvsCtPtr_->value(magUref);
210 const scalar Cp = sink_*UvsCpPtr_->value(magUref);
211
212 if (Cp <= VSMALL || Ct <= VSMALL)
213 {
215 << "Cp and Ct must be greater than zero." << nl
216 << "Cp = " << Cp << ", Ct = " << Ct
217 << exit(FatalError);
218 }
219
220 // Calibrated thrust/power coeffs from power/thrust curves (LSRMTK:Eq. 6)
221 const scalar CtStar = Ct*sqr(magUref/magUdisk);
222 const scalar CpStar = Cp*pow3(magUref/magUdisk);
223
224 // Compute calibrated thrust/power (LSRMTK:Eq. 5)
225 const vector diskDir = this->diskDir();
226 const scalar T = 0.5*rhoRef*diskArea_*magSqr(Udisk & diskDir)*CtStar;
227 const scalar P = 0.5*rhoRef*diskArea_*pow3(mag(Udisk & diskDir))*CpStar;
228
229 for (const label celli : cells_)
230 {
231 Usource[celli] += (cellsV[celli]/totalV*T)*diskDir;
232 }
233
234 if
235 (
236 mesh_.time().timeOutputValue() >= writeFileStart_
237 && mesh_.time().timeOutputValue() <= writeFileEnd_
238 )
239 {
240 Ostream& os = file();
241 writeCurrentTime(os);
242
243 os << Uref << tab << Cp << tab << Ct << tab
244 << Udisk << tab << CpStar << tab << CtStar << tab << T << tab << P
245 << tab << diskDir << endl;
246 }
247}
248
249
250// ************************************************************************* //
enum forceMethodType forceMethod_
The type of the force computation method.
@ FROUDE
"Froude's ideal actuator disk method"
@ VARIABLE_SCALING
"Variable-scaling actuator disk method"
U
Definition pEqn.H:72
const volScalarField & T
const volScalarField & Cp
Definition EEqn.H:7
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const expr V(m.psi().mesh().V())
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
Field< vector > vectorField
Specialisation of Field<T> for vector.
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...
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
volScalarField & alpha