layerAdditionRemoval.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-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Description
27  Cell layer addition/removal mesh modifier
28 
29 \*---------------------------------------------------------------------------*/
30 
31 #include "layerAdditionRemoval.H"
32 #include "polyTopoChanger.H"
33 #include "polyMesh.H"
34 #include "Time.H"
35 #include "primitiveMesh.H"
36 #include "polyTopoChange.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(layerAdditionRemoval, 0);
45  (
46  polyMeshModifier,
47  layerAdditionRemoval,
48  dictionary
49  );
50 }
51 
52 
53 const Foam::scalar Foam::layerAdditionRemoval::addDelta_ = 0.3;
54 const Foam::scalar Foam::layerAdditionRemoval::removeDelta_ = 0.1;
55 
56 
57 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
58 
59 void Foam::layerAdditionRemoval::checkDefinition()
60 {
61  if (!faceZoneID_.active())
62  {
64  << "Master face zone named " << faceZoneID_.name()
65  << " cannot be found."
66  << abort(FatalError);
67  }
68 
69  if
70  (
71  minLayerThickness_ < VSMALL
72  || maxLayerThickness_ < minLayerThickness_
73  )
74  {
76  << "Incorrect layer thickness definition."
77  << abort(FatalError);
78  }
79 
80  label nFaces = topoChanger().mesh().faceZones()[faceZoneID_.index()].size();
81 
82  reduce(nFaces, sumOp<label>());
83 
84  if (nFaces == 0)
85  {
87  << "Face extrusion zone contains no faces. "
88  << "Please check your mesh definition."
89  << abort(FatalError);
90  }
91 
92  if (debug)
93  {
94  Pout<< "Cell layer addition/removal object " << name() << " :" << nl
95  << " faceZoneID: " << faceZoneID_ << endl;
96  }
97 }
98 
99 
100 Foam::scalar Foam::layerAdditionRemoval::readOldThickness
101 (
102  const dictionary& dict
103 )
104 {
105  return dict.lookupOrDefault("oldLayerThickness", -1.0);
106 }
107 
108 
109 void Foam::layerAdditionRemoval::clearAddressing() const
110 {
111  if (pointsPairingPtr_)
112  {
113  if (debug)
114  {
115  Pout<< "layerAdditionRemoval::clearAddressing()" << nl
116  << " clearing pointsPairingPtr_" << endl;
117  }
118 
119  deleteDemandDrivenData(pointsPairingPtr_);
120  }
121 
122  if (facesPairingPtr_)
123  {
124  if (debug)
125  {
126  Pout<< "layerAdditionRemoval::clearAddressing()" << nl
127  << " clearing facesPairingPtr_" << endl;
128  }
129 
130  deleteDemandDrivenData(facesPairingPtr_);
131  }
132 }
133 
134 
135 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
136 
137 Foam::layerAdditionRemoval::layerAdditionRemoval
138 (
139  const word& name,
140  const label index,
141  const polyTopoChanger& ptc,
142  const word& zoneName,
143  const scalar minThickness,
144  const scalar maxThickness,
145  const bool thicknessFromVolume
146 )
147 :
148  polyMeshModifier(name, index, ptc, true),
149  faceZoneID_(zoneName, ptc.mesh().faceZones()),
150  minLayerThickness_(minThickness),
151  maxLayerThickness_(maxThickness),
152  thicknessFromVolume_(thicknessFromVolume),
153  oldLayerThickness_(-1.0),
154  pointsPairingPtr_(nullptr),
155  facesPairingPtr_(nullptr),
156  triggerRemoval_(-1),
157  triggerAddition_(-1)
158 {
159  checkDefinition();
160 }
161 
162 
163 Foam::layerAdditionRemoval::layerAdditionRemoval
164 (
165  const word& name,
166  const dictionary& dict,
167  const label index,
168  const polyTopoChanger& ptc
169 )
170 :
171  polyMeshModifier(name, index, ptc, dict.get<bool>("active")),
172  faceZoneID_(dict.lookup("faceZoneName"), ptc.mesh().faceZones()),
173  minLayerThickness_(dict.get<scalar>("minLayerThickness")),
174  maxLayerThickness_(dict.get<scalar>("maxLayerThickness")),
175  thicknessFromVolume_(dict.lookupOrDefault("thicknessFromVolume", true)),
176  oldLayerThickness_(readOldThickness(dict)),
177  pointsPairingPtr_(nullptr),
178  facesPairingPtr_(nullptr),
179  triggerRemoval_(-1),
180  triggerAddition_(-1)
181 {
182  checkDefinition();
183 }
184 
185 
186 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
187 
189 {
190  clearAddressing();
191 }
192 
193 
194 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
195 
197 {
198  // Protect from multiple calculation in the same time-step
199  if (triggerRemoval_ > -1 || triggerAddition_ > -1)
200  {
201  return true;
202  }
203 
204  // Go through all the cells in the master layer and calculate
205  // approximate layer thickness as the ratio of the cell volume and
206  // face area in the face zone.
207  // Layer addition:
208  // When the max thickness exceeds the threshold, trigger refinement.
209  // Layer removal:
210  // When the min thickness falls below the threshold, trigger removal.
211 
212  const polyMesh& mesh = topoChanger().mesh();
213 
214  const faceZone& fz = mesh.faceZones()[faceZoneID_.index()];
215  const labelList& mc = fz.masterCells();
216 
217  const scalarField& V = mesh.cellVolumes();
218  const vectorField& S = mesh.faceAreas();
219 
220  if (min(V) < -VSMALL)
221  {
223  << "negative cell volume. Error in mesh motion before "
224  << "topological change.\n V: " << V
225  << abort(FatalError);
226  }
227 
228  scalar avgDelta = 0;
229  scalar minDelta = GREAT;
230  scalar maxDelta = 0;
231  label nDelta = 0;
232 
233  if (thicknessFromVolume_)
234  {
235  // Thickness calculated from cell volume/face area
236  forAll(fz, facei)
237  {
238  scalar curDelta = V[mc[facei]]/mag(S[fz[facei]]);
239  avgDelta += curDelta;
240  minDelta = min(minDelta, curDelta);
241  maxDelta = max(maxDelta, curDelta);
242  }
243 
244  nDelta = fz.size();
245  }
246  else
247  {
248  // Thickness calculated from edges on layer
249  const Map<label>& zoneMeshPointMap = fz().meshPointMap();
250 
251  // Edges with only one point on zone
252  forAll(mc, facei)
253  {
254  const cell& cFaces = mesh.cells()[mc[facei]];
255  const edgeList cellEdges(cFaces.edges(mesh.faces()));
256 
257  forAll(cellEdges, i)
258  {
259  const edge& e = cellEdges[i];
260 
261  if (zoneMeshPointMap.found(e[0]))
262  {
263  if (!zoneMeshPointMap.found(e[1]))
264  {
265  scalar curDelta = e.mag(mesh.points());
266  avgDelta += curDelta;
267  nDelta++;
268  minDelta = min(minDelta, curDelta);
269  maxDelta = max(maxDelta, curDelta);
270  }
271  }
272  else
273  {
274  if (zoneMeshPointMap.found(e[1]))
275  {
276  scalar curDelta = e.mag(mesh.points());
277  avgDelta += curDelta;
278  nDelta++;
279  minDelta = min(minDelta, curDelta);
280  maxDelta = max(maxDelta, curDelta);
281  }
282  }
283  }
284  }
285  }
286 
287  reduce(minDelta, minOp<scalar>());
288  reduce(maxDelta, maxOp<scalar>());
289  reduce(avgDelta, sumOp<scalar>());
290  reduce(nDelta, sumOp<label>());
291 
292  avgDelta /= nDelta;
293 
294  if (debug)
295  {
296  Pout<< "bool layerAdditionRemoval::changeTopology() const "
297  << " for object " << name() << " : " << nl
298  << "Layer thickness: min: " << minDelta
299  << " max: " << maxDelta << " avg: " << avgDelta
300  << " old thickness: " << oldLayerThickness_ << nl
301  << "Removal threshold: " << minLayerThickness_
302  << " addition threshold: " << maxLayerThickness_ << endl;
303  }
304 
305  bool topologicalChange = false;
306 
307  // If the thickness is decreasing and crosses the min thickness,
308  // trigger removal
309  if (oldLayerThickness_ < 0)
310  {
311  if (debug)
312  {
313  Pout<< "First step. No addition/removal" << endl;
314  }
315 
316  // No topological changes allowed before first mesh motion
317  oldLayerThickness_ = avgDelta;
318 
319  topologicalChange = false;
320  }
321  else if (avgDelta < oldLayerThickness_)
322  {
323  // Layers moving towards removal
324  if (minDelta < minLayerThickness_)
325  {
326  // Check layer pairing
327  if (setLayerPairing())
328  {
329  // A mesh layer detected. Check that collapse is valid
330  if (validCollapse())
331  {
332  // At this point, info about moving the old mesh
333  // in a way to collapse the cells in the removed
334  // layer is available. Not sure what to do with
335  // it.
336 
337  if (debug)
338  {
339  Pout<< "bool layerAdditionRemoval::changeTopology() "
340  << " const for object " << name() << " : "
341  << "Triggering layer removal" << endl;
342  }
343 
344  triggerRemoval_ = mesh.time().timeIndex();
345 
346  // Old thickness looses meaning.
347  // Set it up to indicate layer removal
348  oldLayerThickness_ = GREAT;
349 
350  topologicalChange = true;
351  }
352  else
353  {
354  // No removal, clear addressing
355  clearAddressing();
356  }
357  }
358  }
359  else
360  {
361  oldLayerThickness_ = avgDelta;
362  }
363  }
364  else
365  {
366  // Layers moving towards addition
367  if (maxDelta > maxLayerThickness_)
368  {
369  if (debug)
370  {
371  Pout<< "bool layerAdditionRemoval::changeTopology() const "
372  << " for object " << name() << " : "
373  << "Triggering layer addition" << endl;
374  }
375 
376  triggerAddition_ = mesh.time().timeIndex();
377 
378  // Old thickness looses meaning.
379  // Set it up to indicate layer removal
380  oldLayerThickness_ = 0;
381 
382  topologicalChange = true;
383  }
384  else
385  {
386  oldLayerThickness_ = avgDelta;
387  }
388  }
389 
390  return topologicalChange;
391 }
392 
393 
395 {
396  // Insert the layer addition/removal instructions
397  // into the topological change
398 
399  if (triggerRemoval_ == topoChanger().mesh().time().timeIndex())
400  {
401  removeCellLayer(ref);
402 
403  // Clear addressing. This also resets the addition/removal data
404  if (debug)
405  {
406  Pout<< "layerAdditionRemoval::setRefinement(polyTopoChange&) "
407  << "for object " << name() << " : "
408  << "Clearing addressing after layer removal" << endl;
409  }
410 
411  triggerRemoval_ = -1;
412  clearAddressing();
413  }
414 
415  if (triggerAddition_ == topoChanger().mesh().time().timeIndex())
416  {
417  addCellLayer(ref);
418 
419  // Clear addressing. This also resets the addition/removal data
420  if (debug)
421  {
422  Pout<< "layerAdditionRemoval::setRefinement(polyTopoChange&) "
423  << "for object " << name() << " : "
424  << "Clearing addressing after layer addition" << endl;
425  }
426 
427  triggerAddition_ = -1;
428  clearAddressing();
429  }
430 }
431 
432 
434 {
435  if (debug)
436  {
437  Pout<< "layerAdditionRemoval::updateMesh(const mapPolyMesh&) "
438  << "for object " << name() << " : "
439  << "Clearing addressing on external request";
440 
441  if (pointsPairingPtr_ || facesPairingPtr_)
442  {
443  Pout<< "Pointers set." << endl;
444  }
445  else
446  {
447  Pout<< "Pointers not set." << endl;
448  }
449  }
450 
451  // Mesh has changed topologically. Update local topological data
452  faceZoneID_.update(topoChanger().mesh().faceZones());
453 
454  clearAddressing();
455 }
456 
457 
459 {
460  if (t < VSMALL || maxLayerThickness_ < t)
461  {
463  << "Incorrect layer thickness definition."
464  << abort(FatalError);
465  }
466 
467  minLayerThickness_ = t;
468 }
469 
470 
472 {
473  if (t < minLayerThickness_)
474  {
476  << "Incorrect layer thickness definition."
477  << abort(FatalError);
478  }
479 
480  maxLayerThickness_ = t;
481 }
482 
483 
485 {
486  os << nl << type() << nl
487  << name()<< nl
488  << faceZoneID_ << nl
489  << minLayerThickness_ << nl
490  << oldLayerThickness_ << nl
491  << maxLayerThickness_ << nl
492  << thicknessFromVolume_ << endl;
493 }
494 
495 
497 {
498  os << nl << name() << nl << token::BEGIN_BLOCK << nl
499  << " type " << type()
501  << " faceZoneName " << faceZoneID_.name()
503  << " minLayerThickness " << minLayerThickness_
505  << " maxLayerThickness " << maxLayerThickness_
507  << " thicknessFromVolume " << thicknessFromVolume_
509  << " oldLayerThickness " << oldLayerThickness_
511  << " active " << active()
513  << token::END_BLOCK << endl;
514 }
515 
516 
517 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::maxOp
Definition: ops.H:223
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1038
ref
rDeltaT ref()
Foam::layerAdditionRemoval::changeTopology
virtual bool changeTopology() const
Check for topology change.
Definition: layerAdditionRemoval.C:196
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::layerAdditionRemoval::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: layerAdditionRemoval.C:433
Foam::polyMeshModifier::topoChanger
const polyTopoChanger & topoChanger() const
Return reference to morph engine.
Definition: polyMeshModifier.C:71
Foam::minOp
Definition: ops.H:224
Foam::polyTopoChanger
List of mesh modifiers defining the mesh dynamics.
Definition: polyTopoChanger.H:67
polyTopoChanger.H
polyTopoChange.H
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:100
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::layerAdditionRemoval::setRefinement
virtual void setRefinement(polyTopoChange &) const
Insert the layer addition/removal instructions.
Definition: layerAdditionRemoval.C:394
Foam::Map< label >
primitiveMesh.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::cell::edges
edgeList edges(const faceUList &f) const
Return cell edges.
Definition: cell.C:122
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1241
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
polyMesh.H
layerAdditionRemoval.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::sumOp
Definition: ops.H:213
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::layerAdditionRemoval::~layerAdditionRemoval
virtual ~layerAdditionRemoval()
Destructor.
Definition: layerAdditionRemoval.C:188
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::layerAdditionRemoval::setMaxLayerThickness
void setMaxLayerThickness(const scalar t) const
Set max layer thickness which triggers removal.
Definition: layerAdditionRemoval.C:471
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
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::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:65
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::layerAdditionRemoval::setMinLayerThickness
void setMinLayerThickness(const scalar t) const
Set min layer thickness which triggers removal.
Definition: layerAdditionRemoval.C:458
Foam::token::END_STATEMENT
End entry [isseparator].
Definition: token.H:116
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::token::END_BLOCK
End block [isseparator].
Definition: token.H:122
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:419
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
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::primitiveMesh::cellVolumes
const scalarField & cellVolumes() const
Definition: primitiveMeshCellCentresAndVols.C:187
Foam::token::BEGIN_BLOCK
Begin block [isseparator].
Definition: token.H:121
Foam::polyMeshModifier
Virtual base class for mesh modifiers.
Definition: polyMeshModifier.H:71
Foam::polyMeshModifier::name
const word & name() const
Return name of this modifier.
Definition: polyMeshModifier.H:153
Time.H
Foam::faceZone::masterCells
const labelList & masterCells() const
Definition: faceZone.C:405
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1063
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::polyTopoChanger::mesh
const polyMesh & mesh() const
Return the mesh reference.
Definition: polyTopoChanger.H:116
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:160
Foam::layerAdditionRemoval::write
virtual void write(Ostream &) const
Write.
Definition: layerAdditionRemoval.C:484
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::TimeState::timeIndex
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:36
Foam::DynamicID::name
const keyType & name() const
Return name.
Definition: DynamicID.H:98
Foam::layerAdditionRemoval::writeDict
virtual void writeDict(Ostream &) const
Write dictionary.
Definition: layerAdditionRemoval.C:496
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::DynamicID::active
bool active() const
Has the zone been found.
Definition: DynamicID.H:116
Foam::DynamicID::index
label index() const
Return index of first matching zone.
Definition: DynamicID.H:110
Foam::primitiveMesh::faceAreas
const vectorField & faceAreas() const
Definition: primitiveMeshFaceCentresAndAreas.C:156