Merge branch 'upstream/OpenFOAM' into master
[freefoam.git] / src / topoChangerFvMesh / mixerFvMesh / mixerFvMesh.C
blobf20c5e8ea52a7b69ef080afa6a8c4560639fdeaa
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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"
28 #include <OpenFOAM/Time.H>
29 #include <meshTools/regionSplit.H>
30 #include <dynamicMesh/slidingInterface.H>
31 #include <OpenFOAM/addToRunTimeSelectionTable.H>
32 #include <OpenFOAM/mapPolyMesh.H>
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 namespace Foam
38     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
50     if
51     (
52         pointZones().size()
53      || faceZones().size()
54      || cellZones().size()
55      || topoChanger_.size()
56     )
57     {
58         Info<< "void mixerFvMesh::addZonesAndModifiers() : "
59             << "Zones and modifiers already present.  Skipping."
60             << endl;
62         return;
63     }
65     Info<< "Time = " << time().timeName() << endl
66         << "Adding zones and modifiers to the mesh" << endl;
68     // Add zones
69     List<pointZone*> pz(1);
71     // Add an empty zone for cut points
73     pz[0] = new pointZone
74     (
75         "cutPointZone",
76         labelList(0),
77         0,
78         pointZones()
79     );
82     // Do face zones for slider
84     List<faceZone*> fz(3);
86     // Inner slider
87     const word innerSliderName(motionDict_.subDict("slider").lookup("inside"));
88     const polyPatch& innerSlider =
89         boundaryMesh()[boundaryMesh().findPatchID(innerSliderName)];
91     labelList isf(innerSlider.size());
93     forAll (isf, i)
94     {
95         isf[i] = innerSlider.start() + i;
96     }
98     fz[0] = new faceZone
99     (
100         "insideSliderZone",
101         isf,
102         boolList(innerSlider.size(), false),
103         0,
104         faceZones()
105     );
107     // Outer slider
108     const word outerSliderName(motionDict_.subDict("slider").lookup("outside"));
109     const polyPatch& outerSlider =
110         boundaryMesh()[boundaryMesh().findPatchID(outerSliderName)];
112     labelList osf(outerSlider.size());
114     forAll (osf, i)
115     {
116         osf[i] = outerSlider.start() + i;
117     }
119     fz[1] = new faceZone
120     (
121         "outsideSliderZone",
122         osf,
123         boolList(outerSlider.size(), false),
124         1,
125         faceZones()
126     );
128     // Add empty zone for cut faces
129     fz[2] = new faceZone
130     (
131         "cutFaceZone",
132         labelList(0),
133         boolList(0, false),
134         2,
135         faceZones()
136     );
138     List<cellZone*> cz(1);
140     // Mark every cell with its topological region
141     regionSplit rs(*this);
143     // Get the region of the cell containing the origin.
144     label originRegion = rs[findNearestCell(cs().origin())];
146     labelList movingCells(nCells());
147     label nMovingCells = 0;
149     forAll(rs, cellI)
150     {
151         if (rs[cellI] == originRegion)
152         {
153             movingCells[nMovingCells] = cellI;
154             nMovingCells++;
155         }
156     }
158     movingCells.setSize(nMovingCells);
159     Info << "Number of cells in the moving region: " << nMovingCells << endl;
161     cz[0] = new cellZone
162     (
163         "movingCells",
164         movingCells,
165         0,
166         cellZones()
167     );
169     Info << "Adding point, face and cell zones" << endl;
170     addZones(pz, fz, cz);
172     // Add a topology modifier
173     Info << "Adding topology modifiers" << endl;
174     topoChanger_.setSize(1);
175     topoChanger_.set
176     (
177         0,
178         new slidingInterface
179         (
180             "mixerSlider",
181             0,
182             topoChanger_,
183             outerSliderName + "Zone",
184             innerSliderName + "Zone",
185             "cutPointZone",
186             "cutFaceZone",
187             outerSliderName,
188             innerSliderName,
189             slidingInterface::INTEGRAL
190         )
191     );
192     topoChanger_.writeOpt() = IOobject::AUTO_WRITE;
194     write();
198 void Foam::mixerFvMesh::calcMovingMasks() const
200     if (debug)
201     {
202         Info<< "void mixerFvMesh::calcMovingMasks() const : "
203             << "Calculating point and cell masks"
204             << endl;
205     }
207     if (movingPointsMaskPtr_)
208     {
209         FatalErrorIn("void mixerFvMesh::calcMovingMasks() const")
210             << "point mask already calculated"
211             << abort(FatalError);
212     }
214     // Set the point mask
215     movingPointsMaskPtr_ = new scalarField(points().size(), 0);
216     scalarField& movingPointsMask = *movingPointsMaskPtr_;
218     const cellList& c = cells();
219     const faceList& f = faces();
221     const labelList& cellAddr =
222         cellZones()[cellZones().findZoneID("movingCells")];
224     forAll (cellAddr, cellI)
225     {
226         const cell& curCell = c[cellAddr[cellI]];
228         forAll (curCell, faceI)
229         {
230             // Mark all the points as moving
231             const face& curFace = f[curCell[faceI]];
233             forAll (curFace, pointI)
234             {
235                 movingPointsMask[curFace[pointI]] = 1;
236             }
237         }
238     }
240     const word innerSliderZoneName
241     (
242         word(motionDict_.subDict("slider").lookup("inside"))
243       + "Zone"
244     );
246     const labelList& innerSliderAddr =
247         faceZones()[faceZones().findZoneID(innerSliderZoneName)];
249     forAll (innerSliderAddr, faceI)
250     {
251         const face& curFace = f[innerSliderAddr[faceI]];
253         forAll (curFace, pointI)
254         {
255             movingPointsMask[curFace[pointI]] = 1;
256         }
257     }
259     const word outerSliderZoneName
260     (
261         word(motionDict_.subDict("slider").lookup("outside"))
262       + "Zone"
263     );
265     const labelList& outerSliderAddr =
266         faceZones()[faceZones().findZoneID(outerSliderZoneName)];
268     forAll (outerSliderAddr, faceI)
269     {
270         const face& curFace = f[outerSliderAddr[faceI]];
272         forAll (curFace, pointI)
273         {
274             movingPointsMask[curFace[pointI]] = 0;
275         }
276     }
280 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
282 // Construct from components
283 Foam::mixerFvMesh::mixerFvMesh
285     const IOobject& io
288     topoChangerFvMesh(io),
289     motionDict_
290     (
291         IOdictionary
292         (
293             IOobject
294             (
295                 "dynamicMeshDict",
296                 time().constant(),
297                 *this,
298                 IOobject::MUST_READ,
299                 IOobject::NO_WRITE
300             )
301         ).subDict(typeName + "Coeffs")
302     ),
303     csPtr_
304     (
305         coordinateSystem::New
306         (
307             "coordinateSystem",
308             motionDict_.subDict("coordinateSystem")
309         )
310     ),
311     rpm_(readScalar(motionDict_.lookup("rpm"))),
312     movingPointsMaskPtr_(NULL)
314     addZonesAndModifiers();
316     Info<< "Mixer mesh:" << nl
317         << "    origin: " << cs().origin() << nl
318         << "    axis: " << cs().axis() << nl
319         << "    rpm: " << rpm_ << endl;
323 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
325 Foam::mixerFvMesh::~mixerFvMesh()
327     deleteDemandDrivenData(movingPointsMaskPtr_);
330 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
332 // Return moving points mask.  Moving points marked with 1
333 const Foam::scalarField& Foam::mixerFvMesh::movingPointsMask() const
335     if (!movingPointsMaskPtr_)
336     {
337         calcMovingMasks();
338     }
340     return *movingPointsMaskPtr_;
344 bool Foam::mixerFvMesh::update()
346      // Rotational speed needs to be converted from rpm
347     movePoints
348     (
349         csPtr_->globalPosition
350         (
351             csPtr_->localPosition(points())
352           + vector(0, rpm_*360.0*time().deltaT().value()/60.0, 0)
353             *movingPointsMask()
354         )
355     );
357     // Make changes. Use inflation (so put new points in topoChangeMap)
358     autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh(true);
360     if (topoChangeMap.valid())
361     {
362         if (debug)
363         {
364             Info << "Mesh topology is changing" << endl;
365         }
367         deleteDemandDrivenData(movingPointsMaskPtr_);
368     }
370     movePoints
371     (
372         csPtr_->globalPosition
373         (
374             csPtr_->localPosition(oldPoints())
375           + vector(0, rpm_*360.0*time().deltaT().value()/60.0, 0)
376             *movingPointsMask()
377         )
378     );
380     return true;
384 // ************************ vim: set sw=4 sts=4 et: ************************ //