Loading...
Searching...
No Matches
movingConeTopoFvMesh.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-2020 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
30#include "Time.H"
31#include "mapPolyMesh.H"
34#include "meshTools.H"
35#include "OFstream.H"
38using namespace Foam::constant::mathematical;
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
42namespace Foam
43{
47 (
51 );
53 (
56 doInit
57 );
58}
59
60
61// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
62
63Foam::tmp<Foam::scalarField> Foam::movingConeTopoFvMesh::vertexMarkup
64(
65 const pointField& p,
66 const scalar curLeft,
67 const scalar curRight
68) const
69{
70 Info<< "Updating vertex markup. curLeft: "
71 << curLeft << " curRight: " << curRight << endl;
72
73 auto tvertexMarkup = tmp<scalarField>::New(p.size());
74 auto& vertexMarkup = tvertexMarkup.ref();
75
76 forAll(p, pI)
77 {
78 if (p[pI].x() < curLeft - SMALL)
79 {
80 vertexMarkup[pI] = -1;
81 }
82 else if (p[pI].x() > curRight + SMALL)
83 {
84 vertexMarkup[pI] = 1;
85 }
86 else
87 {
88 vertexMarkup[pI] = 0;
89 }
90 }
91
92 return tvertexMarkup;
93}
94
95
96void Foam::movingConeTopoFvMesh::addZonesAndModifiers()
97{
98 // Add zones and modifiers for motion action
99
100 if
101 (
102 pointZones().size()
103 || faceZones().size()
104 || cellZones().size()
105 || topoChanger_.size()
106 )
107 {
109 << "Zones and modifiers already present. Skipping."
110 << endl;
111
112 return;
113 }
114
115 Info<< "Time = " << time().timeName() << endl
116 << "Adding zones and modifiers to the mesh" << endl;
117
118 const vectorField& fc = faceCentres();
119 const vectorField& fa = faceAreas();
120
121 labelList zone1(fc.size());
122 boolList flipZone1(fc.size(), false);
123 label nZoneFaces1 = 0;
124
125 labelList zone2(fc.size());
126 boolList flipZone2(fc.size(), false);
127 label nZoneFaces2 = 0;
128
129 forAll(fc, facei)
130 {
131 if
132 (
133 fc[facei].x() > -0.003501
134 && fc[facei].x() < -0.003499
135 )
136 {
137 if (fa[facei].x() < 0)
138 {
139 flipZone1[nZoneFaces1] = true;
140 }
141
142 zone1[nZoneFaces1] = facei;
143 Info<< "face " << facei << " for zone 1. Flip: "
144 << flipZone1[nZoneFaces1] << endl;
145 ++nZoneFaces1;
146 }
147 else if
148 (
149 fc[facei].x() > -0.00701
150 && fc[facei].x() < -0.00699
151 )
152 {
153 zone2[nZoneFaces2] = facei;
154
155 if (fa[facei].x() > 0)
156 {
157 flipZone2[nZoneFaces2] = true;
158 }
159
160 Info<< "face " << facei << " for zone 2. Flip: "
161 << flipZone2[nZoneFaces2] << endl;
162 ++nZoneFaces2;
163 }
164 }
165
166 zone1.setSize(nZoneFaces1);
167 flipZone1.setSize(nZoneFaces1);
168
169 zone2.setSize(nZoneFaces2);
170 flipZone2.setSize(nZoneFaces2);
171
172 Info<< "zone: " << zone1 << endl;
173 Info<< "zone: " << zone2 << endl;
174
175 List<pointZone*> pz(0);
176 List<faceZone*> fz(2);
177 List<cellZone*> cz(0);
178
179 label nFz = 0;
180
181 fz[nFz] =
182 new faceZone
183 (
184 "rightExtrusionFaces",
185 std::move(zone1),
186 std::move(flipZone1),
187 nFz,
188 faceZones()
189 );
190 ++nFz;
191
192 fz[nFz] =
193 new faceZone
194 (
195 "leftExtrusionFaces",
196 std::move(zone2),
197 std::move(flipZone2),
198 nFz,
199 faceZones()
200 );
201 ++nFz;
202
203 fz.setSize(nFz);
204
205 Info<< "Adding mesh zones." << endl;
206 addZones(pz, fz, cz);
207
208
209 // Add layer addition/removal interfaces
210
212 label nMods = 0;
213
214 tm[nMods] =
216 (
217 "right",
218 nMods,
219 topoChanger_,
220 "rightExtrusionFaces",
221 motionDict_.subDict("right").get<scalar>("minThickness"),
222 motionDict_.subDict("right").get<scalar>("maxThickness")
223 );
224 ++nMods;
225
226 tm[nMods] = new layerAdditionRemoval
227 (
228 "left",
229 nMods,
230 topoChanger_,
231 "leftExtrusionFaces",
232 motionDict_.subDict("left").get<scalar>("minThickness"),
233 motionDict_.subDict("left").get<scalar>("maxThickness")
234 );
235 ++nMods;
236 tm.setSize(nMods);
237
238 Info<< "Adding " << nMods << " mesh modifiers" << endl;
239 topoChanger_.addTopologyModifiers(tm);
241 write();
242}
243
244
245// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
246
247Foam::movingConeTopoFvMesh::movingConeTopoFvMesh
248(
249 const IOobject& io,
250 const bool doInit
251)
252:
253 topoChangerFvMesh(io, doInit),
254 motionDict_
255 (
256 IOdictionary::readContents
257 (
259 (
260 "dynamicMeshDict",
261 time().constant(),
262 *this,
263 IOobject::MUST_READ
264 )
265 ).optionalSubDict(typeName + "Coeffs")
266 )
267{
268 if (doInit)
269 {
270 init(false); // do not initialise lower levels
271 }
272}
273
274
275bool Foam::movingConeTopoFvMesh::init(const bool doInit)
276{
277 if (doInit)
278 {
280 }
281
282 motionVelAmplitude_ = motionDict_.get<vector>("motionVelAmplitude");
283 motionVelPeriod_ = motionDict_.get<scalar>("motionVelPeriod");
284 curMotionVel_ =
285 motionVelAmplitude_*sin(time().value()*pi/motionVelPeriod_);
286 leftEdge_ = motionDict_.get<scalar>("leftEdge");
287 curLeft_ = motionDict_.get<scalar>("leftObstacleEdge");
288 curRight_ = motionDict_.get<scalar>("rightObstacleEdge");
289
290 Pout<< "Initial time:" << time().value()
291 << " Initial curMotionVel_:" << curMotionVel_
292 << endl;
293
294 addZonesAndModifiers();
295
296 curLeft_ = average
297 (
298 faceZones()["leftExtrusionFaces"].patch().localPoints()
299 ).x() - SMALL;
300
301 curRight_ = average
302 (
303 faceZones()["rightExtrusionFaces"].patch().localPoints()
304 ).x() + SMALL;
305
306 motionMask_ = vertexMarkup
307 (
308 points(),
309 curLeft_,
310 curRight_
311 );
312
313 // Assume something changed
314 return true;
315}
316
317
318// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
321{}
322
323
324// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
325
327{
328 // Do mesh changes (use inflation - put new points in topoChangeMap)
329 autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh(true);
330
331 // Calculate the new point positions depending on whether the
332 // topological change has happened or not
333 pointField newPoints;
334
335 vector curMotionVel_ =
336 motionVelAmplitude_*sin(time().value()*pi/motionVelPeriod_);
337
338 Pout<< "time:" << time().value() << " curMotionVel_:" << curMotionVel_
339 << " curLeft:" << curLeft_ << " curRight:" << curRight_
340 << endl;
341
342 if (topoChangeMap)
343 {
344 Info<< "Topology change. Calculating motion points" << endl;
345
346 if (topoChangeMap().hasMotionPoints())
347 {
348 Info<< "Topology change. Has premotion points" << endl;
349
350 motionMask_ =
351 vertexMarkup
352 (
353 topoChangeMap().preMotionPoints(),
354 curLeft_,
355 curRight_
356 );
357
358 // Move points inside the motionMask
359 newPoints =
360 topoChangeMap().preMotionPoints()
361 + (
362 pos0(0.5 - mag(motionMask_)) // cells above the body
363 )*curMotionVel_*time().deltaTValue();
364 }
365 else
366 {
367 Info<< "Topology change. Already set mesh points" << endl;
368
369 motionMask_ =
370 vertexMarkup
371 (
372 points(),
373 curLeft_,
374 curRight_
375 );
376
377 // Move points inside the motionMask
378 newPoints =
379 points()
380 + (
381 pos0(0.5 - mag(motionMask_)) // cells above the body
382 )*curMotionVel_*time().deltaTValue();
383 }
384 }
385 else
386 {
387 Info<< "No topology change" << endl;
388 // Set the mesh motion
389 newPoints =
390 points()
391 + (
392 pos0(0.5 - mag(motionMask_)) // cells above the body
393 )*curMotionVel_*time().deltaTValue();
394 }
395
396 // The mesh now contains the cells with zero volume
397 Info << "Executing mesh motion" << endl;
398 movePoints(newPoints);
399
400 // The mesh now has got non-zero volume cells
401
402 curLeft_ = average
403 (
404 faceZones()["leftExtrusionFaces"].patch().localPoints()
405 ).x() - SMALL;
406
407 curRight_ = average
408 (
409 faceZones()["rightExtrusionFaces"].patch().localPoints()
410 ).x() + SMALL;
411
412 return true;
413}
414
415
416// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
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.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
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
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
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
Cell layer addition mesh modifier.
Sample topoChangerFvMesh that moves an object in x direction and introduces/removes layers.
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
virtual bool update()
Update the mesh for both mesh motion and topology change.
virtual ~movingConeTopoFvMesh()
Destructor.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition polyMesh.H:671
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
Abstract base class for a topology changing fvMesh.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
const auto & io
const pointField & points
#define InfoInFunction
Report an information message using Foam::Info.
Mathematical constants.
constexpr scalar pi(M_PI)
Different types of constants.
Namespace for finite-area.
Definition limitHeight.C:30
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
dimensionedScalar pos0(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< bool > boolList
A List of bools.
Definition List.H:60
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
runTime write()
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299