radialActuationDisk.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "radialActuationDisk.H"
27 #include "volFields.H"
28 #include "fvMatrix.H"
29 #include "geometricOneField.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace fv
37 {
41  (
42  fvModel,
44  dictionary,
45  radialActuationDiskSource,
46  "radialActuationDiskSource"
47  );
48 }
49 }
50 
51 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 void Foam::fv::radialActuationDisk::readCoeffs()
54 {
55  coeffs().lookup("coeffs") >> radialCoeffs_;
56 }
57 
58 
59 
60 template<class RhoFieldType>
61 void Foam::fv::radialActuationDisk::
62 addRadialActuationDiskAxialInertialResistance
63 (
64  vectorField& Usource,
65  const labelList& cells,
66  const scalarField& Vcells,
67  const RhoFieldType& rho,
68  const vectorField& U
69 ) const
70 {
71  scalar a = 1.0 - Cp_/Ct_;
72  scalarField Tr(cells.size());
73  const vector uniDiskDir = diskDir_/mag(diskDir_);
74 
75  tensor E(Zero);
76  E.xx() = uniDiskDir.x();
77  E.yy() = uniDiskDir.y();
78  E.zz() = uniDiskDir.z();
79 
80  const Field<vector> zoneCellCentres(mesh().cellCentres(), cells);
81  const Field<scalar> zoneCellVolumes(mesh().cellVolumes(), cells);
82 
83  const vector avgCentre = gSum(zoneCellVolumes*zoneCellCentres)/set_.V();
84  const scalar maxR = gMax(mag(zoneCellCentres - avgCentre));
85 
86  scalar intCoeffs =
87  radialCoeffs_[0]
88  + radialCoeffs_[1]*sqr(maxR)/2.0
89  + radialCoeffs_[2]*pow4(maxR)/3.0;
90 
91  vector upU = vector(vGreat, vGreat, vGreat);
92  scalar upRho = vGreat;
93  if (upstreamCellId_ != -1)
94  {
95  upU = U[upstreamCellId_];
96  upRho = rho[upstreamCellId_];
97  }
98  reduce(upU, minOp<vector>());
99  reduce(upRho, minOp<scalar>());
100 
101  scalar T = 2.0*upRho*diskArea_*mag(upU)*a*(1.0 - a);
102  forAll(cells, i)
103  {
104  scalar r2 = magSqr(mesh().cellCentres()[cells[i]] - avgCentre);
105 
106  Tr[i] =
107  T
108  *(radialCoeffs_[0] + radialCoeffs_[1]*r2 + radialCoeffs_[2]*sqr(r2))
109  /intCoeffs;
110 
111  Usource[cells[i]] += ((Vcells[cells[i]]/set_.V())*Tr[i]*E) & upU;
112  }
113 
114  if (debug)
115  {
116  Info<< "Source name: " << name() << nl
117  << "Average centre: " << avgCentre << nl
118  << "Maximum radius: " << maxR << endl;
119  }
120 }
121 
122 
123 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
124 
126 (
127  const word& name,
128  const word& modelType,
129  const fvMesh& mesh,
130  const dictionary& dict
131 )
132 :
133  actuationDisk(name, modelType, mesh, dict),
134  radialCoeffs_()
135 {
136  readCoeffs();
137 }
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
143 (
144  const volVectorField& U,
145  fvMatrix<vector>& eqn
146 ) const
147 {
148  addRadialActuationDiskAxialInertialResistance
149  (
150  eqn.source(),
151  set_.cells(),
152  mesh().V(),
154  U
155  );
156 }
157 
158 
160 (
161  const volScalarField& rho,
162  const volVectorField& U,
163  fvMatrix<vector>& eqn
164 ) const
165 {
166  addRadialActuationDiskAxialInertialResistance
167  (
168  eqn.source(),
169  set_.cells(),
170  mesh().V(),
171  rho,
172  U
173  );
174 }
175 
176 
178 {
180  {
181  readCoeffs();
182  return true;
183  }
184  else
185  {
186  return false;
187  }
188 }
189 
190 
191 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Field< Type > & source()
Definition: fvMatrix.H:307
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
Actuation disk source.
virtual bool read(const dictionary &dict)
Read dictionary.
Actuation disk source including radial thrust.
virtual void addSup(const volVectorField &U, fvMatrix< vector > &eqn) const
Source term to momentum equation.
radialActuationDisk(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
virtual bool read(const dictionary &dict)
Read dictionary.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
A class for handling words, derived from string.
Definition: word.H:62
const cellShapeList & cells
U
Definition: pEqn.H:72
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
addBackwardCompatibleToRunTimeSelectionTable(fvConstraint, fixedTemperature, dictionary, fixedTemperatureConstraint, "fixedTemperatureConstraint")
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
Type gSum(const FieldField< Field, Type > &f)
List< label > labelList
A List of labels.
Definition: labelList.H:56
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
messageStream Info
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedScalar pow4(const dimensionedScalar &ds)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Type gMax(const FieldField< Field, Type > &f)
labelList fv(nPoints)
dictionary dict