1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2006-7 H. Jasak All rights reserved
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
25 \*---------------------------------------------------------------------------*/
27 #include "mixerGgiFvMesh.H"
29 #include "regionSplit.H"
30 #include "ggiPolyPatch.H"
31 #include "polyPatchID.H"
32 #include "addToRunTimeSelectionTable.H"
33 #include "mapPolyMesh.H"
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 defineTypeNameAndDebug(mixerGgiFvMesh, 0);
40 addToRunTimeSelectionTable(dynamicFvMesh, mixerGgiFvMesh, IOobject);
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 void Foam::mixerGgiFvMesh::addZonesAndModifiers()
48 // Add zones and modifiers for motion action
50 if (cellZones().size() > 0)
52 Info<< "void mixerGgiFvMesh::addZonesAndModifiers() : "
53 << "Zones and modifiers already present. Skipping."
59 Info<< "Time = " << time().timeName() << endl
60 << "Adding zones and modifiers to the mesh" << endl;
63 List<pointZone*> pz(0);
64 List<faceZone*> fz(0);
65 List<cellZone*> cz(1);
67 // Copy the face zones associated with the GGI interfaces
68 if (faceZones().size() > 0)
71 Info << "Copying existing face zones" << endl;
73 fz.setSize(faceZones().size());
75 forAll (faceZones(), i)
77 fz[i] = faceZones()[i].clone(faceZones()).ptr();
81 regionSplit rs(*this);
83 // Get the region of the cell containing the origin.
84 label originRegion = rs[findNearestCell(cs().origin())];
86 labelList movingCells(nCells());
87 label nMovingCells = 0;
91 if (rs[cellI] == originRegion)
93 movingCells[nMovingCells] = cellI;
98 movingCells.setSize(nMovingCells);
99 Info << "Number of cells in the moving region: " << nMovingCells << endl;
109 Info << "Adding point, face and cell zones" << endl;
111 addZones(pz, fz, cz);
119 void Foam::mixerGgiFvMesh::calcMovingMasks() const
123 Info<< "void mixerGgiFvMesh::calcMovingMasks() const : "
124 << "Calculating point and cell masks"
128 if (movingPointsMaskPtr_)
130 FatalErrorIn("void mixerGgiFvMesh::calcMovingMasks() const")
131 << "point mask already calculated"
132 << abort(FatalError);
135 // Set the point mask
136 movingPointsMaskPtr_ = new scalarField(allPoints().size(), 0);
137 scalarField& movingPointsMask = *movingPointsMaskPtr_;
139 const cellList& c = cells();
140 const faceList& f = allFaces();
142 label movingCellsID = cellZones().findZoneID("movingCells");
144 if (movingCellsID < 0)
146 FatalErrorIn("void mixerGgiFvMesh::calcMovingMasks() const")
147 << "Cannot find moving cell zone ID"
148 << abort(FatalError);
151 const labelList& cellAddr = cellZones()[movingCellsID];
153 forAll (cellAddr, cellI)
155 const cell& curCell = c[cellAddr[cellI]];
157 forAll (curCell, faceI)
159 // Mark all the points as moving
160 const face& curFace = f[curCell[faceI]];
162 forAll (curFace, pointI)
164 movingPointsMask[curFace[pointI]] = 1;
169 // Grab the ggi patches on the moving side
170 wordList movingPatches(dict_.subDict("slider").lookup("moving"));
172 forAll (movingPatches, patchI)
174 const label movingSliderID =
175 boundaryMesh().findPatchID(movingPatches[patchI]);
177 if (movingSliderID < 0)
179 FatalErrorIn("void mixerGgiFvMesh::calcMovingMasks() const")
180 << "Moving slider named " << movingPatches[patchI]
181 << " not found. Valid patch names: "
182 << boundaryMesh().names() << abort(FatalError);
185 const ggiPolyPatch& movingGgiPatch =
186 refCast<const ggiPolyPatch>(boundaryMesh()[movingSliderID]);
188 const labelList& movingSliderAddr = movingGgiPatch.zone();
190 forAll (movingSliderAddr, faceI)
192 const face& curFace = f[movingSliderAddr[faceI]];
194 forAll (curFace, pointI)
196 movingPointsMask[curFace[pointI]] = 1;
201 // Grab the ggi patches on the static side
202 wordList staticPatches(dict_.subDict("slider").lookup("static"));
204 forAll (staticPatches, patchI)
206 const label staticSliderID =
207 boundaryMesh().findPatchID(staticPatches[patchI]);
209 if (staticSliderID < 0)
211 FatalErrorIn("void mixerGgiFvMesh::calcMovingMasks() const")
212 << "Static slider named " << staticPatches[patchI]
213 << " not found. Valid patch names: "
214 << boundaryMesh().names() << abort(FatalError);
217 const ggiPolyPatch& staticGgiPatch =
218 refCast<const ggiPolyPatch>(boundaryMesh()[staticSliderID]);
220 const labelList& staticSliderAddr = staticGgiPatch.zone();
222 forAll (staticSliderAddr, faceI)
224 const face& curFace = f[staticSliderAddr[faceI]];
226 forAll (curFace, pointI)
228 movingPointsMask[curFace[pointI]] = 0;
235 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
237 // Construct from components
238 Foam::mixerGgiFvMesh::mixerGgiFvMesh
256 ).subDict(typeName + "Coeffs")
261 dict_.subDict("coordinateSystem")
263 rpm_(readScalar(dict_.lookup("rpm"))),
264 movingPointsMaskPtr_(NULL)
266 // Make sure the coordinate system does not operate in degrees
267 // Bug fix, HJ, 3/Oct/2011
268 if (!cs_.inDegrees())
270 WarningIn("mixerGgiFvMesh::mixerGgiFvMesh(const IOobject& io)")
271 << "Mixer coordinate system is set to operate in radians. "
272 << "Changing to rad for correct calculation of angular velocity."
274 << "To remove this message please add entry" << nl << nl
275 << "inDegrees true;" << nl << nl
276 << "to the specification of the coordinate system"
279 cs_.inDegrees() = true;
282 addZonesAndModifiers();
284 Info<< "Mixer mesh:" << nl
285 << " origin: " << cs().origin() << nl
286 << " axis : " << cs().axis() << nl
287 << " rpm : " << rpm_ << endl;
291 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
293 Foam::mixerGgiFvMesh::~mixerGgiFvMesh()
295 deleteDemandDrivenData(movingPointsMaskPtr_);
299 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
301 // Return moving points mask. Moving points marked with 1
302 const Foam::scalarField& Foam::mixerGgiFvMesh::movingPointsMask() const
304 if (!movingPointsMaskPtr_)
309 return *movingPointsMaskPtr_;
313 bool Foam::mixerGgiFvMesh::update()
315 // Rotational speed needs to be converted from rpm
320 cs_.localPosition(allPoints())
321 + vector(0, rpm_*360.0*time().deltaT().value()/60.0, 0)
326 // The mesh is not morphing
331 // ************************************************************************* //