Forward compatibility: flex
[foam-extend-3.2.git] / src / dynamicMesh / topoChangerFvMesh / multiMixerFvMesh / mixerRotor.C
blobf66493a06ede712b1647f7b34749bd3efdc58007
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "mixerRotor.H"
27 #include "regionSplit.H"
28 #include "polyTopoChanger.H"
29 #include "slidingInterface.H"
30 #include "foamTime.H"
32 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
34 void Foam::mixerRotor::addZones
36     DynamicList<pointZone*>& pz,
37     DynamicList<faceZone*>& fz,
38     DynamicList<cellZone*>& cz,
39     const regionSplit& rs
42     // Get the region of the cell containing the origin.
43     label originRegion = rs[mesh_.findNearestCell(rotatingRegionMarker_)];
45     labelList movingCells(mesh_.nCells());
46     label nMovingCells = 0;
48     forAll(rs, cellI)
49     {
50         if (rs[cellI] == originRegion)
51         {
52             movingCells[nMovingCells] = cellI;
53             nMovingCells++;
54         }
55     }
57     movingCells.setSize(nMovingCells);
58     Info<< "Number of moving cells for " << name_ << ": "
59         << nMovingCells << endl;
61     cz.append
62     (
63         new cellZone
64         (
65             "movingCellsZone" + name_,
66             movingCells,
67             cz.size(),
68             mesh_.cellZones()
69         )
70     );
72     if (useTopoSliding_)
73     {
74         // Add an empty zone for cut points
75         pz.append
76         (
77             new pointZone
78             (
79                 "cutPointZone" + name_,
80                 labelList(0),
81                 pz.size(),
82                 mesh_.pointZones()
83             )
84         );
87         // Do face zones for slider
89         // Moving slider
90         label movingSliderIndex =
91             mesh_.boundaryMesh().findPatchID(movingSliderName_);
93         if (movingSliderIndex < 0)
94         {
95             FatalErrorIn("void mixerRotor::addZones(...) const")
96                 << "Moving slider patch not found in boundary"
97                 << abort(FatalError);
98         }
100         label staticSliderIndex =
101             mesh_.boundaryMesh().findPatchID(staticSliderName_);
103         if (staticSliderIndex < 0)
104         {
105             FatalErrorIn("void mixerRotor::addZones(...) const")
106                 << "Static slider patch not found in boundary"
107                 << abort(FatalError);
108         }
110         const polyPatch& movingSlider =
111             mesh_.boundaryMesh()[movingSliderIndex];
113         labelList isf(movingSlider.size());
115         forAll (isf, i)
116         {
117             isf[i] = movingSlider.start() + i;
118         }
120         fz.append
121         (
122             new faceZone
123             (
124                 movingSliderName_ + "Zone" + name_,
125                 isf,
126                 boolList(movingSlider.size(), false),
127                 fz.size(),
128                 mesh_.faceZones()
129             )
130         );
132         // Static slider
133         const polyPatch& staticSlider =
134             mesh_.boundaryMesh()[staticSliderIndex];
136         labelList osf(staticSlider.size());
138         forAll (osf, i)
139         {
140             osf[i] = staticSlider.start() + i;
141         }
143         fz.append
144         (
145             new faceZone
146             (
147                 staticSliderName_ + "Zone" + name_,
148                 osf,
149                 boolList(staticSlider.size(), false),
150                 fz.size(),
151                 mesh_.faceZones()
152             )
153         );
155         // Add empty zone for cut faces
156         fz.append
157         (
158             new faceZone
159             (
160                 "cutFaceZone" + name_,
161                 labelList(0),
162                 boolList(0, false),
163                 fz.size(),
164                 mesh_.faceZones()
165             )
166         );
167     }
171 void Foam::mixerRotor::addModifiers
173     polyTopoChanger& tc,
174     label& nextI
177     // Add a topology modifier
178     if (useTopoSliding())
179     {
180         Info << "Adding topology modifier for rotor " << name_ << endl;
182         tc.set
183         (
184             nextI,
185             new slidingInterface
186             (
187                 "mixerSlider"  + name_,
188                 nextI,
189                 tc,
190                 staticSliderName_ + "Zone"  + name_,
191                 movingSliderName_ + "Zone" + name_,
192                 "cutPointZone" + name_,
193                 "cutFaceZone" + name_,
194                 staticSliderName_,
195                 movingSliderName_,
196                 slidingInterface::INTEGRAL,   // Edge matching algorithm
197                 attachDetach_,                // Attach-detach action
198                 intersection::VISIBLE         // Projection algorithm
199             )
200         );
202         nextI++;
203     }
207 void Foam::mixerRotor::calcMovingMask() const
209     if (movingPointsMaskPtr_)
210     {
211         FatalErrorIn("void mixerRotor::calcMovingMask() const")
212             << "point mask already calculated"
213             << abort(FatalError);
214     }
216     // Set the point mask
217     movingPointsMaskPtr_ = new scalarField(mesh_.allPoints().size(), 0);
218     scalarField& movingPointsMask = *movingPointsMaskPtr_;
220     const cellList& c = mesh_.cells();
221     const faceList& f = mesh_.allFaces();
223     const labelList& cellAddr = mesh_.cellZones()
224         [mesh_.cellZones().findZoneID("movingCellsZone" + name_)];
226     forAll (cellAddr, cellI)
227     {
228         const cell& curCell = c[cellAddr[cellI]];
230         forAll (curCell, faceI)
231         {
232             // Mark all the points as moving
233             const face& curFace = f[curCell[faceI]];
235             forAll (curFace, pointI)
236             {
237                 movingPointsMask[curFace[pointI]] = 1;
238             }
239         }
240     }
242     // Attempt to enforce motion on sliders if zones exist
243     const label msI =
244         mesh_.faceZones().findZoneID(movingSliderName_ + "Zone" + name_);
246     if (msI > -1)
247     {
248         const labelList& movingSliderAddr = mesh_.faceZones()[msI];
250         forAll (movingSliderAddr, faceI)
251         {
252             const face& curFace = f[movingSliderAddr[faceI]];
254             forAll (curFace, pointI)
255             {
256                 movingPointsMask[curFace[pointI]] = 1;
257             }
258         }
259     }
261     const label ssI =
262         mesh_.faceZones().findZoneID(staticSliderName_ + "Zone" + name_);
264     if (ssI > -1)
265     {
266         const labelList& staticSliderAddr = mesh_.faceZones()[ssI];
268         forAll (staticSliderAddr, faceI)
269         {
270             const face& curFace = f[staticSliderAddr[faceI]];
272             forAll (curFace, pointI)
273             {
274                 movingPointsMask[curFace[pointI]] = 0;
275             }
276         }
277     }
281 // Return moving points mask.  Moving points marked with 1
282 const Foam::scalarField& Foam::mixerRotor::movingPointsMask() const
284     if (!movingPointsMaskPtr_)
285     {
286         calcMovingMask();
287     }
289     return *movingPointsMaskPtr_;
293 void Foam::mixerRotor::clearPointMask()
295     deleteDemandDrivenData(movingPointsMaskPtr_);
299 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
301 Foam::mixerRotor::mixerRotor
303     const word& name,
304     const polyMesh& mesh,
305     const dictionary& dict
308     name_(name),
309     mesh_(mesh),
310     cs_
311     (
312         "coordinateSystem",
313         dict.subDict("coordinateSystem")
314     ),
315     rpm_(readScalar(dict.lookup("rpm"))),
316     movingSliderName_(dict.lookup("movingPatch")),
317     staticSliderName_(dict.lookup("staticPatch")),
318     rotatingRegionMarker_
319     (
320         dict.lookupOrDefault<point>("rotatingRegionMarker", cs_.origin())
321     ),
322     invertMotionMask_
323     (
324         dict.lookupOrDefault<bool>("invertMotionMask", false)
325     ),
326     useTopoSliding_(dict.lookup("useTopoSliding")),
327     attachDetach_(dict.lookupOrDefault<bool>("attachDetach", true)),
328     movingPointsMaskPtr_(NULL)
330     // Make sure the coordinate system does not operate in degrees
331     // Bug fix, HJ, 3/Oct/2011
332     if (!cs_.inDegrees())
333     {
334         WarningIn
335         (
336             "mixerRotor::mixerRotor\n"
337             "(\n"
338             "    const word& name,\n"
339             "    const polyMesh& mesh,\n"
340             "    const dictionary& dict\n"
341             ")"
342         )   << "Mixer coordinate system is set to operate in radians.  "
343             << "Changing to rad for correct calculation of angular velocity."
344             << nl
345             << "To remove this message please add entry" << nl << nl
346             << "inDegrees true;" << nl << nl
347             << "to the specification of the coordinate system"
348             << endl;
350         cs_.inDegrees() = true;
351     }
353     Info<< "Rotor " << name << ":" << nl
354         << "    origin      : " << cs().origin() << nl
355         << "    axis        : " << cs().axis() << nl
356         << "    rpm         : " << rpm_ << nl
357         << "    invert mask : " << invertMotionMask_ << nl
358         << "    topo sliding: " << useTopoSliding_ << endl;
360     if (useTopoSliding_)
361     {
362         Info<< "    attach-detach: " << attachDetach_ << endl;
363     }
367 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
369 Foam::mixerRotor::~mixerRotor()
371     clearPointMask();
375 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
377 Foam::tmp<Foam::vectorField> Foam::mixerRotor::pointMotion() const
379     // Rotational speed needs to be converted from rpm
380     scalarField mpm = movingPointsMask();
382     if (invertMotionMask_)
383     {
384         Info << "Inverting motion mask" << endl;
385         mpm = 1 - mpm;
386     }
388     return cs_.globalPosition
389     (
390         // Motion vector in cylindrical coordinate system (x theta z)
391         cs_.localPosition(mesh_.allPoints())
392       + vector(0, rpm_*360.0*mesh_.time().deltaT().value()/60.0, 0)*mpm
393     ) - mesh_.allPoints();
397 void Foam::mixerRotor::updateTopology()
399     clearPointMask();
403 // ************************************************************************* //