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 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 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(mixerFvMesh, 0);
42  addToRunTimeSelectionTable(topoChangerFvMesh, mixerFvMesh, IOobject);
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::mixerFvMesh::addZonesAndModifiers()
49 {
50  // Add zones and modifiers for motion action
51 
52  if
53  (
54  pointZones().size()
55  || faceZones().size()
56  || cellZones().size()
57  || topoChanger_.size()
58  )
59  {
61  << "Zones and modifiers already present. Skipping."
62  << endl;
63 
64  return;
65  }
66 
67  Info<< "Time = " << time().timeName() << endl
68  << "Adding zones and modifiers to the mesh" << endl;
69 
70  // Add zones
71  List<pointZone*> pz(1);
72 
73  // An empty zone for cut points
74  pz[0] = new pointZone("cutPointZone", 0, pointZones());
75 
76 
77  // Do face zones for slider
78 
79  List<faceZone*> fz(3);
80 
81  // Inner slider
82  const word innerSliderName
83  (
84  motionDict_.subDict("slider").get<word>("inside")
85  );
86  const polyPatch& innerSlider = boundaryMesh()[innerSliderName];
87 
88  fz[0] = new faceZone
89  (
90  "insideSliderZone",
91  identity(innerSlider.size(), innerSlider.start()),
92  false, // none are flipped
93  0,
94  faceZones()
95  );
96 
97  // Outer slider
98  const word outerSliderName
99  (
100  motionDict_.subDict("slider").get<word>("outside")
101  );
102  const polyPatch& outerSlider = boundaryMesh()[outerSliderName];
103 
104  fz[1] = new faceZone
105  (
106  "outsideSliderZone",
107  identity(outerSlider.size(), outerSlider.start()),
108  false, // none are flipped
109  1,
110  faceZones()
111  );
112 
113  // An empty zone for cut faces
114  fz[2] = new faceZone("cutFaceZone", 2, faceZones());
115 
116  List<cellZone*> cz(1);
117 
118  // Mark every cell with its topological region
119  regionSplit rs(*this);
120 
121  // Get the region of the cell containing the origin.
122  const label originRegion = rs[findNearestCell(csys_.origin())];
123 
124  labelList movingCells(nCells());
125  label nMovingCells = 0;
126 
127  forAll(rs, celli)
128  {
129  if (rs[celli] == originRegion)
130  {
131  movingCells[nMovingCells] = celli;
132  ++nMovingCells;
133  }
134  }
135 
136  movingCells.resize(nMovingCells);
137  Info<< "Number of cells in the moving region: " << nMovingCells << endl;
138 
139  cz[0] = new cellZone
140  (
141  "movingCells",
142  std::move(movingCells),
143  0,
144  cellZones()
145  );
146 
147  Info<< "Adding point, face and cell zones" << endl;
148  addZones(pz, fz, cz);
149 
150  // Add a topology modifier
151  Info<< "Adding topology modifiers" << endl;
154  (
155  0,
156  new slidingInterface
157  (
158  "mixerSlider",
159  0,
160  topoChanger_,
161  outerSliderName + "Zone",
162  innerSliderName + "Zone",
163  "cutPointZone",
164  "cutFaceZone",
165  outerSliderName,
166  innerSliderName,
168  )
169  );
171 
172  write();
173 }
174 
175 
176 void Foam::mixerFvMesh::calcMovingMasks() const
177 {
178  if (debug)
179  {
181  << "Calculating point and cell masks"
182  << endl;
183  }
184 
185  if (movingPointsMaskPtr_)
186  {
188  << "point mask already calculated"
189  << abort(FatalError);
190  }
191 
192  // Set the point mask
193  movingPointsMaskPtr_ = new scalarField(points().size(), Zero);
194  scalarField& movingPointsMask = *movingPointsMaskPtr_;
195 
196  const cellList& c = cells();
197  const faceList& f = faces();
198 
199  const labelList& cellAddr = cellZones()["movingCells"];
200 
201  for (const label celli : cellAddr)
202  {
203  const cell& curCell = c[celli];
204 
205  for (const label facei : curCell)
206  {
207  // Mark all the points as moving
208  const face& curFace = f[facei];
209 
210  forAll(curFace, pointi)
211  {
212  movingPointsMask[curFace[pointi]] = 1;
213  }
214  }
215  }
216 
217  const word innerSliderZoneName
218  (
219  motionDict_.subDict("slider").get<word>("inside") + "Zone"
220  );
221 
222  const labelList& innerSliderAddr = faceZones()[innerSliderZoneName];
223 
224  for (const label facei : innerSliderAddr)
225  {
226  const face& curFace = f[facei];
227 
228  forAll(curFace, pointi)
229  {
230  movingPointsMask[curFace[pointi]] = 1;
231  }
232  }
233 
234  const word outerSliderZoneName
235  (
236  motionDict_.subDict("slider").get<word>("outside") + "Zone"
237  );
238 
239  const labelList& outerSliderAddr = faceZones()[outerSliderZoneName];
240 
241  for (const label facei : outerSliderAddr)
242  {
243  const face& curFace = f[facei];
244 
245  forAll(curFace, pointi)
246  {
247  movingPointsMask[curFace[pointi]] = 0;
248  }
249  }
250 }
251 
252 
253 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
254 
255 Foam::mixerFvMesh::mixerFvMesh
256 (
257  const IOobject& io
258 )
259 :
260  topoChangerFvMesh(io),
261  motionDict_
262  (
264  (
265  IOobject
266  (
267  "dynamicMeshDict",
268  time().constant(),
269  *this,
272  false
273  )
274  ).optionalSubDict(typeName + "Coeffs")
275  ),
276  csys_(),
277  rpm_(motionDict_.get<scalar>("rpm")),
278  movingPointsMaskPtr_(nullptr)
279 {
280  if (motionDict_.found(coordinateSystem::typeName_()))
281  {
282  // New() for access to indirect (global) coordSystem.
283  static_cast<coordinateSystem&>(csys_) =
284  *coordinateSystem::New(*this, motionDict_);
285  }
286  else
287  {
288  csys_ = coordSystem::cylindrical(motionDict_);
289  }
290 
291  addZonesAndModifiers();
292 
293  Info<< "Mixer mesh:" << nl
294  << " origin: " << csys_.origin() << nl
295  << " axis: " << csys_.e3() << nl
296  << " rpm: " << rpm_ << endl;
297 }
298 
299 
300 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
301 
303 {
304  deleteDemandDrivenData(movingPointsMaskPtr_);
305 }
306 
307 
308 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
309 
310 // Return moving points mask. Moving points marked with 1
311 const Foam::scalarField& Foam::mixerFvMesh::movingPointsMask() const
312 {
313  if (!movingPointsMaskPtr_)
314  {
315  calcMovingMasks();
316  }
317 
318  return *movingPointsMaskPtr_;
319 }
320 
321 
323 {
324  // The tangential sweep (radians)
325  const vector theta(0, rpmToRads(rpm_)*time().deltaTValue(), 0);
326 
327  movePoints
328  (
329  csys_.globalPosition
330  (
331  csys_.localPosition(points())
332  + theta
333  *movingPointsMask()
334  )
335  );
336 
337  // Make changes. Use inflation (so put new points in topoChangeMap)
338  autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh(true);
339 
340  if (topoChangeMap.valid())
341  {
342  if (debug)
343  {
344  InfoInFunction << "Mesh topology is changing" << endl;
345  }
346 
347  deleteDemandDrivenData(movingPointsMaskPtr_);
348  }
349 
350  movePoints
351  (
352  csys_.globalPosition
353  (
354  csys_.localPosition(oldPoints())
355  + theta
356  *movingPointsMask()
357  )
358  );
359 
360  return true;
361 }
362 
363 
364 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::HashTable< regIOobject * >::size
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:316
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::fvMesh::write
virtual bool write(const bool valid=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:939
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
mapPolyMesh.H
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:764
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:483
unitConversion.H
Unit conversion functions.
Foam::autoPtr::valid
bool valid() const noexcept
True if the managed pointer is non-null.
Definition: autoPtrI.H:107
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:477
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::coordinateSystem::origin
virtual const point & origin() const
Return origin.
Definition: coordinateSystem.H:459
Foam::primitiveMesh::findNearestCell
label findNearestCell(const point &location) const
Find the cell with the nearest cell centre to location.
Definition: primitiveMeshFindCell.C:88
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::label
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:62
Foam::Field< scalar >
Foam::IOobject::writeOpt
writeOption writeOpt() const
The write option.
Definition: IOobjectI.H:153
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::polyMesh::pointZones
const pointZoneMesh & pointZones() const
Return point zone mesh.
Definition: polyMesh.H:471
Foam::coordinateSystem::New
static autoPtr< coordinateSystem > New(word modelType, const objectRegistry &obr, const dictionary &dict)
Definition: coordinateSystemNew.C:82
mixerFvMesh.H
regionSplit.H
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:523
Foam::cellList
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:47
Foam::PtrList::setSize
void setSize(const label newLen)
Same as resize()
Definition: PtrListI.H:108
Foam::coordSystem::cylindrical
A cylindrical coordinate system (r-theta-z). The coordinate system angle theta is always in radians.
Definition: cylindricalCS.H:71
Foam::FatalError
error FatalError
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::PtrList::set
const T * set(const label i) const
Return const pointer to element (if set) or nullptr.
Definition: PtrListI.H:143
Foam::mixerFvMesh::update
virtual bool update()
Update the mesh for both mesh motion and topology change.
Definition: mixerFvMesh.C:322
Foam::rpmToRads
constexpr scalar rpmToRads(const scalar rpm) noexcept
Conversion from revolutions/minute to radians/sec.
Definition: unitConversion.H:73
Time.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
slidingInterface.H
Foam::nl
constexpr char nl
Definition: Ostream.H:372
f
labelList f(nPoints)
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::Vector< scalar >
Foam::topoChangerFvMesh
Abstract base class for a topology changing fvMesh.
Definition: topoChangerFvMesh.H:52
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::IOobject::MUST_READ_IF_MODIFIED
Definition: IOobject.H:121
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::mixerFvMesh::~mixerFvMesh
virtual ~mixerFvMesh()
Destructor.
Definition: mixerFvMesh.C:302
constant
constant condensation/saturation model.
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::slidingInterface::INTEGRAL
Definition: slidingInterface.H:84
Foam::polyMesh::addZones
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
Definition: polyMesh.C:968
Foam::coordinateSystem
Base class for coordinate system specification, the default coordinate system type is cartesian .
Definition: coordinateSystem.H:122
Foam::topoChangerFvMesh::topoChanger_
polyTopoChanger topoChanger_
Definition: topoChangerFvMesh.H:69