actuationDisk.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 "actuationDisk.H"
27 #include "fvMesh.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  actuationDiskSource,
46  "actuationDiskSource"
47  );
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
54 void Foam::fv::actuationDisk::readCoeffs()
55 {
56  phaseName_ = coeffs().lookupOrDefault<word>("phase", word::null);
57 
58  UName_ =
59  coeffs().lookupOrDefault<word>
60  (
61  "U",
63  );
64 
65  diskDir_ = coeffs().lookup<vector>("diskDir");
66  if (mag(diskDir_) < vSmall)
67  {
69  << "disk direction vector is approximately zero"
70  << exit(FatalIOError);
71  }
72 
73  Cp_ = coeffs().lookup<scalar>("Cp");
74  Ct_ = coeffs().lookup<scalar>("Ct");
75  if (Cp_ <= vSmall || Ct_ <= vSmall)
76  {
78  << "Cp and Ct must be greater than zero"
79  << exit(FatalIOError);
80  }
81 
82  diskArea_ = coeffs().lookup<scalar>("diskArea");
83  if (magSqr(diskArea_) <= vSmall)
84  {
86  << "diskArea is approximately zero"
87  << exit(FatalIOError);
88  }
89 
90  upstreamPoint_ = coeffs().lookup<point>("upstreamPoint");
92  if (returnReduce(upstreamCellId_, maxOp<label>()) == -1)
93  {
95  << "upstream location " << upstreamPoint_ << " not found in mesh"
96  << exit(FatalIOError);
97  }
98 }
99 
100 
101 template<class AlphaFieldType, class RhoFieldType>
102 void Foam::fv::actuationDisk::addActuationDiskAxialInertialResistance
103 (
104  vectorField& Usource,
105  const labelList& cells,
106  const scalarField& Vcells,
107  const AlphaFieldType& alpha,
108  const RhoFieldType& rho,
109  const vectorField& U
110 ) const
111 {
112  const scalar a = 1 - Cp_/Ct_;
113  const vector dHat(diskDir_/mag(diskDir_));
114 
115  scalar dHatUo(vGreat);
116  if (upstreamCellId_ != -1)
117  {
118  dHatUo = dHat & U[upstreamCellId_];
119  }
120  reduce(dHatUo, minOp<scalar>());
121 
122  const vector T = 2*diskArea_*sqr(dHatUo)*a*(1 - a)*dHat;
123 
124  forAll(cells, i)
125  {
126  Usource[cells[i]] +=
127  (alpha[cells[i]]*rho[cells[i]]*(Vcells[cells[i]]/set_.V()))*T;
128  }
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
133 
135 (
136  const word& name,
137  const word& modelType,
138  const fvMesh& mesh,
139  const dictionary& dict
140 )
141 :
142  fvModel(name, modelType, mesh, dict),
143  set_(mesh, coeffs()),
144  phaseName_(word::null),
145  UName_(word::null),
146  diskDir_(vector::uniform(NaN)),
147  Cp_(NaN),
148  Ct_(NaN),
149  diskArea_(NaN),
150  upstreamPoint_(vector::uniform(NaN)),
151  upstreamCellId_(-1)
152 {
153  readCoeffs();
154 }
155 
156 
157 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158 
160 {
161  return wordList(1, UName_);
162 }
163 
164 
166 (
167  const volVectorField& U,
168  fvMatrix<vector>& eqn
169 ) const
170 {
171  addActuationDiskAxialInertialResistance
172  (
173  eqn.source(),
174  set_.cells(),
175  mesh().V(),
178  U
179  );
180 }
181 
182 
184 (
185  const volScalarField& rho,
186  const volVectorField& U,
187  fvMatrix<vector>& eqn
188 ) const
189 {
190  addActuationDiskAxialInertialResistance
191  (
192  eqn.source(),
193  set_.cells(),
194  mesh().V(),
196  rho,
197  U
198  );
199 }
200 
201 
203 (
204  const volScalarField& alpha,
205  const volScalarField& rho,
206  const volVectorField& U,
207  fvMatrix<vector>& eqn
208 ) const
209 {
210  addActuationDiskAxialInertialResistance
211  (
212  eqn.source(),
213  set_.cells(),
214  mesh().V(),
215  alpha,
216  rho,
217  U
218  );
219 }
220 
221 
223 {
224  set_.movePoints();
225  return true;
226 }
227 
228 
230 {
231  set_.topoChange(map);
232 }
233 
234 
236 {
237  set_.mapMesh(map);
238 }
239 
240 
242 {
243  set_.distribute(map);
244 }
245 
246 
248 {
249  if (fvModel::read(dict))
250  {
251  set_.read(coeffs());
252  readCoeffs();
253  return true;
254  }
255  else
256  {
257  return false;
258  }
259 }
260 
261 
262 // ************************************************************************* //
#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.
static word groupName(Name name, const word &group)
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
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
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:199
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:53
Actuation disk source.
virtual bool movePoints()
Update for mesh motion.
virtual void addSup(const volVectorField &U, fvMatrix< vector > &eqn) const
Source term to momentum equation.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
scalar Cp_
Power coefficient.
word UName_
Name of the velocity field.
word phaseName_
The name of the phase to which this fvModel applies.
vector diskDir_
Disk area normal.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read dictionary.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
point upstreamPoint_
Upstream point sample.
actuationDisk(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
label upstreamCellId_
Upstream cell ID.
scalar diskArea_
Disk area.
scalar Ct_
Thrust coefficient.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1727
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
const cellShapeList & cells
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
addBackwardCompatibleToRunTimeSelectionTable(fvConstraint, fixedTemperature, dictionary, fixedTemperatureConstraint, "fixedTemperatureConstraint")
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< word > wordList
A List of words.
Definition: fileName.H:54
List< label > labelList
A List of labels.
Definition: labelList.H:56
dimensionedSymmTensor sqr(const dimensionedVector &dv)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
vector point
Point is a vector.
Definition: point.H:41
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)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
labelList fv(nPoints)
dictionary dict