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  );
53  (
54  topoChangerFvMesh,
55  movingConeTopoFvMesh,
56  doInit
57  );
58 }
59 
60 
61 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
62 
63 Foam::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  tmp<scalarField> tvertexMarkup(new scalarField(p.size()));
74  scalarField& 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 
96 void 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] & vector(1, 0, 0)) < 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] & vector(1, 0, 0)) > 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 
211  List<polyMeshModifier*> tm(2);
212  label nMods = 0;
213 
214  tm[nMods] =
215  new layerAdditionRemoval
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);
240 
241  write();
242 }
243 
244 
245 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
246 
247 Foam::movingConeTopoFvMesh::movingConeTopoFvMesh
248 (
249  const IOobject& io,
250  const bool doInit
251 )
252 :
253  topoChangerFvMesh(io, doInit),
254  motionDict_
255  (
257  (
258  IOobject
259  (
260  "dynamicMeshDict",
261  time().constant(),
262  *this,
265  false
266  )
267  ).optionalSubDict(typeName + "Coeffs")
268  )
269 {
270  if (doInit)
271  {
272  init(false); // do not initialise lower levels
273  }
274 }
275 
276 
277 bool Foam::movingConeTopoFvMesh::init(const bool doInit)
278 {
279  if (doInit)
280  {
281  topoChangerFvMesh::init(doInit);
282  }
283 
284  motionVelAmplitude_ = motionDict_.get<vector>("motionVelAmplitude");
285  motionVelPeriod_ = motionDict_.get<scalar>("motionVelPeriod");
286  curMotionVel_ =
287  motionVelAmplitude_*sin(time().value()*pi/motionVelPeriod_);
288  leftEdge_ = motionDict_.get<scalar>("leftEdge");
289  curLeft_ = motionDict_.get<scalar>("leftObstacleEdge");
290  curRight_ = motionDict_.get<scalar>("rightObstacleEdge");
291 
292  Pout<< "Initial time:" << time().value()
293  << " Initial curMotionVel_:" << curMotionVel_
294  << endl;
295 
296  addZonesAndModifiers();
297 
298  curLeft_ = average
299  (
300  faceZones()
301  [
302  faceZones().findZoneID("leftExtrusionFaces")
303  ]().localPoints()
304  ).x() - SMALL;
305 
306  curRight_ = average
307  (
308  faceZones()
309  [
310  faceZones().findZoneID("rightExtrusionFaces")
311  ]().localPoints()
312  ).x() + SMALL;
313 
314  motionMask_ = vertexMarkup
315  (
316  points(),
317  curLeft_,
318  curRight_
319  );
320 
321  // Assume something changed
322  return true;
323 }
324 
325 
326 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
327 
329 {}
330 
331 
332 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
333 
335 {
336  // Do mesh changes (use inflation - put new points in topoChangeMap)
337  autoPtr<mapPolyMesh> topoChangeMap = topoChanger_.changeMesh(true);
338 
339  // Calculate the new point positions depending on whether the
340  // topological change has happened or not
341  pointField newPoints;
342 
343  vector curMotionVel_ =
344  motionVelAmplitude_*sin(time().value()*pi/motionVelPeriod_);
345 
346  Pout<< "time:" << time().value() << " curMotionVel_:" << curMotionVel_
347  << " curLeft:" << curLeft_ << " curRight:" << curRight_
348  << endl;
349 
350  if (topoChangeMap)
351  {
352  Info<< "Topology change. Calculating motion points" << endl;
353 
354  if (topoChangeMap().hasMotionPoints())
355  {
356  Info<< "Topology change. Has premotion points" << endl;
357 
358  motionMask_ =
359  vertexMarkup
360  (
361  topoChangeMap().preMotionPoints(),
362  curLeft_,
363  curRight_
364  );
365 
366  // Move points inside the motionMask
367  newPoints =
368  topoChangeMap().preMotionPoints()
369  + (
370  pos0(0.5 - mag(motionMask_)) // cells above the body
371  )*curMotionVel_*time().deltaT().value();
372  }
373  else
374  {
375  Info<< "Topology change. Already set mesh points" << endl;
376 
377  motionMask_ =
378  vertexMarkup
379  (
380  points(),
381  curLeft_,
382  curRight_
383  );
384 
385  // Move points inside the motionMask
386  newPoints =
387  points()
388  + (
389  pos0(0.5 - mag(motionMask_)) // cells above the body
390  )*curMotionVel_*time().deltaT().value();
391  }
392  }
393  else
394  {
395  Info<< "No topology change" << endl;
396  // Set the mesh motion
397  newPoints =
398  points()
399  + (
400  pos0(0.5 - mag(motionMask_)) // cells above the body
401  )*curMotionVel_*time().deltaT().value();
402  }
403 
404  // The mesh now contains the cells with zero volume
405  Info << "Executing mesh motion" << endl;
406  movePoints(newPoints);
407 
408  // The mesh now has got non-zero volume cells
409 
410  curLeft_ = average
411  (
412  faceZones()
413  [
414  faceZones().findZoneID("leftExtrusionFaces")
415  ]().localPoints()
416  ).x() - SMALL;
417 
418  curRight_ = average
419  (
420  faceZones()
421  [
422  faceZones().findZoneID("rightExtrusionFaces")
423  ]().localPoints()
424  ).x() + SMALL;
425 
426  return true;
427 }
428 
429 
430 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
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::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:169
meshTools.H
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:350
mathematicalConstants.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::movingConeTopoFvMesh::~movingConeTopoFvMesh
virtual ~movingConeTopoFvMesh()
Destructor.
Definition: movingConeTopoFvMesh.C:328
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
mapPolyMesh.H
movingConeTopoFvMesh.H
Foam::dynamicFvMesh::init
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: dynamicFvMesh.C:90
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:65
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::movingConeTopoFvMesh::update
virtual bool update()
Update the mesh for both mesh motion and topology change.
Definition: movingConeTopoFvMesh.C:334
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::Field< vector >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
init
mesh init(true)
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::movingConeTopoFvMesh::init
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: movingConeTopoFvMesh.C:277
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:53
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::IOobject::MUST_READ_IF_MODIFIED
Definition: IOobject.H:186
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:36
constant
constant condensation/saturation model.
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:328