Loading...
Searching...
No Matches
rotorDiskSource.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) 2018-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
29#include "rotorDiskSource.H"
31#include "trimModel.H"
32#include "fvMatrices.H"
33#include "geometricOneField.H"
34#include "syncTools.H"
35#include "unitConversion.H"
36
37using namespace Foam::constant;
39// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
40
41namespace Foam
42{
43 namespace fv
44 {
47 }
48}
49
50
51const Foam::Enum
52<
54>
56({
57 { geometryModeType::gmAuto, "auto" },
58 { geometryModeType::gmSpecified, "specified" },
59});
60
61
62const Foam::Enum
63<
65>
67({
68 { inletFlowType::ifFixed, "fixed" },
69 { inletFlowType::ifSurfaceNormal, "surfaceNormal" },
70 { inletFlowType::ifLocal, "local" },
71});
72
73
74// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
75
77{
78 // Set inflow type
79 switch (selectionMode())
80 {
81 case smCellSet:
82 case smCellZone:
83 case smAll:
84 {
85 // Set the profile ID for each blade section
87 switch (inletFlow_)
88 {
89 case ifFixed:
90 {
91 coeffs_.readEntry("inletVelocity", inletVelocity_);
92 break;
93 }
94 case ifSurfaceNormal:
95 {
96 scalar UIn(coeffs_.get<scalar>("inletNormalVelocity"));
97 inletVelocity_ = -coordSys_.e3()*UIn;
98 break;
99 }
100 case ifLocal:
101 {
102 break;
103 }
104 default:
105 {
107 << "Unknown inlet velocity type" << abort(FatalError);
108 }
109 }
110
111 break;
112 }
113 default:
114 {
116 << "Source cannot be used with '"
118 << "' mode. Please use one of: " << nl
122 << exit(FatalError);
123 }
124 }
125}
126
127
129{
130 area_ = 0.0;
131
132 static const scalar tol = 0.8;
133
134 const label nInternalFaces = mesh_.nInternalFaces();
135 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
136 const vectorField& Sf = mesh_.Sf();
137 const scalarField& magSf = mesh_.magSf();
138
139 vector n = Zero;
140
141 // Calculate cell addressing for selected cells
142 labelList cellAddr(mesh_.nCells(), -1);
143 labelUIndList(cellAddr, cells_) = identity(cells_.size());
144 labelList nbrFaceCellAddr(mesh_.nBoundaryFaces(), -1);
145 forAll(pbm, patchi)
146 {
147 const polyPatch& pp = pbm[patchi];
148
149 if (pp.coupled())
150 {
151 forAll(pp, i)
152 {
153 label facei = pp.start() + i;
154 label nbrFacei = facei - nInternalFaces;
155 label own = mesh_.faceOwner()[facei];
156 nbrFaceCellAddr[nbrFacei] = cellAddr[own];
157 }
158 }
159 }
160
161 // Correct for parallel running
162 syncTools::swapBoundaryFaceList(mesh_, nbrFaceCellAddr);
163
164 // Add internal field contributions
165 for (label facei = 0; facei < nInternalFaces; facei++)
166 {
167 const label own = cellAddr[mesh_.faceOwner()[facei]];
168 const label nbr = cellAddr[mesh_.faceNeighbour()[facei]];
169
170 if ((own != -1) && (nbr == -1))
171 {
172 vector nf = Sf[facei]/magSf[facei];
173
174 if ((nf & axis) > tol)
175 {
176 area_[own] += magSf[facei];
177 n += Sf[facei];
178 }
179 }
180 else if ((own == -1) && (nbr != -1))
181 {
182 vector nf = Sf[facei]/magSf[facei];
183
184 if ((-nf & axis) > tol)
185 {
186 area_[nbr] += magSf[facei];
187 n -= Sf[facei];
188 }
189 }
190 }
191
192
193 // Add boundary contributions
194 forAll(pbm, patchi)
195 {
196 const polyPatch& pp = pbm[patchi];
197 const vectorField& Sfp = mesh_.Sf().boundaryField()[patchi];
198 const scalarField& magSfp = mesh_.magSf().boundaryField()[patchi];
199
200 if (pp.coupled())
201 {
202 forAll(pp, j)
203 {
204 const label facei = pp.start() + j;
205 const label own = cellAddr[mesh_.faceOwner()[facei]];
206 const label nbr = nbrFaceCellAddr[facei - nInternalFaces];
207 const vector nf = Sfp[j]/magSfp[j];
208
209 if ((own != -1) && (nbr == -1) && ((nf & axis) > tol))
210 {
211 area_[own] += magSfp[j];
212 n += Sfp[j];
213 }
214 }
215 }
216 else
217 {
218 forAll(pp, j)
219 {
220 const label facei = pp.start() + j;
221 const label own = cellAddr[mesh_.faceOwner()[facei]];
222 const vector nf = Sfp[j]/magSfp[j];
223
224 if ((own != -1) && ((nf & axis) > tol))
225 {
226 area_[own] += magSfp[j];
227 n += Sfp[j];
228 }
229 }
230 }
231 }
232
233 if (correct)
234 {
235 reduce(n, sumOp<vector>());
236 axis = n/mag(n);
237 }
238
239 if (debug)
240 {
241 auto tarea = volScalarField::New
242 (
243 IOobject::scopedName(name_, "area"),
245 mesh_,
247 );
248 auto& area = tarea.ref();
249
250 UIndirectList<scalar>(area.primitiveFieldRef(), cells_) = area_;
251
252 Info<< type() << ": " << name_
253 << " writing field " << area.name() << endl;
254
255 area.write();
256 }
257}
258
259
261{
262 // Construct the local rotor coordinate system
263 vector origin(Zero);
264 vector axis(Zero);
265 vector refDir(Zero);
266
267 geometryModeType gm =
268 geometryModeTypeNames_.get("geometryMode", coeffs_);
269
270 switch (gm)
271 {
272 case gmAuto:
273 {
274 // Determine rotation origin (cell volume weighted)
275 scalar sumV = 0.0;
276 const scalarField& V = mesh_.V();
277 const vectorField& C = mesh_.C();
278 forAll(cells_, i)
279 {
280 const label celli = cells_[i];
281 sumV += V[celli];
282 origin += V[celli]*C[celli];
283 }
284 reduce(origin, sumOp<vector>());
285 reduce(sumV, sumOp<scalar>());
286 origin /= sumV;
287
288 // Determine first radial vector
289 vector dx1(Zero);
290 scalar magR = -GREAT;
291 forAll(cells_, i)
292 {
293 const label celli = cells_[i];
294 vector test = C[celli] - origin;
295 if (mag(test) > magR)
296 {
297 dx1 = test;
298 magR = mag(test);
299 }
300 }
302 magR = mag(dx1);
303
304 // Determine second radial vector and cross to determine axis
305 forAll(cells_, i)
306 {
307 const label celli = cells_[i];
308 vector dx2 = C[celli] - origin;
309 if (mag(dx2) > 0.5*magR)
310 {
311 axis = dx1 ^ dx2;
312 if (mag(axis) > SMALL)
313 {
314 break;
315 }
316 }
317 }
319 axis.normalise();
320
321 // Correct the axis direction using a point above the rotor
322 {
323 vector pointAbove(coeffs_.get<vector>("pointAbove"));
324 vector dir = pointAbove - origin;
325 dir.normalise();
326 if ((dir & axis) < 0)
327 {
328 axis *= -1.0;
329 }
330 }
331
332 coeffs_.readEntry("refDirection", refDir);
333
334 // Set the face areas and apply correction to calculated axis
335 // e.g. if cellZone is more than a single layer in thickness
336 setFaceArea(axis, true);
337
338 break;
339 }
340 case gmSpecified:
341 {
342 coeffs_.readEntry("origin", origin);
343 coeffs_.readEntry("axis", axis);
344 coeffs_.readEntry("refDirection", refDir);
345
346 setFaceArea(axis, false);
347
348 break;
349 }
350 default:
351 {
353 << "Unknown geometryMode " << geometryModeTypeNames_[gm]
354 << ". Available geometry modes include "
355 << geometryModeTypeNames_
356 << exit(FatalError);
357 }
358 }
359
360 coordSys_ = coordSystem::cylindrical(origin, axis, refDir);
361
362 const scalar sumArea = gSum(area_);
363 const scalar diameter = Foam::sqrt(4.0*sumArea/mathematical::pi);
364 Info<< " Rotor gometry:" << nl
365 << " - disk diameter = " << diameter << nl
366 << " - disk area = " << sumArea << nl
367 << " - origin = " << coordSys_.origin() << nl
368 << " - r-axis = " << coordSys_.e1() << nl
369 << " - psi-axis = " << coordSys_.e2() << nl
370 << " - z-axis = " << coordSys_.e3() << endl;
371}
372
373
375{
376 const pointUIndList cc(mesh_.C(), cells_);
377
378 // Optional: for later transform(), invTransform()
380
381 forAll(cells_, i)
382 {
383 if (area_[i] > ROOTVSMALL)
384 {
385 // Position in (planar) rotor coordinate system
386 x_[i] = coordSys_.localPosition(cc[i]);
387
388 // Cache max radius
389 rMax_ = max(rMax_, x_[i].x());
390
391 // Swept angle relative to rDir axis [radians] in range 0 -> 2*pi
392 scalar psi = x_[i].y();
393
394 // Blade flap angle [radians]
395 scalar beta =
396 flap_.beta0 - flap_.beta1c*cos(psi) - flap_.beta2s*sin(psi);
397
398 // Determine rotation tensor to convert from planar system into the
399 // rotor cone system
400 scalar c = cos(beta);
401 scalar s = sin(beta);
402 Rcone_[i] = tensor(c, 0, -s, 0, 1, 0, s, 0, c);
403 }
404 }
405}
406
407
409(
410 const volVectorField& U
411) const
412{
413 switch (inletFlow_)
414 {
415 case ifFixed:
416 case ifSurfaceNormal:
417 {
418 return tmp<vectorField>::New(mesh_.nCells(), inletVelocity_);
419
420 break;
421 }
422 case ifLocal:
423 {
424 return U.primitiveField();
425
426 break;
427 }
428 default:
429 {
431 << "Unknown inlet flow specification" << abort(FatalError);
432 }
433 }
435 return tmp<vectorField>::New(mesh_.nCells(), Zero);
436}
437
438
439// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
440
442(
443 const word& name,
444 const word& modelType,
445 const dictionary& dict,
446 const fvMesh& mesh
447
448)
449:
450 fv::cellSetOption(name, modelType, dict, mesh),
451 rhoRef_(1.0),
452 omega_(0.0),
453 nBlades_(0),
454 inletFlow_(ifLocal),
455 inletVelocity_(Zero),
456 tipEffect_(1.0),
457 flap_(),
458 x_(cells_.size(), Zero),
459 Rcone_(cells_.size(), I),
460 area_(cells_.size(), Zero),
461 coordSys_(),
462 rMax_(0.0),
463 trim_(trimModel::New(*this, coeffs_)),
464 blade_(coeffs_.subDict("blade")),
465 profiles_(coeffs_.subDict("profiles"))
467 read(dict);
468}
469
470
471// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
472
474(
475 fvMatrix<vector>& eqn,
476 const label fieldi
477)
478{
479 auto tforce = volVectorField::New
480 (
481 IOobject::scopedName(name_, "rotorForce"),
483 mesh_,
485 );
486 auto& force = tforce.ref();
487
488 // Read the reference density for incompressible flow
489 coeffs_.readEntry("rhoRef", rhoRef_);
490
491 const vectorField Uin(inflowVelocity(eqn.psi()));
492 trim_->correct(Uin, force);
493 calculate(geometricOneField(), Uin, trim_->thetag(), force);
494
495 // Add source to rhs of eqn
496 eqn -= force;
497
498 if (mesh_.time().writeTime())
499 {
500 force.write();
501 }
502}
503
504
506(
507 const volScalarField& rho,
508 fvMatrix<vector>& eqn,
509 const label fieldi
510)
511{
512 auto tforce = volVectorField::New
513 (
514 IOobject::scopedName(name_, "rotorForce"),
516 mesh_,
518 );
519 auto& force = tforce.ref();
520
521 const vectorField Uin(inflowVelocity(eqn.psi()));
522 trim_->correct(rho, Uin, force);
523 calculate(rho, Uin, trim_->thetag(), force);
524
525 // Add source to rhs of eqn
526 eqn -= force;
527
528 if (mesh_.time().writeTime())
529 {
530 force.write();
531 }
532}
533
534
536{
538 {
539 coeffs_.readEntry("fields", fieldNames_);
541
542 // Read coordinate system/geometry invariant properties
543 omega_ = rpmToRads(coeffs_.get<scalar>("rpm"));
544
545 coeffs_.readEntry("nBlades", nBlades_);
546
547 inletFlowTypeNames_.readEntry("inletFlowType", coeffs_, inletFlow_);
548
549 coeffs_.readEntry("tipEffect", tipEffect_);
550
551 const dictionary& flapCoeffs(coeffs_.subDict("flapCoeffs"));
552 flap_.beta0 = degToRad(flapCoeffs.get<scalar>("beta0"));
553 flap_.beta1c = degToRad(flapCoeffs.get<scalar>("beta1c"));
554 flap_.beta2s = degToRad(flapCoeffs.get<scalar>("beta2s"));
555
556
557 // Create coordinate system
558 createCoordinateSystem();
559
560 // Read coordinate system dependent properties
561 checkData();
562
563 constructGeometry();
564
565 trim_->read(coeffs_);
566
567 if (debug)
568 {
569 writeField("thetag", trim_->thetag()(), true);
570 writeField("faceArea", area_, true);
571 }
572
573 return true;
574 }
575
576 return false;
577}
578
579
580// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
Graphite solid properties.
Definition C.H:49
C() noexcept
Default construct.
Definition C.C:36
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
A List with indirect addressing. Like IndirectList but does not store addressing.
Vector< Cmpt > & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the vector by its magnitude.
Definition VectorI.H:114
const List< label > & profileID() const
Return const access to the profile ID list.
Definition bladeModel.C:139
const List< word > & profileName() const
Return const access to the profile name list.
Definition bladeModel.C:133
A cylindrical coordinate system (r-theta-z). The coordinate system angle theta is always in radians.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
const dimensionSet & dimensions() const noexcept
Definition fvMatrix.H:530
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Intermediate abstract class for handling cell-set options for the derived fvOptions.
scalar V() const noexcept
Return const access to the total cell volume.
labelList cells_
Set of cells to apply source to.
virtual bool read(const dictionary &dict)
Read source dictionary.
cellSetOption(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
static const Enum< selectionModeType > selectionModeTypeNames_
List of selection mode type names.
selectionModeType selectionMode() const noexcept
Return the cell selection mode.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition fvOption.H:124
const word & name() const noexcept
Return const access to the source name.
Definition fvOptionI.H:24
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
static autoPtr< option > New(const word &name, const dictionary &dict, const fvMesh &mesh)
Return a reference to the selected fvOption model.
Definition fvOption.C:78
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition fvOption.H:157
dictionary coeffs_
Dictionary containing source coefficients.
Definition fvOption.H:152
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition fvOption.C:41
virtual void correct(volScalarField &field)
Definition fvOption.C:314
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition fvOptionI.H:30
const word name_
Source name.
Definition fvOption.H:132
Applies cell-based momentum sources on velocity (i.e. U) within a specified cylindrical region to app...
scalar omega_
Rotational speed [rad/s].
rotorDiskSource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
static const Enum< geometryModeType > geometryModeTypeNames_
Names for geometryModeType.
void calculate(const RhoFieldType &rho, const vectorField &U, const scalarField &thetag, vectorField &force, const bool divideVolume=true, const bool output=true) const
Calculate forces.
scalar rhoRef_
Reference density for incompressible case.
void setFaceArea(vector &axis, const bool correct)
Set the face areas per cell, and optionally correct the rotor axis.
void checkData()
Check data.
List< scalar > area_
Area [m2].
scalar tipEffect_
Tip effect [0-1].
inletFlowType inletFlow_
Inlet flow type.
void constructGeometry()
Construct geometry.
void writeField(const word &name, const List< Type > &values, const bool writeNow=false) const
Helper function to write rotor values.
coordSystem::cylindrical coordSys_
Rotor local cylindrical coordinate system (r-theta-z).
scalar rMax_
Maximum radius.
virtual bool read(const dictionary &dict)
Read source dictionary.
List< point > x_
Cell centre positions in local rotor frame.
bladeModel blade_
Blade data.
tmp< vectorField > inflowVelocity(const volVectorField &U) const
Return the inlet flow field.
void createCoordinateSystem()
Create the coordinate system.
flapData flap_
Blade flap coefficients [rad/s].
vector inletVelocity_
Inlet velocity for specified inflow.
inletFlowType
Options for the inlet flow type specification.
List< tensor > Rcone_
Rotation tensor for flap angle.
label nBlades_
Number of blades.
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Add explicit contribution to momentum equation.
autoPtr< trimModel > trim_
Trim model.
profileModelList profiles_
Profile data.
static const Enum< inletFlowType > inletFlowTypeNames_
Names for inletFlowType.
geometryModeType
Options for the geometry type specification.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
void connectBlades(const List< word > &names, List< label > &addr) const
Set blade->profile addressing.
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const bool parRun=UPstream::parRun())
Swap coupled boundary face values. Uses eqOp.
Definition syncTools.H:524
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
Base class for trim models for handling blade characteristics and thrust-torque relations.
Definition trimModel.H:108
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
thermo correct()
const volScalarField & psi
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
A special matrix type and solver, designed for finite volume solutions of scalar equations.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
constexpr scalar pi(M_PI)
Different types of constants.
Namespace for handling debugging switches.
Definition debug.C:45
const wordList area
Standard area field types (scalar, vector, tensor, etc).
Namespace for finite-volume.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Type gSum(const FieldField< Field, Type > &f)
GeometricField< vector, fvPatchField, volMesh > volVectorField
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
UIndirectList< label > labelUIndList
UIndirectList of labels.
const dimensionSet dimArea(sqr(dimLength))
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
static const Identity< scalar > I
Definition Identity.H:100
UIndirectList< point > pointUIndList
UIndirectList of point.
Tensor< scalar > tensor
Definition symmTensor.H:57
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
label checkData(const objectRegistry &obr, const instantList &timeDirs, wordList &objectNames, const fileName &local=fileName::null)
Check if fields are good to use (available at all times).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr scalar rpmToRads() noexcept
Multiplication factor for revolutions/minute to radians/sec.
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensionedScalar cos(const dimensionedScalar &ds)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
fvOptions correct(rho)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.