mixerFvMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2018-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "mixerFvMesh.H"
30 #include "Time.H"
31 #include "regionSplit.H"
32 #include "slidingInterface.H"
34 #include "mapPolyMesh.H"
35 #include "unitConversion.H"
36 #include "demandDrivenData.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(mixerFvMesh, 0);
43  addToRunTimeSelectionTable(topoChangerFvMesh, mixerFvMesh, IOobject);
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::mixerFvMesh::addZonesAndModifiers()
50 {
51  // Add zones and modifiers for motion action
52 
53  if
54  (
55  pointZones().size()
56  || faceZones().size()
57  || cellZones().size()
58  || topoChanger_.size()
59  )
60  {
62  << "Zones and modifiers already present. Skipping."
63  << endl;
64 
65  return;
66  }
67 
68  Info<< "Time = " << time().timeName() << endl
69  << "Adding zones and modifiers to the mesh" << endl;
70 
71  // Add zones
72  List<pointZone*> pz(1);
73 
74  // An empty zone for cut points
75  pz[0] = new pointZone("cutPointZone", 0, pointZones());
76 
77 
78  // Do face zones for slider
79 
80  List<faceZone*> fz(3);
81 
82  // Inner slider
83  const word innerSliderName
84  (
85  motionDict_.subDict("slider").get<word>("inside")
86  );
87  const polyPatch& innerSlider = boundaryMesh()[innerSliderName];
88 
89  fz[0] = new faceZone
90  (
91  "insideSliderZone",
92  identity(innerSlider.range()),
93  false, // none are flipped
94  0,
95  faceZones()
96  );
97 
98  // Outer slider
99  const word outerSliderName
100  (
101  motionDict_.subDict("slider").get<word>("outside")
102  );
103  const polyPatch& outerSlider = boundaryMesh()[outerSliderName];
104 
105  fz[1] = new faceZone
106  (
107  "outsideSliderZone",
108  identity(outerSlider.range()),
109  false, // none are flipped
110  1,
111  faceZones()
112  );
113 
114  // An empty zone for cut faces
115  fz[2] = new faceZone("cutFaceZone", 2, faceZones());
116 
117  List<cellZone*> cz(1);
118 
119  // Mark every cell with its topological region
120  regionSplit rs(*this);
121 
122  // Get the region of the cell containing the origin.
123  const label originRegion = rs[findNearestCell(csys_.origin())];
124 
125  labelList movingCells(nCells());
126  label nMovingCells = 0;
127 
128  forAll(rs, celli)
129  {
130  if (rs[celli] == originRegion)
131  {
132  movingCells[nMovingCells] = celli;
133  ++nMovingCells;
134  }
135  }
136 
137  movingCells.resize(nMovingCells);
138  Info<< "Number of cells in the moving region: " << nMovingCells << endl;
139 
140  cz[0] = new cellZone
141  (
142  "movingCells",
143  std::move(movingCells),
144  0,
145  cellZones()
146  );
147 
148  Info<< "Adding point, face and cell zones" << endl;
149  addZones(pz, fz, cz);
150 
151  // Add a topology modifier
152  Info<< "Adding topology modifiers" << endl;
155  (
156  0,
157  new slidingInterface
158  (
159  "mixerSlider",
160  0,
161  topoChanger_,
162  outerSliderName + "Zone",
163  innerSliderName + "Zone",
164  "cutPointZone",
165  "cutFaceZone",
166  outerSliderName,
167  innerSliderName,
169  )
170  );
172 
173  write();
174 }
175 
176 
177 void Foam::mixerFvMesh::calcMovingMasks() const
178 {
179  DebugInFunction << "Calculating point and cell masks" << endl;
180 
181  if (movingPointsMaskPtr_)
182  {
184  << "point mask already calculated"
185  << abort(FatalError);
186  }
187 
188  // Set the point mask
189  movingPointsMaskPtr_ = new scalarField(points().size(), Zero);
190  scalarField& movingPointsMask = *movingPointsMaskPtr_;
191 
192  const cellList& c = cells();
193  const faceList& f = faces();
194 
195  const labelList& cellAddr = cellZones()["movingCells"];
196 
197  for (const label celli : cellAddr)
198  {
199  const cell& curCell = c[celli];
200 
201  for (const label facei : curCell)
202  {
203  // Mark all the points as moving
204  const face& curFace = f[facei];
205 
206  forAll(curFace, pointi)
207  {
208  movingPointsMask[curFace[pointi]] = 1;
209  }
210  }
211  }
212 
213  const word innerSliderZoneName
214  (
215  motionDict_.subDict("slider").get<word>("inside") + "Zone"
216  );
217 
218  const labelList& innerSliderAddr = faceZones()[innerSliderZoneName];
219 
220  for (const label facei : innerSliderAddr)
221  {
222  const face& curFace = f[facei];
223 
224  forAll(curFace, pointi)
225  {
226  movingPointsMask[curFace[pointi]] = 1;
227  }
228  }
229 
230  const word outerSliderZoneName
231  (
232  motionDict_.subDict("slider").get<word>("outside") + "Zone"
233  );
234 
235  const labelList& outerSliderAddr = faceZones()[outerSliderZoneName];
236 
237  for (const label facei : outerSliderAddr)
238  {
239  const face& curFace = f[facei];
240 
241  forAll(curFace, pointi)
242  {
243  movingPointsMask[curFace[pointi]] = 0;
244  }
245  }
246 }
247 
248 
249 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
250 
251 Foam::mixerFvMesh::mixerFvMesh
252 (
253  const IOobject& io
254 )
255 :
256  topoChangerFvMesh(io),
257  motionDict_
258  (
259  IOdictionary::readContents
260  (
261  IOobject
262  (
263  "dynamicMeshDict",
264  time().constant(),
265  *this,
266  IOobject::MUST_READ
267  )
268  ).optionalSubDict(typeName + "Coeffs")
269  ),
270  csys_(),
271  rpm_(motionDict_.get<scalar>("rpm")),
272  movingPointsMaskPtr_(nullptr)
273 {
274  // New() for access to indirect (global) coordSystem.
275 
276  auto csysPtr = coordinateSystem::NewIfPresent(*this, dict);
277 
278  if (csysPtr)
279  {
280  static_cast<coordinateSystem&>(csys_) = csysPtr();
281  }
282  else
283  {
284  csys_ = coordSystem::cylindrical(motionDict_);
285  }
286 
287  addZonesAndModifiers();
288 
289  Info<< "Mixer mesh:" << nl
290  << " origin: " << csys_.origin() << nl
291  << " axis: " << csys_.e3() << nl
292  << " rpm: " << rpm_ << endl;
293 }
294 
295 
296 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
297 
299 {
300  deleteDemandDrivenData(movingPointsMaskPtr_);
301 }
302 
303 
304 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
305 
306 // Return moving points mask. Moving points marked with 1
307 const Foam::scalarField& Foam::mixerFvMesh::movingPointsMask() const
308 {
309  if (!movingPointsMaskPtr_)
310  {
311  calcMovingMasks();
312  }
313 
314  return *movingPointsMaskPtr_;
315 }
316 
317 
319 {
320  // The tangential sweep (radians)
321  const vector theta(0, rpmToRads(rpm_)*time().deltaTValue(), 0);
322 
323  movePoints
324  (
325  csys_.globalPosition
326  (
327  csys_.localPosition(points())
328  + theta
329  *movingPointsMask()
330  )
331  );
332 
333  // Make changes. Use inflation (so put new points in topoChangeMap)
334  autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh(true);
335 
336  if (topoChangeMap)
337  {
338  DebugInFunction << "Mesh topology is changing" << nl;
339 
340  deleteDemandDrivenData(movingPointsMaskPtr_);
341  }
342 
343  movePoints
344  (
345  csys_.globalPosition
346  (
347  csys_.localPosition(oldPoints())
348  + theta
349  *movingPointsMask()
350  )
351  );
352 
353  return true;
354 }
355 
356 
357 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
List< cell > cellList
List of cell.
Definition: cellListFwd.H:39
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
constexpr scalar rpmToRads(const scalar rpm) noexcept
Conversion from revolutions/minute to radians/sec.
Unit conversion functions.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
virtual const vector e3() const
The local Cartesian z-axis in global coordinates.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
Macros for easy insertion into run-time selection tables.
virtual const point & origin() const
Return origin.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
label size() const noexcept
The number of elements in table.
Definition: HashTable.H:342
writeOption writeOpt() const noexcept
Get the write option.
const cellShapeList & cells
const pointField & points
virtual ~mixerFvMesh()
Destructor.
Definition: mixerFvMesh.C:291
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
#define DebugInFunction
Report an information message using Foam::Info.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
void addZones(PtrList< pointZone > &&pz, PtrList< faceZone > &&fz, PtrList< cellZone > &&cz)
Add mesh zones.
Definition: polyMesh.C:1009
Vector< scalar > vector
Definition: vector.H:57
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1102
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:670
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
void setSize(const label newLen)
Same as resize()
Definition: PtrList.H:337
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition: polyMesh.H:662
Template functions to aid in the implementation of demand driven data.
label nCells() const noexcept
Number of mesh cells.
const dimensionedScalar c
Speed of light in a vacuum.
label findNearestCell(const point &location) const
Find the cell with the nearest cell centre to location.
Automatically write from objectRegistry::writeObject()
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:678
messageStream Info
Information stream (stdout output on master, null elsewhere)
List< label > labelList
A List of labels.
Definition: List.H:62
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
polyTopoChanger topoChanger_
virtual bool update()
Update the mesh for both mesh motion and topology change.
Definition: mixerFvMesh.C:311
void deleteDemandDrivenData(DataPtr &dataPtr)
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
#define InfoInFunction
Report an information message using Foam::Info.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127