Loading...
Searching...
No Matches
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-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
43}
44
45
46// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47
48void Foam::mixerFvMesh::addZonesAndModifiers()
49{
50 // Add zones and modifiers for motion action
51
52 if
53 (
55 || faceZones().size()
56 || cellZones().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.range()),
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.range()),
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;
152 topoChanger_.setSize(1);
153 topoChanger_.set
154 (
155 0,
156 new slidingInterface
157 (
158 "mixerSlider",
159 0,
161 outerSliderName + "Zone",
162 innerSliderName + "Zone",
163 "cutPointZone",
164 "cutFaceZone",
165 outerSliderName,
166 innerSliderName,
168 )
169 );
171
172 write();
173}
174
175
176void Foam::mixerFvMesh::calcMovingMasks() const
177{
178 DebugInFunction << "Calculating point and cell masks" << endl;
179
180 if (movingPointsMaskPtr_)
181 {
183 << "point mask already calculated"
184 << abort(FatalError);
185 }
186
187 // Set the point mask
188 movingPointsMaskPtr_ = std::make_unique<scalarField>(points().size(), Zero);
189 auto& movingPointsMask = *movingPointsMaskPtr_;
190
191 const cellList& c = cells();
192 const faceList& f = faces();
193
194 const labelList& cellAddr = cellZones()["movingCells"];
195
196 for (const label celli : cellAddr)
197 {
198 const cell& curCell = c[celli];
199
200 for (const label facei : curCell)
201 {
202 // Mark all the points as moving
203 const face& curFace = f[facei];
204
205 forAll(curFace, pointi)
206 {
207 movingPointsMask[curFace[pointi]] = 1;
208 }
209 }
210 }
211
212 const word innerSliderZoneName
213 (
214 motionDict_.subDict("slider").get<word>("inside") + "Zone"
215 );
216
217 const labelList& innerSliderAddr = faceZones()[innerSliderZoneName];
218
219 for (const label facei : innerSliderAddr)
220 {
221 const face& curFace = f[facei];
222
223 forAll(curFace, pointi)
224 {
225 movingPointsMask[curFace[pointi]] = 1;
226 }
227 }
228
229 const word outerSliderZoneName
230 (
231 motionDict_.subDict("slider").get<word>("outside") + "Zone"
232 );
233
234 const labelList& outerSliderAddr = faceZones()[outerSliderZoneName];
235
236 for (const label facei : outerSliderAddr)
237 {
238 const face& curFace = f[facei];
239
240 forAll(curFace, pointi)
241 {
242 movingPointsMask[curFace[pointi]] = 0;
244 }
245}
246
247
248// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
249
250Foam::mixerFvMesh::mixerFvMesh
251(
252 const IOobject& io
253)
254:
256 motionDict_
257 (
258 IOdictionary::readContents
259 (
261 (
262 "dynamicMeshDict",
263 time().constant(),
264 *this,
265 IOobject::MUST_READ
266 )
267 ).optionalSubDict(typeName + "Coeffs")
268 ),
269 rpm_(motionDict_.get<scalar>("rpm"))
270{
271 // New() for access to indirect (global) coordSystem.
272
273 auto csysPtr = coordinateSystem::NewIfPresent(*this, dict);
274
275 if (csysPtr)
276 {
277 static_cast<coordinateSystem&>(csys_) = *csysPtr;
278 }
279 else
280 {
281 csys_ = coordSystem::cylindrical(motionDict_);
282 }
283
284 addZonesAndModifiers();
285
286 Info<< "Mixer mesh:" << nl
287 << " origin: " << csys_.origin() << nl
288 << " axis: " << csys_.e3() << nl
289 << " rpm: " << rpm_ << endl;
290}
291
292
293// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
294
296{
297 movingPointsMaskPtr_.reset(nullptr);
298}
299
300
301// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
302
303// Return moving points mask. Moving points marked with 1
304const Foam::scalarField& Foam::mixerFvMesh::movingPointsMask() const
305{
306 if (!movingPointsMaskPtr_)
307 {
308 calcMovingMasks();
309 }
310
311 return *movingPointsMaskPtr_;
312}
313
314
316{
317 // The tangential sweep (radians)
318 const vector theta(0, rpmToRads(rpm_)*time().deltaTValue(), 0);
319
320 movePoints
321 (
322 csys_.globalPosition
323 (
324 csys_.localPosition(points())
325 + theta
326 *movingPointsMask()
327 )
328 );
329
330 // Make changes. Use inflation (so put new points in topoChangeMap)
331 autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh(true);
332
333 if (topoChangeMap)
334 {
335 DebugInFunction << "Mesh topology is changing" << nl;
336
337 movingPointsMaskPtr_.reset(nullptr);
338 }
339
340 movePoints
341 (
342 csys_.globalPosition
343 (
344 csys_.localPosition(oldPoints())
345 + theta
346 *movingPointsMask()
347 )
348 );
349
350 return true;
351}
352
353
354// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
label size() const noexcept
Definition HashTable.H:358
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
static dictionary readContents(const IOobject &io)
Read and return contents, testing for "dictionary" type. The IOobject will not be registered.
@ MUST_READ
Reading required.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
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
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
A cylindrical coordinate system (r-theta-z). The coordinate system angle theta is always in radians.
Base class for coordinate system specification, the default coordinate system type is cartesian .
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 dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition dictionary.C:560
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition fvMesh.C:884
A rotating slider mesh.
Definition mixerFvMesh.H:50
virtual bool update()
Update the mesh for both mesh motion and topology change.
virtual ~mixerFvMesh()
Destructor.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition polyMesh.H:671
virtual const pointField & oldPoints() const
Return old points (mesh motion).
Definition polyMesh.C:1113
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition polyMesh.H:679
void addZones(PtrList< pointZone > &&pz, PtrList< faceZone > &&fz, PtrList< cellZone > &&cz)
Add mesh zones.
Definition polyMesh.C:994
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition polyMesh.H:663
label findNearestCell(const point &location) const
Find the cell with the nearest cell centre to location.
label nCells() const noexcept
Number of mesh cells.
Abstract base class for a topology changing fvMesh.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
const pointField & points
const cellShapeList & cells
#define DebugInFunction
Report an information message using Foam::Info.
#define InfoInFunction
Report an information message using Foam::Info.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
const dimensionedScalar c
Speed of light in a vacuum.
Different types of constants.
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr scalar rpmToRads() noexcept
Multiplication factor for revolutions/minute to radians/sec.
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
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...
errorManip< error > abort(error &err)
Definition errorManip.H:139
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
runTime write()
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.