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 -------------------------------------------------------------------------------
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 "fvMesh.H"
32 #include "fvMatrix.H"
33 #include "volFields.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 template<class AlphaFieldType, class RhoFieldType>
38 void Foam::fv::actuationDiskSource::calc
39 (
40  const AlphaFieldType& alpha,
41  const RhoFieldType& rho,
42  fvMatrix<vector>& eqn
43 )
44 {
45  switch (forceMethod_)
46  {
47  case forceMethodType::FROUDE:
48  {
49  calcFroudeMethod(alpha, rho, eqn);
50  break;
51  }
52 
53  case forceMethodType::VARIABLE_SCALING:
54  {
55  calcVariableScalingMethod(alpha, rho, eqn);
56  break;
57  }
58 
59  default:
60  break;
61  }
62 }
63 
64 
65 template<class AlphaFieldType, class RhoFieldType>
66 void 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 scalar a = 1.0 - Cp/Ct;
114  const scalar T = 2.0*rhoRef*diskArea_*magSqr(Uref & diskDir_)*a*(1 - a);
115 
116  for (const label celli : cells_)
117  {
118  Usource[celli] += ((cellsV[celli]/V())*T)*diskDir_;
119  }
120 
121  if
122  (
123  mesh_.time().timeOutputValue() >= writeFileStart_
124  && mesh_.time().timeOutputValue() <= writeFileEnd_
125  )
126  {
127  Ostream& os = file();
128  writeCurrentTime(os);
129 
130  os << Uref << tab << Cp << tab << Ct << tab << a << tab << T
131  << endl;
132  }
133 }
134 
135 
136 template<class AlphaFieldType, class RhoFieldType>
137 void Foam::fv::actuationDiskSource::calcVariableScalingMethod
138 (
139  const AlphaFieldType& alpha,
140  const RhoFieldType& rho,
141  fvMatrix<vector>& eqn
142 )
143 {
144  const vectorField& U = eqn.psi();
145  vectorField& Usource = eqn.source();
146  const scalarField& cellsV = mesh_.V();
147 
148  // Monitor and average monitor-region U and rho
149  vector Uref(Zero);
150  scalar rhoRef = 0.0;
151  label szMonitorCells = monitorCells_.size();
152 
153  for (const label celli : monitorCells_)
154  {
155  Uref += U[celli];
156  rhoRef = rhoRef + rho[celli];
157  }
158  reduce(Uref, sumOp<vector>());
159  reduce(rhoRef, sumOp<scalar>());
160  reduce(szMonitorCells, sumOp<label>());
161 
162  if (szMonitorCells == 0)
163  {
165  << "No cell is available for incoming velocity monitoring."
166  << exit(FatalError);
167  }
168 
169  Uref /= szMonitorCells;
170  const scalar magUref = mag(Uref);
171  rhoRef /= szMonitorCells;
172 
173  // Monitor and average U and rho on actuator disk
174  vector Udisk(Zero);
175  scalar rhoDisk = 0.0;
176  scalar totalV = 0.0;
177 
178  for (const label celli : cells_)
179  {
180  Udisk += U[celli]*cellsV[celli];
181  rhoDisk += rho[celli]*cellsV[celli];
182  totalV += cellsV[celli];
183  }
184  reduce(Udisk, sumOp<vector>());
185  reduce(rhoDisk, sumOp<scalar>());
186  reduce(totalV, sumOp<scalar>());
187 
188  if (totalV < SMALL)
189  {
191  << "No cell in the actuator disk."
192  << exit(FatalError);
193  }
194 
195  Udisk /= totalV;
196  const scalar magUdisk = mag(Udisk);
197  rhoDisk /= totalV;
198 
199  if (mag(Udisk) < SMALL)
200  {
202  << "Velocity spatial-averaged on actuator disk is zero." << nl
203  << "Please check if the initial U field is zero."
204  << exit(FatalError);
205  }
206 
207  // Interpolated thrust/power coeffs from power/thrust curves
208  const scalar Ct = sink_*UvsCtPtr_->value(magUref);
209  const scalar Cp = sink_*UvsCpPtr_->value(magUref);
210 
211  if (Cp <= VSMALL || Ct <= VSMALL)
212  {
214  << "Cp and Ct must be greater than zero." << nl
215  << "Cp = " << Cp << ", Ct = " << Ct
216  << exit(FatalError);
217  }
218 
219  // Calibrated thrust/power coeffs from power/thrust curves (LSRMTK:Eq. 6)
220  const scalar CtStar = Ct*sqr(magUref/magUdisk);
221  const scalar CpStar = Cp*pow3(magUref/magUdisk);
222 
223  // Compute calibrated thrust/power (LSRMTK:Eq. 5)
224  const scalar T = 0.5*rhoRef*diskArea_*magSqr(Udisk & diskDir_)*CtStar;
225  const scalar P = 0.5*rhoRef*diskArea_*pow3(mag(Udisk & diskDir_))*CpStar;
226 
227  for (const label celli : cells_)
228  {
229  Usource[celli] += (cellsV[celli]/totalV*T)*diskDir_;
230  }
231 
232  if
233  (
234  mesh_.time().timeOutputValue() >= writeFileStart_
235  && mesh_.time().timeOutputValue() <= writeFileEnd_
236  )
237  {
238  Ostream& os = file();
239  writeCurrentTime(os);
240 
241  os << Uref << tab << Cp << tab << Ct << tab
242  << Udisk << tab << CpStar << tab << CtStar << tab << T << tab << P
243  << endl;
244  }
245 }
246 
247 
248 // ************************************************************************* //
volFields.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
fvMatrix.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
rho
rho
Definition: readInitialConditions.H:88
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::fv::actuationDiskSource::forceMethod_
enum forceMethodType forceMethod_
The type of the force computation method.
Definition: actuationDiskSource.H:398
Foam::FatalError
error FatalError
os
OBJstream os(runTime.globalPath()/outputName)
fvMesh.H
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::tab
constexpr char tab
Definition: Ostream.H:403
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Cp
const volScalarField & Cp
Definition: EEqn.H:7
actuationDiskSource.H