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 -------------------------------------------------------------------------------
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 "movingConeTopoFvMesh.H"
30 #include "Time.H"
31 #include "mapPolyMesh.H"
32 #include "layerAdditionRemoval.H"
34 #include "meshTools.H"
35 #include "OFstream.H"
36 #include "mathematicalConstants.H"
37 
38 using namespace Foam::constant::mathematical;
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44  defineTypeNameAndDebug(movingConeTopoFvMesh, 0);
45 
47  (
48  topoChangerFvMesh,
49  movingConeTopoFvMesh,
50  IOobject
51  );
52 }
53 
54 
55 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
56 
57 Foam::tmp<Foam::scalarField> Foam::movingConeTopoFvMesh::vertexMarkup
58 (
59  const pointField& p,
60  const scalar curLeft,
61  const scalar curRight
62 ) const
63 {
64  Info<< "Updating vertex markup. curLeft: "
65  << curLeft << " curRight: " << curRight << endl;
66 
67  tmp<scalarField> tvertexMarkup(new scalarField(p.size()));
68  scalarField& vertexMarkup = tvertexMarkup.ref();
69 
70  forAll(p, pI)
71  {
72  if (p[pI].x() < curLeft - SMALL)
73  {
74  vertexMarkup[pI] = -1;
75  }
76  else if (p[pI].x() > curRight + SMALL)
77  {
78  vertexMarkup[pI] = 1;
79  }
80  else
81  {
82  vertexMarkup[pI] = 0;
83  }
84  }
85 
86  return tvertexMarkup;
87 }
88 
89 
90 void Foam::movingConeTopoFvMesh::addZonesAndModifiers()
91 {
92  // Add zones and modifiers for motion action
93 
94  if
95  (
96  pointZones().size()
97  || faceZones().size()
98  || cellZones().size()
99  || topoChanger_.size()
100  )
101  {
103  << "Zones and modifiers already present. Skipping."
104  << endl;
105 
106  return;
107  }
108 
109  Info<< "Time = " << time().timeName() << endl
110  << "Adding zones and modifiers to the mesh" << endl;
111 
112  const vectorField& fc = faceCentres();
113  const vectorField& fa = faceAreas();
114 
115  labelList zone1(fc.size());
116  boolList flipZone1(fc.size(), false);
117  label nZoneFaces1 = 0;
118 
119  labelList zone2(fc.size());
120  boolList flipZone2(fc.size(), false);
121  label nZoneFaces2 = 0;
122 
123  forAll(fc, facei)
124  {
125  if
126  (
127  fc[facei].x() > -0.003501
128  && fc[facei].x() < -0.003499
129  )
130  {
131  if ((fa[facei] & vector(1, 0, 0)) < 0)
132  {
133  flipZone1[nZoneFaces1] = true;
134  }
135 
136  zone1[nZoneFaces1] = facei;
137  Info<< "face " << facei << " for zone 1. Flip: "
138  << flipZone1[nZoneFaces1] << endl;
139  ++nZoneFaces1;
140  }
141  else if
142  (
143  fc[facei].x() > -0.00701
144  && fc[facei].x() < -0.00699
145  )
146  {
147  zone2[nZoneFaces2] = facei;
148 
149  if ((fa[facei] & vector(1, 0, 0)) > 0)
150  {
151  flipZone2[nZoneFaces2] = true;
152  }
153 
154  Info<< "face " << facei << " for zone 2. Flip: "
155  << flipZone2[nZoneFaces2] << endl;
156  ++nZoneFaces2;
157  }
158  }
159 
160  zone1.setSize(nZoneFaces1);
161  flipZone1.setSize(nZoneFaces1);
162 
163  zone2.setSize(nZoneFaces2);
164  flipZone2.setSize(nZoneFaces2);
165 
166  Info<< "zone: " << zone1 << endl;
167  Info<< "zone: " << zone2 << endl;
168 
169  List<pointZone*> pz(0);
170  List<faceZone*> fz(2);
171  List<cellZone*> cz(0);
172 
173  label nFz = 0;
174 
175  fz[nFz] =
176  new faceZone
177  (
178  "rightExtrusionFaces",
179  std::move(zone1),
180  std::move(flipZone1),
181  nFz,
182  faceZones()
183  );
184  ++nFz;
185 
186  fz[nFz] =
187  new faceZone
188  (
189  "leftExtrusionFaces",
190  std::move(zone2),
191  std::move(flipZone2),
192  nFz,
193  faceZones()
194  );
195  ++nFz;
196 
197  fz.setSize(nFz);
198 
199  Info<< "Adding mesh zones." << endl;
200  addZones(pz, fz, cz);
201 
202 
203  // Add layer addition/removal interfaces
204 
205  List<polyMeshModifier*> tm(2);
206  label nMods = 0;
207 
208  tm[nMods] =
209  new layerAdditionRemoval
210  (
211  "right",
212  nMods,
213  topoChanger_,
214  "rightExtrusionFaces",
215  motionDict_.subDict("right").get<scalar>("minThickness"),
216  motionDict_.subDict("right").get<scalar>("maxThickness")
217  );
218  ++nMods;
219 
220  tm[nMods] = new layerAdditionRemoval
221  (
222  "left",
223  nMods,
224  topoChanger_,
225  "leftExtrusionFaces",
226  motionDict_.subDict("left").get<scalar>("minThickness"),
227  motionDict_.subDict("left").get<scalar>("maxThickness")
228  );
229  ++nMods;
230  tm.setSize(nMods);
231 
232  Info<< "Adding " << nMods << " mesh modifiers" << endl;
233  topoChanger_.addTopologyModifiers(tm);
234 
235  write();
236 }
237 
238 
239 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
240 
241 Foam::movingConeTopoFvMesh::movingConeTopoFvMesh(const IOobject& io)
242 :
243  topoChangerFvMesh(io),
244  motionDict_
245  (
247  (
248  IOobject
249  (
250  "dynamicMeshDict",
251  time().constant(),
252  *this,
253  IOobject::MUST_READ_IF_MODIFIED,
254  IOobject::NO_WRITE,
255  false
256  )
257  ).optionalSubDict(typeName + "Coeffs")
258  ),
259  motionVelAmplitude_(motionDict_.get<vector>("motionVelAmplitude")),
260  motionVelPeriod_(motionDict_.get<scalar>("motionVelPeriod")),
261  curMotionVel_
262  (
263  motionVelAmplitude_*sin(time().value()*pi/motionVelPeriod_)
264  ),
265  leftEdge_(motionDict_.get<scalar>("leftEdge")),
266  curLeft_(motionDict_.get<scalar>("leftObstacleEdge")),
267  curRight_(motionDict_.get<scalar>("rightObstacleEdge"))
268 {
269  Pout<< "Initial time:" << time().value()
270  << " Initial curMotionVel_:" << curMotionVel_
271  << endl;
272 
273  addZonesAndModifiers();
274 
275  curLeft_ = average
276  (
277  faceZones()
278  [
279  faceZones().findZoneID("leftExtrusionFaces")
280  ]().localPoints()
281  ).x() - SMALL;
282 
283  curRight_ = average
284  (
285  faceZones()
286  [
287  faceZones().findZoneID("rightExtrusionFaces")
288  ]().localPoints()
289  ).x() + SMALL;
290 
291  motionMask_ = vertexMarkup
292  (
293  points(),
294  curLeft_,
295  curRight_
296  );
297 }
298 
299 
300 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
301 
303 {}
304 
305 
306 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
307 
309 {
310  // Do mesh changes (use inflation - put new points in topoChangeMap)
311  autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh(true);
312 
313  // Calculate the new point positions depending on whether the
314  // topological change has happened or not
315  pointField newPoints;
316 
317  vector curMotionVel_ =
318  motionVelAmplitude_*sin(time().value()*pi/motionVelPeriod_);
319 
320  Pout<< "time:" << time().value() << " curMotionVel_:" << curMotionVel_
321  << " curLeft:" << curLeft_ << " curRight:" << curRight_
322  << endl;
323 
324  if (topoChangeMap.valid())
325  {
326  Info<< "Topology change. Calculating motion points" << endl;
327 
328  if (topoChangeMap().hasMotionPoints())
329  {
330  Info<< "Topology change. Has premotion points" << endl;
331 
332  motionMask_ =
333  vertexMarkup
334  (
335  topoChangeMap().preMotionPoints(),
336  curLeft_,
337  curRight_
338  );
339 
340  // Move points inside the motionMask
341  newPoints =
342  topoChangeMap().preMotionPoints()
343  + (
344  pos0(0.5 - mag(motionMask_)) // cells above the body
345  )*curMotionVel_*time().deltaT().value();
346  }
347  else
348  {
349  Info<< "Topology change. Already set mesh points" << endl;
350 
351  motionMask_ =
352  vertexMarkup
353  (
354  points(),
355  curLeft_,
356  curRight_
357  );
358 
359  // Move points inside the motionMask
360  newPoints =
361  points()
362  + (
363  pos0(0.5 - mag(motionMask_)) // cells above the body
364  )*curMotionVel_*time().deltaT().value();
365  }
366  }
367  else
368  {
369  Info<< "No topology change" << endl;
370  // Set the mesh motion
371  newPoints =
372  points()
373  + (
374  pos0(0.5 - mag(motionMask_)) // cells above the body
375  )*curMotionVel_*time().deltaT().value();
376  }
377 
378  // The mesh now contains the cells with zero volume
379  Info << "Executing mesh motion" << endl;
380  movePoints(newPoints);
381 
382  // The mesh now has got non-zero volume cells
383 
384  curLeft_ = average
385  (
386  faceZones()
387  [
388  faceZones().findZoneID("leftExtrusionFaces")
389  ]().localPoints()
390  ).x() - SMALL;
391 
392  curRight_ = average
393  (
394  faceZones()
395  [
396  faceZones().findZoneID("rightExtrusionFaces")
397  ]().localPoints()
398  ).x() + SMALL;
399 
400  return true;
401 }
402 
403 
404 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1038
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
meshTools.H
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:320
mathematicalConstants.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::movingConeTopoFvMesh::~movingConeTopoFvMesh
virtual ~movingConeTopoFvMesh()
Destructor.
Definition: movingConeTopoFvMesh.C:302
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
mapPolyMesh.H
movingConeTopoFvMesh.H
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:69
Foam::autoPtr::valid
bool valid() const noexcept
True if the managed pointer is non-null.
Definition: autoPtr.H:148
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::movingConeTopoFvMesh::update
virtual bool update()
Update the mesh for both mesh motion and topology change.
Definition: movingConeTopoFvMesh.C:308
layerAdditionRemoval.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
OFstream.H
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:477
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Time.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::constant::mathematical
Mathematical constants.
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::Vector< scalar >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::topoChangerFvMesh
Abstract base class for a topology changing fvMesh.
Definition: topoChangerFvMesh.H:52
points
const pointField & points
Definition: gmvOutputHeader.H:1
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:35
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:248
constant
constant condensation/saturation model.
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:328