1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
26 Hrvoje Jasak, Wikki Ltd. All rights reserved.
27 Fethi Tekin, All rights reserved.
28 Oliver Borm, All rights reserved.
30 \*---------------------------------------------------------------------------*/
32 #include "turboFvMesh.H"
34 #include "addToRunTimeSelectionTable.H"
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 defineTypeNameAndDebug(turboFvMesh, 0);
42 addToRunTimeSelectionTable(dynamicFvMesh, turboFvMesh, IOobject);
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 void Foam::turboFvMesh::calcMovingPoints() const
52 Info<< "void turboFvMesh::calcMovingMasks() const : "
53 << "Calculating point and cell masks"
59 FatalErrorIn("void turboFvMesh::calcMovingMasks() const")
60 << "point mask already calculated"
64 // Retrieve the cell zone Names
65 const wordList cellZoneNames = cellZones().names();
68 movingPointsPtr_ = new vectorField(allPoints().size(), vector::zero);
70 vectorField& movingPoints = *movingPointsPtr_;
72 const cellList& c = cells();
73 const faceList& f = allFaces();
77 forAll (cellZoneNames,cellZoneI)
79 const labelList& cellAddr =
80 cellZones()[cellZones().findZoneID(cellZoneNames[cellZoneI])];
82 if (dict_.subDict("rpm").found(cellZoneNames[cellZoneI]))
86 dict_.subDict("rpm").lookup(cellZoneNames[cellZoneI])
89 Info<< "Moving Cell Zone Name: " << cellZoneNames[cellZoneI]
90 << " rpm: " << rpm << endl;
92 forAll (cellAddr, cellI)
94 const cell& curCell = c[cellAddr[cellI]];
96 forAll (curCell, faceI)
98 const face& curFace = f[curCell[faceI]];
100 forAll (curFace, pointI)
102 // The rotation data is saved within the cell data. For
103 // non-rotating regions rpm is zero, so mesh movement
104 // is also zero. The conversion of rotational speed
106 // Note: deltaT changes during the run: moved to
107 // turboFvMesh::update(). HJ, 14/Oct/2010
108 movingPoints[curFace[pointI]] =
109 vector(0, rpm/60.0*360.0, 0);
116 // Retrieve the face zone Names
117 const wordList faceZoneNames = faceZones().names();
119 // This is not bullet proof, as one could specify the wrong name of a
120 // faceZone, which is then not selected. The solver will crash after the
121 // first meshMotion iteration.
122 forAll (faceZoneNames, faceZoneI)
124 if (dict_.subDict("slider").found(faceZoneNames[faceZoneI]))
128 dict_.subDict("slider").lookup(faceZoneNames[faceZoneI])
131 Info<< "Moving Face Zone Name: " << faceZoneNames[faceZoneI]
132 << " rpm: " << rpm << endl;
134 faceZoneID zoneID(faceZoneNames[faceZoneI], faceZones());
136 const labelList& movingSliderAddr = faceZones()[zoneID.index()];
138 forAll (movingSliderAddr, faceI)
140 const face& curFace = f[movingSliderAddr[faceI]];
142 forAll (curFace, pointI)
144 movingPoints[curFace[pointI]] =
145 vector( 0, rpm/60.0*360.0, 0);
153 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
155 // Construct from components
156 Foam::turboFvMesh::turboFvMesh
174 ).subDict(typeName + "Coeffs")
179 dict_.subDict("coordinateSystem")
181 movingPointsPtr_(NULL)
183 // Make sure the coordinate system does not operate in degrees
184 // Bug fix, HJ, 3/Oct/2011
185 if (!cs_.inDegrees())
187 WarningIn("turboFvMesh::turboFvMesh(const IOobject& io)")
188 << "Mixer coordinate system is set to operate in radians. "
189 << "Changing to rad for correct calculation of angular velocity."
191 << "To remove this message please add entry" << nl << nl
192 << "inDegrees true;" << nl << nl
193 << "to the specification of the coordinate system"
196 cs_.inDegrees() = true;
200 Info<< "Turbomachine Mixer mesh:" << nl
201 << " origin: " << cs().origin() << nl
202 << " axis : " << cs().axis() << endl;
206 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
208 Foam::turboFvMesh::~turboFvMesh()
210 deleteDemandDrivenData(movingPointsPtr_);
214 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
216 // Return moving points mask
217 const Foam::vectorField& Foam::turboFvMesh::movingPoints() const
219 if (!movingPointsPtr_)
224 return *movingPointsPtr_;
228 bool Foam::turboFvMesh::update()
234 cs_.localPosition(allPoints())
235 + movingPoints()*time().deltaT().value()
239 // The mesh is not morphing
244 // ************************************************************************* //