1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend 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 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
25 Hrvoje Jasak, Wikki Ltd. All rights reserved.
26 Fethi Tekin, All rights reserved.
27 Oliver Borm, All rights reserved.
29 \*---------------------------------------------------------------------------*/
31 #include "turboFvMesh.H"
33 #include "addToRunTimeSelectionTable.H"
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 defineTypeNameAndDebug(turboFvMesh, 0);
41 addToRunTimeSelectionTable(dynamicFvMesh, turboFvMesh, IOobject);
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 void Foam::turboFvMesh::calcMovingPoints() const
51 Info<< "void turboFvMesh::calcMovingMasks() const : "
52 << "Calculating point and cell masks"
58 FatalErrorIn("void turboFvMesh::calcMovingMasks() const")
59 << "point mask already calculated"
63 // Retrieve the cell zone Names
64 const wordList cellZoneNames = cellZones().names();
67 movingPointsPtr_ = new vectorField(allPoints().size(), vector::zero);
69 vectorField& movingPoints = *movingPointsPtr_;
71 const cellList& c = cells();
72 const faceList& f = allFaces();
76 forAll (cellZoneNames,cellZoneI)
78 const labelList& cellAddr =
79 cellZones()[cellZones().findZoneID(cellZoneNames[cellZoneI])];
81 if (dict_.subDict("rpm").found(cellZoneNames[cellZoneI]))
85 dict_.subDict("rpm").lookup(cellZoneNames[cellZoneI])
88 Info<< "Moving Cell Zone Name: " << cellZoneNames[cellZoneI]
89 << " rpm: " << rpm << endl;
91 forAll (cellAddr, cellI)
93 const cell& curCell = c[cellAddr[cellI]];
95 forAll (curCell, faceI)
97 const face& curFace = f[curCell[faceI]];
99 forAll (curFace, pointI)
101 // The rotation data is saved within the cell data. For
102 // non-rotating regions rpm is zero, so mesh movement
103 // is also zero. The conversion of rotational speed
105 // Note: deltaT changes during the run: moved to
106 // turboFvMesh::update(). HJ, 14/Oct/2010
107 movingPoints[curFace[pointI]] =
108 vector(0, rpm/60.0*360.0, 0);
115 // Retrieve the face zone Names
116 const wordList faceZoneNames = faceZones().names();
118 // This is not bullet proof, as one could specify the wrong name of a
119 // faceZone, which is then not selected. The solver will crash after the
120 // first meshMotion iteration.
121 forAll (faceZoneNames, faceZoneI)
123 if (dict_.subDict("slider").found(faceZoneNames[faceZoneI]))
127 dict_.subDict("slider").lookup(faceZoneNames[faceZoneI])
130 Info<< "Moving Face Zone Name: " << faceZoneNames[faceZoneI]
131 << " rpm: " << rpm << endl;
133 faceZoneID zoneID(faceZoneNames[faceZoneI], faceZones());
135 const labelList& movingSliderAddr = faceZones()[zoneID.index()];
137 forAll (movingSliderAddr, faceI)
139 const face& curFace = f[movingSliderAddr[faceI]];
141 forAll (curFace, pointI)
143 movingPoints[curFace[pointI]] =
144 vector( 0, rpm/60.0*360.0, 0);
152 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
154 // Construct from components
155 Foam::turboFvMesh::turboFvMesh
173 ).subDict(typeName + "Coeffs")
178 dict_.subDict("coordinateSystem")
180 movingPointsPtr_(NULL)
182 // Make sure the coordinate system does not operate in degrees
183 // Bug fix, HJ, 3/Oct/2011
184 if (!cs_.inDegrees())
186 WarningIn("turboFvMesh::turboFvMesh(const IOobject& io)")
187 << "Mixer coordinate system is set to operate in radians. "
188 << "Changing to rad for correct calculation of angular velocity."
190 << "To remove this message please add entry" << nl << nl
191 << "inDegrees true;" << nl << nl
192 << "to the specification of the coordinate system"
195 cs_.inDegrees() = true;
199 Info<< "Turbomachine Mixer mesh:" << nl
200 << " origin: " << cs().origin() << nl
201 << " axis : " << cs().axis() << endl;
205 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
207 Foam::turboFvMesh::~turboFvMesh()
209 deleteDemandDrivenData(movingPointsPtr_);
213 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
215 // Return moving points mask
216 const Foam::vectorField& Foam::turboFvMesh::movingPoints() const
218 if (!movingPointsPtr_)
223 return *movingPointsPtr_;
227 bool Foam::turboFvMesh::update()
233 cs_.localPosition(allPoints())
234 + movingPoints()*time().deltaT().value()
238 // The mesh is not morphing, but flux re-calculation is required
244 // ************************************************************************* //