fixedTrim.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) 2012-2021 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 "fixedTrim.H"
28 #include "unitConversion.H"
29 #include "mathematicalConstants.H"
30 
31 using namespace Foam::constant;
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(fixedTrim, 0);
38 
39  addToRunTimeSelectionTable(trimModel, fixedTrim, dictionary);
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const fv::rotorDiskSource& rotor,
48  const dictionary& dict
49 )
50 :
51  trimModel(rotor, dict, typeName),
52  thetag_(rotor.set().cells().size(), 0.0)
53 {
54  read(dict);
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
59 
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
67 {
68  trimModel::read(dict);
69 
70  scalar theta0 = degToRad(coeffs_.lookup<scalar>("theta0"));
71  scalar theta1c = degToRad(coeffs_.lookup<scalar>("theta1c"));
72  scalar theta1s = degToRad(coeffs_.lookup<scalar>("theta1s"));
73 
74  const List<point>& x = rotor_.x();
75  forAll(thetag_, i)
76  {
77  scalar psi = x[i].y();
78  thetag_[i] = theta0 + theta1c*cos(psi) + theta1s*sin(psi);
79  }
80 }
81 
82 
84 {
85  return tmp<scalarField>(thetag_);
86 }
87 
88 
90 (
91  const vectorField& U,
92  vectorField& force
93 ) const
94 {}
95 
96 
98 (
99  const volScalarField rho,
100  const vectorField& U,
101  vectorField& force
102 ) const
103 {}
104 
105 
106 // ************************************************************************* //
Collection of constants.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual void correct(const vectorField &U, vectorField &force) const
Correct the model.
Definition: fixedTrim.C:90
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
virtual void read(const dictionary &dict)
Read.
Definition: trimModel.C:61
Unit conversion functions.
Trim model base class.
Definition: trimModel.H:50
const fvCellSet & set() const
Macros for easy insertion into run-time selection tables.
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
Cell based momentum source which approximates the mean effects of rotor forces on a cylindrical regio...
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedScalar cos(const dimensionedScalar &ds)
const volScalarField & psi
fixedTrim(const fv::rotorDiskSource &rotor, const dictionary &dict)
Constructor.
Definition: fixedTrim.C:46
dimensionedScalar sin(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
void read(const dictionary &dict)
Read.
Definition: fixedTrim.C:66
virtual tmp< scalarField > thetag() const
Return the geometric angle of attack [rad].
Definition: fixedTrim.C:83
defineTypeNameAndDebug(combustionModel, 0)
A class for managing temporary objects.
Definition: PtrList.H:53
virtual ~fixedTrim()
Destructor.
Definition: fixedTrim.C:60
Namespace for OpenFOAM.