1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd.
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
13 the Free Software Foundation, either version 3 of the License, or
14 (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "mixerFvMesh.H"
28 #include "regionSplit.H"
29 #include "slidingInterface.H"
30 #include "addToRunTimeSelectionTable.H"
31 #include "mapPolyMesh.H"
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 defineTypeNameAndDebug(mixerFvMesh, 0);
39 addToRunTimeSelectionTable(topoChangerFvMesh, mixerFvMesh, IOobject);
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 void Foam::mixerFvMesh::addZonesAndModifiers()
47 // Add zones and modifiers for motion action
54 || topoChanger_.size()
57 Info<< "void mixerFvMesh::addZonesAndModifiers() : "
58 << "Zones and modifiers already present. Skipping."
64 Info<< "Time = " << time().timeName() << endl
65 << "Adding zones and modifiers to the mesh" << endl;
68 List<pointZone*> pz(1);
70 // Add an empty zone for cut points
81 // Do face zones for slider
83 List<faceZone*> fz(3);
86 const word innerSliderName(motionDict_.subDict("slider").lookup("inside"));
87 const polyPatch& innerSlider =
88 boundaryMesh()[boundaryMesh().findPatchID(innerSliderName)];
90 labelList isf(innerSlider.size());
94 isf[i] = innerSlider.start() + i;
101 boolList(innerSlider.size(), false),
107 const word outerSliderName(motionDict_.subDict("slider").lookup("outside"));
108 const polyPatch& outerSlider =
109 boundaryMesh()[boundaryMesh().findPatchID(outerSliderName)];
111 labelList osf(outerSlider.size());
115 osf[i] = outerSlider.start() + i;
122 boolList(outerSlider.size(), false),
127 // Add empty zone for cut faces
137 List<cellZone*> cz(1);
139 // Mark every cell with its topological region
140 regionSplit rs(*this);
142 // Get the region of the cell containing the origin.
143 label originRegion = rs[findNearestCell(cs().origin())];
145 labelList movingCells(nCells());
146 label nMovingCells = 0;
150 if (rs[cellI] == originRegion)
152 movingCells[nMovingCells] = cellI;
157 movingCells.setSize(nMovingCells);
158 Info << "Number of cells in the moving region: " << nMovingCells << endl;
168 Info << "Adding point, face and cell zones" << endl;
169 addZones(pz, fz, cz);
171 // Add a topology modifier
172 Info << "Adding topology modifiers" << endl;
173 topoChanger_.setSize(1);
182 outerSliderName + "Zone",
183 innerSliderName + "Zone",
188 slidingInterface::INTEGRAL
191 topoChanger_.writeOpt() = IOobject::AUTO_WRITE;
197 void Foam::mixerFvMesh::calcMovingMasks() const
201 Info<< "void mixerFvMesh::calcMovingMasks() const : "
202 << "Calculating point and cell masks"
206 if (movingPointsMaskPtr_)
208 FatalErrorIn("void mixerFvMesh::calcMovingMasks() const")
209 << "point mask already calculated"
210 << abort(FatalError);
213 // Set the point mask
214 movingPointsMaskPtr_ = new scalarField(points().size(), 0);
215 scalarField& movingPointsMask = *movingPointsMaskPtr_;
217 const cellList& c = cells();
218 const faceList& f = faces();
220 const labelList& cellAddr =
221 cellZones()[cellZones().findZoneID("movingCells")];
223 forAll (cellAddr, cellI)
225 const cell& curCell = c[cellAddr[cellI]];
227 forAll (curCell, faceI)
229 // Mark all the points as moving
230 const face& curFace = f[curCell[faceI]];
232 forAll (curFace, pointI)
234 movingPointsMask[curFace[pointI]] = 1;
239 const word innerSliderZoneName
241 word(motionDict_.subDict("slider").lookup("inside"))
245 const labelList& innerSliderAddr =
246 faceZones()[faceZones().findZoneID(innerSliderZoneName)];
248 forAll (innerSliderAddr, faceI)
250 const face& curFace = f[innerSliderAddr[faceI]];
252 forAll (curFace, pointI)
254 movingPointsMask[curFace[pointI]] = 1;
258 const word outerSliderZoneName
260 word(motionDict_.subDict("slider").lookup("outside"))
264 const labelList& outerSliderAddr =
265 faceZones()[faceZones().findZoneID(outerSliderZoneName)];
267 forAll (outerSliderAddr, faceI)
269 const face& curFace = f[outerSliderAddr[faceI]];
271 forAll (curFace, pointI)
273 movingPointsMask[curFace[pointI]] = 0;
279 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
281 // Construct from components
282 Foam::mixerFvMesh::mixerFvMesh
287 topoChangerFvMesh(io),
300 ).subDict(typeName + "Coeffs")
304 coordinateSystem::New
307 motionDict_.subDict("coordinateSystem")
310 rpm_(readScalar(motionDict_.lookup("rpm"))),
311 movingPointsMaskPtr_(NULL)
313 addZonesAndModifiers();
315 Info<< "Mixer mesh:" << nl
316 << " origin: " << cs().origin() << nl
317 << " axis: " << cs().axis() << nl
318 << " rpm: " << rpm_ << endl;
322 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
324 Foam::mixerFvMesh::~mixerFvMesh()
326 deleteDemandDrivenData(movingPointsMaskPtr_);
329 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
331 // Return moving points mask. Moving points marked with 1
332 const Foam::scalarField& Foam::mixerFvMesh::movingPointsMask() const
334 if (!movingPointsMaskPtr_)
339 return *movingPointsMaskPtr_;
343 bool Foam::mixerFvMesh::update()
345 // Rotational speed needs to be converted from rpm
348 csPtr_->globalPosition
350 csPtr_->localPosition(points())
351 + vector(0, rpm_*360.0*time().deltaT().value()/60.0, 0)
356 // Make changes. Use inflation (so put new points in topoChangeMap)
357 autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh(true);
359 if (topoChangeMap.valid())
363 Info << "Mesh topology is changing" << endl;
366 deleteDemandDrivenData(movingPointsMaskPtr_);
371 csPtr_->globalPosition
373 csPtr_->localPosition(oldPoints())
374 + vector(0, rpm_*360.0*time().deltaT().value()/60.0, 0)
383 // ************************************************************************* //