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-2022 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 {
38  defineTypeNameAndDebug(actuationDiskSource, 0);
40  (
41  fvModel,
42  actuationDiskSource,
43  dictionary
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",
59  IOobject::groupName("U", phaseName_)
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");
88  upstreamCellId_ = mesh().findCell(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 dictionary& dict,
105  const fvMesh& mesh
106 )
107 :
108  fvModel(name, modelType, dict, mesh),
109  set_(coeffs(), mesh),
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  if (set_.V() > vSmall)
142  {
143  addActuationDiskAxialInertialResistance
144  (
145  Usource,
146  set_.cells(),
147  cellsV,
150  U
151  );
152  }
153 }
154 
155 
157 (
158  const volScalarField& rho,
159  fvMatrix<vector>& eqn,
160  const word& fieldName
161 ) const
162 {
163  const scalarField& cellsV = mesh().V();
164  vectorField& Usource = eqn.source();
165  const vectorField& U = eqn.psi();
166 
167  if (set_.V() > vSmall)
168  {
169  addActuationDiskAxialInertialResistance
170  (
171  Usource,
172  set_.cells(),
173  cellsV,
175  rho,
176  U
177  );
178  }
179 }
180 
181 
183 (
184  const volScalarField& alpha,
185  const volScalarField& rho,
186  fvMatrix<vector>& eqn,
187  const word& fieldName
188 ) const
189 {
190  const scalarField& cellsV = mesh().V();
191  vectorField& Usource = eqn.source();
192  const vectorField& U = eqn.psi();
193 
194  if (set_.V() > vSmall)
195  {
196  addActuationDiskAxialInertialResistance
197  (
198  Usource,
199  set_.cells(),
200  cellsV,
201  alpha,
202  rho,
203  U
204  );
205  }
206 }
207 
208 
210 {
211  set_.movePoints();
212  return true;
213 }
214 
215 
217 {
218  set_.topoChange(map);
219 }
220 
221 
223 {
224  set_.mapMesh(map);
225 }
226 
227 
229 (
230  const polyDistributionMap& map
231 )
232 {
233  set_.distribute(map);
234 }
235 
236 
238 {
239  if (fvModel::read(dict))
240  {
241  set_.read(coeffs());
242  readCoeffs();
243  return true;
244  }
245  else
246  {
247  return false;
248  }
249 }
250 
251 
252 // ************************************************************************* //
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
U
Definition: pEqn.H:72
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:290
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:164
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
actuationDiskSource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Finite volume model abstract base class.
Definition: fvModel.H:57
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual void addSup(fvMatrix< vector > &eqn, const word &fieldName) const
Source term to momentum equation.
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:59
labelList fv(nPoints)
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
static word groupName(Name name, const word &group)
static const word null
An empty word.
Definition: word.H:77
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Field< Type > & source()
Definition: fvMatrix.H:300
virtual bool read(const dictionary &dict)
Read dictionary.
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
List< word > wordList
A List of words.
Definition: fileName.H:54
vector point
Point is a vector.
Definition: point.H:41
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
static Vector< scalar > uniform(const scalar &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
dimensioned< scalar > mag(const dimensioned< Type > &)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual bool movePoints()
Update for mesh motion.
Namespace for OpenFOAM.
IOerror FatalIOError