actuationDiskSource.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-2023 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 "actuationDiskSource.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 {
40  (
41  fvModel,
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::fv::actuationDiskSource::readCoeffs()
52 {
53  phaseName_ = coeffs().lookupOrDefault<word>("phase", word::null);
54 
55  UName_ =
56  coeffs().lookupOrDefault<word>
57  (
58  "U",
60  );
61 
62  diskDir_ = coeffs().lookup<vector>("diskDir");
63  if (mag(diskDir_) < vSmall)
64  {
66  << "disk direction vector is approximately zero"
67  << exit(FatalIOError);
68  }
69 
70  Cp_ = coeffs().lookup<scalar>("Cp");
71  Ct_ = coeffs().lookup<scalar>("Ct");
72  if (Cp_ <= vSmall || Ct_ <= vSmall)
73  {
75  << "Cp and Ct must be greater than zero"
76  << exit(FatalIOError);
77  }
78 
79  diskArea_ = coeffs().lookup<scalar>("diskArea");
80  if (magSqr(diskArea_) <= vSmall)
81  {
83  << "diskArea is approximately zero"
84  << exit(FatalIOError);
85  }
86 
87  upstreamPoint_ = coeffs().lookup<point>("upstreamPoint");
89  if (returnReduce(upstreamCellId_, maxOp<label>()) == -1)
90  {
92  << "upstream location " << upstreamPoint_ << " not found in mesh"
93  << exit(FatalIOError);
94  }
95 }
96 
97 
98 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
99 
101 (
102  const word& name,
103  const word& modelType,
104  const fvMesh& mesh,
105  const dictionary& dict
106 )
107 :
108  fvModel(name, modelType, mesh, dict),
109  set_(mesh, coeffs()),
110  phaseName_(word::null),
111  UName_(word::null),
112  diskDir_(vector::uniform(NaN)),
113  Cp_(NaN),
114  Ct_(NaN),
115  diskArea_(NaN),
116  upstreamPoint_(vector::uniform(NaN)),
117  upstreamCellId_(-1)
118 {
119  readCoeffs();
120 }
121 
122 
123 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
124 
126 {
127  return wordList(1, UName_);
128 }
129 
130 
132 (
133  fvMatrix<vector>& eqn,
134  const word& fieldName
135 ) const
136 {
137  const scalarField& cellsV = mesh().V();
138  vectorField& Usource = eqn.source();
139  const vectorField& U = eqn.psi();
140 
141  addActuationDiskAxialInertialResistance
142  (
143  Usource,
144  set_.cells(),
145  cellsV,
148  U
149  );
150 }
151 
152 
154 (
155  const volScalarField& rho,
156  fvMatrix<vector>& eqn,
157  const word& fieldName
158 ) const
159 {
160  const scalarField& cellsV = mesh().V();
161  vectorField& Usource = eqn.source();
162  const vectorField& U = eqn.psi();
163 
164  addActuationDiskAxialInertialResistance
165  (
166  Usource,
167  set_.cells(),
168  cellsV,
170  rho,
171  U
172  );
173 }
174 
175 
177 (
178  const volScalarField& alpha,
179  const volScalarField& rho,
180  fvMatrix<vector>& eqn,
181  const word& fieldName
182 ) const
183 {
184  const scalarField& cellsV = mesh().V();
185  vectorField& Usource = eqn.source();
186  const vectorField& U = eqn.psi();
187 
188  addActuationDiskAxialInertialResistance
189  (
190  Usource,
191  set_.cells(),
192  cellsV,
193  alpha,
194  rho,
195  U
196  );
197 }
198 
199 
201 {
202  set_.movePoints();
203  return true;
204 }
205 
206 
208 {
209  set_.topoChange(map);
210 }
211 
212 
214 {
215  set_.mapMesh(map);
216 }
217 
218 
220 (
221  const polyDistributionMap& map
222 )
223 {
224  set_.distribute(map);
225 }
226 
227 
229 {
230  if (fvModel::read(dict))
231  {
232  set_.read(coeffs());
233  readCoeffs();
234  return true;
235  }
236  else
237  {
238  return false;
239  }
240 }
241 
242 
243 // ************************************************************************* //
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:160
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
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
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:40
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:187
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
virtual bool movePoints()
Update for mesh motion.
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 addSup(fvMatrix< vector > &eqn, const word &fieldName) const
Source term to momentum equation.
actuationDiskSource(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
point upstreamPoint_
Upstream point sample.
label upstreamCellId_
Upstream cell ID.
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:1775
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:318
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)
static Type NaN()
Return a primitive with all components set to NaN.
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
vector point
Point is a vector.
Definition: point.H:41
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
IOerror FatalIOError
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensioned< scalar > magSqr(const dimensioned< Type > &)
labelList fv(nPoints)
dictionary dict