radialActuationDiskSourceTemplates.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 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
30 #include "volFields.H"
31 #include "fvMatrix.H"
32 #include "fvm.H"
33 #include "mathematicalConstants.H"
34 
35 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
36 
37 template<class RhoFieldType>
38 void Foam::fv::radialActuationDiskSource::
39 addRadialActuationDiskAxialInertialResistance
40 (
41  vectorField& Usource,
42  const labelList& cells,
43  const scalarField& Vcells,
44  const RhoFieldType& rho,
45  const vectorField& U
46 )
47 {
48  scalarField Tr(cells.size());
49 
50  tensor E(Zero);
51  const vector diskDir = this->diskDir();
52  E.diag(diskDir);
53 
54  const Field<vector> zoneCellCentres(mesh().cellCentres(), cells);
55  const Field<scalar> zoneCellVolumes(mesh().cellVolumes(), cells);
56 
57  const vector avgCentre = gSum(zoneCellVolumes*zoneCellCentres)/V();
58  const scalar maxR = gMax(mag(zoneCellCentres - avgCentre));
59 
60  const scalar intCoeffs =
61  radialCoeffs_[0]
62  + radialCoeffs_[1]*sqr(maxR)/2.0
63  + radialCoeffs_[2]*pow4(maxR)/3.0;
64 
65  if (mag(intCoeffs) < VSMALL)
66  {
68  << "Radial distribution coefficients lead to zero polynomial." << nl
69  << "radialCoeffs = " << radialCoeffs_
70  << exit(FatalError);
71  }
72 
73  // Compute upstream U and rho, spatial-averaged over monitor-region
74  vector Uref(Zero);
75  scalar rhoRef = 0.0;
76  label szMonitorCells = monitorCells_.size();
77 
78  for (const label celli : monitorCells_)
79  {
80  Uref += U[celli];
81  rhoRef = rhoRef + rho[celli];
82  }
83  reduce(Uref, sumOp<vector>());
84  reduce(rhoRef, sumOp<scalar>());
85  reduce(szMonitorCells, sumOp<label>());
86 
87  if (szMonitorCells == 0)
88  {
90  << "No cell is available for incoming velocity monitoring."
91  << exit(FatalError);
92  }
93 
94  Uref /= szMonitorCells;
95  rhoRef /= szMonitorCells;
96 
97  const scalar Ct = sink_*UvsCtPtr_->value(mag(Uref));
98  const scalar Cp = sink_*UvsCpPtr_->value(mag(Uref));
99 
100  if (Cp <= VSMALL || Ct <= VSMALL)
101  {
103  << "Cp and Ct must be greater than zero." << nl
104  << "Cp = " << Cp << ", Ct = " << Ct
105  << exit(FatalError);
106  }
107 
108  const scalar a = 1.0 - Cp/Ct;
109  const scalar T = 2.0*rhoRef*diskArea_*mag(Uref)*a*(1.0 - a);
110 
111  forAll(cells, i)
112  {
113  const scalar r2 = magSqr(mesh().cellCentres()[cells[i]] - avgCentre);
114 
115  Tr[i] =
116  T
117  *(radialCoeffs_[0] + radialCoeffs_[1]*r2 + radialCoeffs_[2]*sqr(r2))
118  /intCoeffs;
119 
120  Usource[cells[i]] += ((Vcells[cells[i]]/V_)*Tr[i]*E) & Uref;
121  }
122 
123  if
124  (
127  )
128  {
129  Ostream& os = file();
131 
132  os << Uref << tab << Cp << tab << Ct << tab << a << tab << T << tab
133  << endl;
134  }
135 
136  if (debug)
137  {
138  Info<< "Source name: " << name() << nl
139  << "Average centre: " << avgCentre << nl
140  << "Maximum radius: " << maxR << endl;
141  }
142 }
143 
144 
145 // ************************************************************************* //
label sink_
Flag for body forces to act as a source (false) or a sink (true)
scalar V_
Sum of cell volumes.
virtual OFstream & file()
Return access to the file (if only 1)
Definition: writeFile.C:264
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition: fvOptionI.H:30
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...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
scalar writeFileEnd_
End time for file output.
const labelList & cells() const noexcept
Return const access to the cell selection.
autoPtr< Function1< scalar > > UvsCpPtr_
Velocity vs power coefficients.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Tensor< scalar > tensor
Definition: symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
autoPtr< Function1< scalar > > UvsCtPtr_
Velocity vs thrust coefficients.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
scalar V() const noexcept
Return const access to the total cell volume.
Type gSum(const FieldField< Field, Type > &f)
const cellShapeList & cells
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
vector diskDir() const
Normal disk direction, evaluated at timeOutputValue.
const word & name() const noexcept
Return const access to the source name.
Definition: fvOptionI.H:24
Vector< scalar > vector
Definition: vector.H:57
const volScalarField & Cp
Definition: EEqn.H:7
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
virtual void writeCurrentTime(Ostream &os) const
Write the current time to stream.
Definition: writeFile.C:349
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
U
Definition: pEqn.H:72
scalar timeOutputValue() const
Return the current user-time value. (ie, after applying any timeToUserTime() conversion) ...
Definition: TimeStateI.H:24
dimensionedScalar pow4(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
labelList monitorCells_
Set of cells whereat the incoming velocity is monitored.
scalar diskArea_
Actuator disk planar surface area [m2].
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
scalar writeFileStart_
Start time for file output.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127