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
25 \*---------------------------------------------------------------------------*/
27 #include "mixerFvMesh.H"
29 #include "regionSplit.H"
30 #include "slidingInterface.H"
31 #include "mapPolyMesh.H"
33 #include "addToRunTimeSelectionTable.H"
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 defineTypeNameAndDebug(mixerFvMesh, 0);
40 addToRunTimeSelectionTable(topoChangerFvMesh, mixerFvMesh, IOobject);
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 void Foam::mixerFvMesh::addZonesAndModifiers()
48 // Add zones and modifiers for motion action
52 pointZones().size() > 0
53 || faceZones().size() > 0
54 || cellZones().size() > 0
57 Info<< "void mixerFvMesh::addZonesAndModifiers() : "
58 << "Zones and modifiers already present. Skipping."
61 if (topoChanger_.size() == 0)
65 "void mixerFvMesh::addZonesAndModifiers()"
66 ) << "Mesh modifiers not read properly"
73 Info<< "Time = " << time().timeName() << endl
74 << "Adding zones and modifiers to the mesh" << endl;
77 List<pointZone*> pz(1);
79 // Add an empty zone for cut points
90 // Do face zones for slider
92 List<faceZone*> fz(3);
95 const word movingSliderName(dict_.subDict("slider").lookup("moving"));
96 label movingSliderIndex = boundaryMesh().findPatchID(movingSliderName);
98 if (movingSliderIndex < 0)
100 FatalErrorIn("void mixerFvMesh::addZonesAndModifiers() const")
101 << "Moving slider patch not found in boundary"
102 << abort(FatalError);
105 const word staticSliderName(dict_.subDict("slider").lookup("static"));
106 label staticSliderIndex = boundaryMesh().findPatchID(staticSliderName);
108 if (staticSliderIndex < 0)
110 FatalErrorIn("void mixerFvMesh::addZonesAndModifiers() const")
111 << "Static slider patch not found in boundary"
112 << abort(FatalError);
116 const polyPatch& movingSlider = boundaryMesh()[movingSliderIndex];
118 labelList isf(movingSlider.size());
122 isf[i] = movingSlider.start() + i;
127 movingSliderName + "Zone",
129 boolList(movingSlider.size(), false),
135 const polyPatch& staticSlider = boundaryMesh()[staticSliderIndex];
137 labelList osf(staticSlider.size());
141 osf[i] = staticSlider.start() + i;
146 staticSliderName + "Zone",
148 boolList(staticSlider.size(), false),
153 // Add empty zone for cut faces
163 List<cellZone*> cz(1);
165 // Mark every cell with its topological region
166 regionSplit rs(*this);
168 // Get the region of the cell containing the origin.
169 label originRegion = rs[findNearestCell(rotatingRegionMarker_)];
171 labelList movingCells(nCells());
172 label nMovingCells = 0;
176 if (rs[cellI] == originRegion)
178 movingCells[nMovingCells] = cellI;
183 movingCells.setSize(nMovingCells);
184 Info << "Number of cells in the moving region: " << nMovingCells << endl;
194 Info << "Adding point, face and cell zones" << endl;
195 addZones(pz, fz, cz);
197 // Add a topology modifier
198 Info << "Adding topology modifiers" << endl;
199 topoChanger_.setSize(1);
208 staticSliderName + "Zone",
209 movingSliderName + "Zone",
214 slidingInterface::INTEGRAL, // Edge matching algorithm
215 attachDetach_, // Attach-detach action
216 intersection::VISIBLE // Projection algorithm
220 // Write mesh and modifiers
221 topoChanger_.writeOpt() = IOobject::AUTO_WRITE;
222 topoChanger_.write();
227 bool Foam::mixerFvMesh::attached() const
229 return refCast<const slidingInterface>(topoChanger_[0]).attached();
233 void Foam::mixerFvMesh::calcMovingMask() const
237 Info<< "void mixerFvMesh::calcMovingMask() const : "
238 << "Calculating point and cell masks"
242 if (movingPointsMaskPtr_)
244 FatalErrorIn("void mixerFvMesh::calcMovingMask() const")
245 << "point mask already calculated"
246 << abort(FatalError);
249 // Set the point mask
250 movingPointsMaskPtr_ = new scalarField(allPoints().size(), 0);
251 scalarField& movingPointsMask = *movingPointsMaskPtr_;
253 const cellList& c = cells();
254 const faceList& f = allFaces();
256 const labelList& cellAddr =
257 cellZones()[cellZones().findZoneID("movingCells")];
259 forAll (cellAddr, cellI)
261 const cell& curCell = c[cellAddr[cellI]];
263 forAll (curCell, faceI)
265 // Mark all the points as moving
266 const face& curFace = f[curCell[faceI]];
268 forAll (curFace, pointI)
270 movingPointsMask[curFace[pointI]] = 1;
275 const word movingSliderZoneName
277 word(dict_.subDict("slider").lookup("moving"))
281 const labelList& movingSliderAddr =
282 faceZones()[faceZones().findZoneID(movingSliderZoneName)];
284 forAll (movingSliderAddr, faceI)
286 const face& curFace = f[movingSliderAddr[faceI]];
288 forAll (curFace, pointI)
290 movingPointsMask[curFace[pointI]] = 1;
294 const word staticSliderZoneName
296 word(dict_.subDict("slider").lookup("static"))
300 const labelList& staticSliderAddr =
301 faceZones()[faceZones().findZoneID(staticSliderZoneName)];
303 forAll (staticSliderAddr, faceI)
305 const face& curFace = f[staticSliderAddr[faceI]];
307 forAll (curFace, pointI)
309 movingPointsMask[curFace[pointI]] = 0;
315 // Return moving points mask. Moving points marked with 1
316 const Foam::scalarField& Foam::mixerFvMesh::movingPointsMask() const
318 if (!movingPointsMaskPtr_)
323 return *movingPointsMaskPtr_;
327 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
329 // Construct from components
330 Foam::mixerFvMesh::mixerFvMesh
335 topoChangerFvMesh(io),
348 ).subDict(typeName + "Coeffs")
352 coordinateSystem::New
355 dict_.subDict("coordinateSystem")
358 rpm_(readScalar(dict_.lookup("rpm"))),
359 rotatingRegionMarker_
361 dict_.lookupOrDefault<point>("rotatingRegionMarker", csPtr_->origin())
363 attachDetach_(dict_.lookupOrDefault<bool>("attachDetach", false)),
364 movingPointsMaskPtr_(NULL)
366 Info<< "Rotating region marker point: " << rotatingRegionMarker_
367 << " Attach-detach action = " << attachDetach_ << endl;
369 addZonesAndModifiers();
371 Info<< "Mixer mesh" << nl
372 << " origin : " << cs().origin() << nl
373 << " axis : " << cs().axis() << nl
374 << " rpm : " << rpm_ << nl
375 << " marker : " << rotatingRegionMarker_ << nl
376 << " attach-detach: " << attachDetach_
381 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
383 Foam::mixerFvMesh::~mixerFvMesh()
385 deleteDemandDrivenData(movingPointsMaskPtr_);
389 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
391 bool Foam::mixerFvMesh::update()
393 if (attached() && attachDetach_)
396 Info << "Detaching rotors" << endl;
397 autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh();
399 deleteDemandDrivenData(movingPointsMaskPtr_);
403 pointField oldPointsNew = allPoints();
405 // Move points. Rotational speed needs to be converted from rpm
408 csPtr_->globalPosition
410 csPtr_->localPosition(allPoints())
411 + vector(0, rpm_*360.0*time().deltaT().value()/60.0, 0)
416 autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh();
418 if (topoChangeMap->morphing())
420 Info << "Attaching rotors" << endl;
422 deleteDemandDrivenData(movingPointsMaskPtr_);
424 // Move the sliding interface points to correct position
425 pointField mappedOldPointsNew(allPoints().size());
426 mappedOldPointsNew.map(oldPointsNew, topoChangeMap->pointMap());
428 movePoints(mappedOldPointsNew);
432 // Move the sliding interface points to correct position
433 movePoints(topoChangeMap->preMotionPoints());
436 return topoChangeMap->morphing();
440 // ************************************************************************* //