mixerFvMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "mixerFvMesh.H"
27 #include "Time.H"
28 #include "regionSplit.H"
29 #include "slidingInterface.H"
31 #include "mapPolyMesh.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(mixerFvMesh, 0);
38 
39  addToRunTimeSelectionTable(topoChangerFvMesh, mixerFvMesh, IOobject);
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 void Foam::mixerFvMesh::addZonesAndModifiers()
46 {
47  // Add zones and modifiers for motion action
48 
49  if
50  (
51  pointZones().size()
52  || faceZones().size()
53  || cellZones().size()
54  || topoChanger_.size()
55  )
56  {
58  << "Zones and modifiers already present. Skipping."
59  << endl;
60 
61  return;
62  }
63 
64  Info<< "Time = " << time().timeName() << endl
65  << "Adding zones and modifiers to the mesh" << endl;
66 
67  // Add zones
68  List<pointZone*> pz(1);
69 
70  // Add an empty zone for cut points
71 
72  pz[0] = new pointZone
73  (
74  "cutPointZone",
75  labelList(0),
76  0,
77  pointZones()
78  );
79 
80 
81  // Do face zones for slider
82 
83  List<faceZone*> fz(3);
84 
85  // Inner slider
86  const word innerSliderName(motionDict_.subDict("slider").lookup("inside"));
87  const polyPatch& innerSlider = boundaryMesh()[innerSliderName];
88 
89  labelList isf(innerSlider.size());
90 
91  forAll(isf, i)
92  {
93  isf[i] = innerSlider.start() + i;
94  }
95 
96  fz[0] = new faceZone
97  (
98  "insideSliderZone",
99  isf,
100  boolList(innerSlider.size(), false),
101  0,
102  faceZones()
103  );
104 
105  // Outer slider
106  const word outerSliderName(motionDict_.subDict("slider").lookup("outside"));
107  const polyPatch& outerSlider = boundaryMesh()[outerSliderName];
108 
109  labelList osf(outerSlider.size());
110 
111  forAll(osf, i)
112  {
113  osf[i] = outerSlider.start() + i;
114  }
115 
116  fz[1] = new faceZone
117  (
118  "outsideSliderZone",
119  osf,
120  boolList(outerSlider.size(), false),
121  1,
122  faceZones()
123  );
124 
125  // Add empty zone for cut faces
126  fz[2] = new faceZone
127  (
128  "cutFaceZone",
129  labelList(0),
130  boolList(0, false),
131  2,
132  faceZones()
133  );
134 
135  List<cellZone*> cz(1);
136 
137  // Mark every cell with its topological region
138  regionSplit rs(*this);
139 
140  // Get the region of the cell containing the origin.
141  label originRegion = rs[findNearestCell(cs().origin())];
142 
143  labelList movingCells(nCells());
144  label nMovingCells = 0;
145 
146  forAll(rs, celli)
147  {
148  if (rs[celli] == originRegion)
149  {
150  movingCells[nMovingCells] = celli;
151  nMovingCells++;
152  }
153  }
154 
155  movingCells.setSize(nMovingCells);
156  Info<< "Number of cells in the moving region: " << nMovingCells << endl;
157 
158  cz[0] = new cellZone
159  (
160  "movingCells",
161  movingCells,
162  0,
163  cellZones()
164  );
165 
166  Info<< "Adding point, face and cell zones" << endl;
167  addZones(pz, fz, cz);
168 
169  // Add a topology modifier
170  Info<< "Adding topology modifiers" << endl;
171  topoChanger_.setSize(1);
172  topoChanger_.set
173  (
174  0,
175  new slidingInterface
176  (
177  "mixerSlider",
178  0,
179  topoChanger_,
180  outerSliderName + "Zone",
181  innerSliderName + "Zone",
182  "cutPointZone",
183  "cutFaceZone",
184  outerSliderName,
185  innerSliderName,
187  )
188  );
189  topoChanger_.writeOpt() = IOobject::AUTO_WRITE;
190 
191  write();
192 }
193 
194 
195 void Foam::mixerFvMesh::calcMovingMasks() const
196 {
197  if (debug)
198  {
200  << "Calculating point and cell masks"
201  << endl;
202  }
203 
204  if (movingPointsMaskPtr_)
205  {
207  << "point mask already calculated"
208  << abort(FatalError);
209  }
210 
211  // Set the point mask
212  movingPointsMaskPtr_ = new scalarField(points().size(), 0);
213  scalarField& movingPointsMask = *movingPointsMaskPtr_;
214 
215  const cellList& c = cells();
216  const faceList& f = faces();
217 
218  const labelList& cellAddr = cellZones()["movingCells"];
219 
220  forAll(cellAddr, celli)
221  {
222  const cell& curCell = c[cellAddr[celli]];
223 
224  forAll(curCell, facei)
225  {
226  // Mark all the points as moving
227  const face& curFace = f[curCell[facei]];
228 
229  forAll(curFace, pointi)
230  {
231  movingPointsMask[curFace[pointi]] = 1;
232  }
233  }
234  }
235 
236  const word innerSliderZoneName
237  (
238  word(motionDict_.subDict("slider").lookup("inside"))
239  + "Zone"
240  );
241 
242  const labelList& innerSliderAddr = faceZones()[innerSliderZoneName];
243 
244  forAll(innerSliderAddr, facei)
245  {
246  const face& curFace = f[innerSliderAddr[facei]];
247 
248  forAll(curFace, pointi)
249  {
250  movingPointsMask[curFace[pointi]] = 1;
251  }
252  }
253 
254  const word outerSliderZoneName
255  (
256  word(motionDict_.subDict("slider").lookup("outside"))
257  + "Zone"
258  );
259 
260  const labelList& outerSliderAddr = faceZones()[outerSliderZoneName];
261 
262  forAll(outerSliderAddr, facei)
263  {
264  const face& curFace = f[outerSliderAddr[facei]];
265 
266  forAll(curFace, pointi)
267  {
268  movingPointsMask[curFace[pointi]] = 0;
269  }
270  }
271 }
272 
273 
274 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
275 
276 // Construct from components
277 Foam::mixerFvMesh::mixerFvMesh
278 (
279  const IOobject& io
280 )
281 :
282  topoChangerFvMesh(io),
283  motionDict_
284  (
286  (
287  IOobject
288  (
289  "dynamicMeshDict",
290  time().constant(),
291  *this,
294  false
295  )
296  ).optionalSubDict(typeName + "Coeffs")
297  ),
298  csPtr_
299  (
301  (
302  "coordinateSystem",
303  motionDict_.subDict("coordinateSystem")
304  )
305  ),
306  rpm_(readScalar(motionDict_.lookup("rpm"))),
307  movingPointsMaskPtr_(nullptr)
308 {
309  addZonesAndModifiers();
310 
311  Info<< "Mixer mesh:" << nl
312  << " origin: " << cs().origin() << nl
313  << " axis: " << cs().axis() << nl
314  << " rpm: " << rpm_ << endl;
315 }
316 
317 
318 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
319 
321 {
322  deleteDemandDrivenData(movingPointsMaskPtr_);
323 }
324 
325 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
326 
327 // Return moving points mask. Moving points marked with 1
328 const Foam::scalarField& Foam::mixerFvMesh::movingPointsMask() const
329 {
330  if (!movingPointsMaskPtr_)
331  {
332  calcMovingMasks();
333  }
334 
335  return *movingPointsMaskPtr_;
336 }
337 
338 
340 {
341  // Rotational speed needs to be converted from rpm
342  movePoints
343  (
344  csPtr_->globalPosition
345  (
346  csPtr_->localPosition(points())
347  + vector(0, rpm_*360.0*time().deltaTValue()/60.0, 0)
348  *movingPointsMask()
349  )
350  );
351 
352  // Make changes. Use inflation (so put new points in topoChangeMap)
353  autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh(true);
354 
355  if (topoChangeMap.valid())
356  {
357  if (debug)
358  {
359  InfoInFunction << "Mesh topology is changing" << endl;
360  }
361 
362  deleteDemandDrivenData(movingPointsMaskPtr_);
363  }
364 
365  movePoints
366  (
367  csPtr_->globalPosition
368  (
369  csPtr_->localPosition(oldPoints())
370  + vector(0, rpm_*360.0*time().deltaTValue()/60.0, 0)
371  *movingPointsMask()
372  )
373  );
374 
375  return true;
376 }
377 
378 
379 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Macros for easy insertion into run-time selection tables.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
const cellShapeList & cells
const pointField & points
virtual ~mixerFvMesh()
Destructor.
Definition: mixerFvMesh.C:320
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Abstract base class for a topology changing fvMesh.
List< label > labelList
A List of labels.
Definition: labelList.H:56
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
errorManip< error > abort(error &err)
Definition: errorManip.H:131
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:265
defineTypeNameAndDebug(combustionModel, 0)
static autoPtr< coordinateSystem > New(const objectRegistry &obr, const dictionary &dict)
Select constructed from dictionary and objectRegistry.
messageStream Info
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual bool update()
Update the mesh for both mesh motion and topology change.
Definition: mixerFvMesh.C:339
void deleteDemandDrivenData(DataPtr &dataPtr)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
List< cell > cellList
list of cells
Definition: cellList.H:42
Namespace for OpenFOAM.
#define InfoInFunction
Report an information message using Foam::Info.