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
49 if (cellZones().size() > 0)
51 Info<< "void mixerFvMesh::addZonesAndModifiers() : "
52 << "Zones and modifiers already present. Skipping."
55 // Check definition of the modifier
58 pointZones().size() > 0
59 || faceZones().size() > 0
62 if (topoChanger_.size() == 0)
66 "void mixerFvMesh::addZonesAndModifiers()"
67 ) << "Mesh modifiers not read properly. "
68 << "pointZones = " << pointZones().size()
69 << " faceZones = " << faceZones().size()
77 // Find patches and check sizes
80 const word movingSliderName(dict_.subDict("slider").lookup("moving"));
81 label movingSliderIndex = boundaryMesh().findPatchID(movingSliderName);
83 if (movingSliderIndex < 0)
85 FatalErrorIn("void mixerFvMesh::addZonesAndModifiers() const")
86 << "Moving slider patch not found in boundary"
90 const word staticSliderName(dict_.subDict("slider").lookup("static"));
91 label staticSliderIndex = boundaryMesh().findPatchID(staticSliderName);
93 if (staticSliderIndex < 0)
95 FatalErrorIn("void mixerFvMesh::addZonesAndModifiers() const")
96 << "Static slider patch not found in boundary"
101 const polyPatch& movingSlider = boundaryMesh()[movingSliderIndex];
102 const polyPatch& staticSlider = boundaryMesh()[staticSliderIndex];
107 bool addSlider = false;
109 if (movingSlider.size() > 0 && staticSlider.size() > 0)
116 Info<< "Time = " << time().timeName() << endl
117 << "Adding zones and modifiers to the mesh" << endl;
121 // Add an empty zone for cut points
123 pz[0] = new pointZone
132 // Do face zones for slider
135 labelList isf(movingSlider.size());
139 isf[i] = movingSlider.start() + i;
144 movingSliderName + "Zone",
146 boolList(movingSlider.size(), false),
152 labelList osf(staticSlider.size());
156 osf[i] = staticSlider.start() + i;
161 staticSliderName + "Zone",
163 boolList(staticSlider.size(), false),
168 // Add empty zone for cut faces
179 // Add cell zone even if slider does not exist
181 List<cellZone*> cz(1);
183 // Mark every cell with its topological region
184 regionSplit rs(*this);
186 // Get the region of the cell containing the origin.
187 label originRegion = rs[findNearestCell(rotatingRegionMarker_)];
189 labelList movingCells(nCells());
190 label nMovingCells = 0;
194 if (rs[cellI] == originRegion)
196 movingCells[nMovingCells] = cellI;
201 movingCells.setSize(nMovingCells);
202 Info << "Number of cells in the moving region: " << nMovingCells << endl;
212 Info << "Adding point, face and cell zones" << endl;
213 addZones(pz, fz, cz);
217 // Add a topology modifier
218 Info << "Adding topology modifiers" << endl;
219 topoChanger_.setSize(1);
229 staticSliderName + "Zone",
230 movingSliderName + "Zone",
235 slidingInterface::INTEGRAL, // Edge matching algorithm
236 attachDetach_, // Attach-detach action
237 intersection::VISIBLE // Projection algorithm
242 // Write mesh and modifiers
243 topoChanger_.writeOpt() = IOobject::AUTO_WRITE;
244 topoChanger_.write();
249 bool Foam::mixerFvMesh::attached() const
253 if (topoChanger_.size() > 0)
255 result = refCast<const slidingInterface>(topoChanger_[0]).attached();
258 reduce(result, orOp<bool>());
264 void Foam::mixerFvMesh::calcMovingMask() const
268 Info<< "void mixerFvMesh::calcMovingMask() const : "
269 << "Calculating point and cell masks"
273 if (movingPointsMaskPtr_)
275 FatalErrorIn("void mixerFvMesh::calcMovingMask() const")
276 << "point mask already calculated"
277 << abort(FatalError);
280 // Set the point mask
281 movingPointsMaskPtr_ = new scalarField(allPoints().size(), 0);
282 scalarField& movingPointsMask = *movingPointsMaskPtr_;
284 const cellList& c = cells();
285 const faceList& f = allFaces();
287 const labelList& cellAddr =
288 cellZones()[cellZones().findZoneID("movingCells")];
290 forAll (cellAddr, cellI)
292 const cell& curCell = c[cellAddr[cellI]];
294 forAll (curCell, faceI)
296 // Mark all the points as moving
297 const face& curFace = f[curCell[faceI]];
299 forAll (curFace, pointI)
301 movingPointsMask[curFace[pointI]] = 1;
306 if (topoChanger_.size() > 0)
308 // Topo changer present. Use zones for motion
309 const word movingSliderZoneName
311 word(dict_.subDict("slider").lookup("moving")) + "Zone"
314 const labelList& movingSliderAddr =
315 faceZones()[faceZones().findZoneID(movingSliderZoneName)];
317 forAll (movingSliderAddr, faceI)
319 const face& curFace = f[movingSliderAddr[faceI]];
321 forAll (curFace, pointI)
323 movingPointsMask[curFace[pointI]] = 1;
327 const word staticSliderZoneName
329 word(dict_.subDict("slider").lookup("static")) + "Zone"
332 const labelList& staticSliderAddr =
333 faceZones()[faceZones().findZoneID(staticSliderZoneName)];
335 forAll (staticSliderAddr, faceI)
337 const face& curFace = f[staticSliderAddr[faceI]];
339 forAll (curFace, pointI)
341 movingPointsMask[curFace[pointI]] = 0;
348 // Return moving points mask. Moving points marked with 1
349 const Foam::scalarField& Foam::mixerFvMesh::movingPointsMask() const
351 if (!movingPointsMaskPtr_)
356 return *movingPointsMaskPtr_;
360 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
362 // Construct from components
363 Foam::mixerFvMesh::mixerFvMesh
368 topoChangerFvMesh(io),
381 ).subDict(typeName + "Coeffs")
386 dict_.subDict("coordinateSystem")
388 rpm_(readScalar(dict_.lookup("rpm"))),
389 rotatingRegionMarker_
391 dict_.lookupOrDefault<point>("rotatingRegionMarker", cs_.origin())
393 attachDetach_(dict_.lookupOrDefault<bool>("attachDetach", false)),
394 movingPointsMaskPtr_(NULL)
396 // Make sure the coordinate system does not operate in degrees
397 // Bug fix, HJ, 3/Oct/2011
398 if (!cs_.inDegrees())
400 WarningIn("mixerFvMesh::mixerFvMesh(const IOobject& io)")
401 << "Mixer coordinate system is set to operate in radians. "
402 << "Changing to rad for correct calculation of angular velocity."
404 << "To remove this message please add entry" << nl << nl
405 << "inDegrees true;" << nl << nl
406 << "to the specification of the coordinate system"
409 cs_.inDegrees() = true;
412 Info<< "Rotating region marker point: " << rotatingRegionMarker_
413 << " Attach-detach action = " << attachDetach_ << endl;
415 addZonesAndModifiers();
417 Info<< "Mixer mesh" << nl
418 << " origin : " << cs().origin() << nl
419 << " axis : " << cs().axis() << nl
420 << " rpm : " << rpm_ << nl
421 << " marker : " << rotatingRegionMarker_ << nl
422 << " attach-detach: " << attachDetach_
427 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
429 Foam::mixerFvMesh::~mixerFvMesh()
431 deleteDemandDrivenData(movingPointsMaskPtr_);
435 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
437 bool Foam::mixerFvMesh::update()
439 if (attached() && attachDetach_)
442 Info << "Detaching rotors" << endl;
443 autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh();
445 deleteDemandDrivenData(movingPointsMaskPtr_);
449 pointField oldPointsNew = allPoints();
451 // Move points. Rotational speed needs to be converted from rpm
456 cs_.localPosition(allPoints())
457 + vector(0, rpm_*360.0*time().deltaT().value()/60.0, 0)
462 autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh();
464 bool localMorphing = topoChangeMap->morphing();
465 bool globalMorphing = localMorphing;
467 reduce(globalMorphing, orOp<bool>());
471 Info << "Attaching rotors" << endl;
473 deleteDemandDrivenData(movingPointsMaskPtr_);
475 // Move the sliding interface points to correct position
478 pointField mappedOldPointsNew(allPoints().size());
479 mappedOldPointsNew.map(oldPointsNew, topoChangeMap->pointMap());
481 movePoints(mappedOldPointsNew);
486 // Move the sliding interface points to correct position
487 movePoints(topoChangeMap->preMotionPoints());
491 pointField newPoints = allPoints();
492 movePoints(oldPointsNew);
497 // Move the sliding interface points to correct position
498 movePoints(newPoints);
502 return topoChangeMap->morphing();
506 // ************************************************************************* //