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 | Copyright (C) 2011-2014 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25  Cell layer addition/removal mesh modifier
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "layerAdditionRemoval.H"
30 #include "polyTopoChanger.H"
31 #include "polyMesh.H"
32 #include "Time.H"
33 #include "primitiveMesh.H"
34 #include "polyTopoChange.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(layerAdditionRemoval, 0);
43  (
44  polyMeshModifier,
45  layerAdditionRemoval,
46  dictionary
47  );
48 }
49 
50 
51 const Foam::scalar Foam::layerAdditionRemoval::addDelta_ = 0.3;
52 const Foam::scalar Foam::layerAdditionRemoval::removeDelta_ = 0.1;
53 
54 
55 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
56 
57 void Foam::layerAdditionRemoval::checkDefinition()
58 {
59  if (!faceZoneID_.active())
60  {
61  FatalErrorIn("void Foam::layerAdditionRemoval::checkDefinition()")
62  << "Master face zone named " << faceZoneID_.name()
63  << " cannot be found."
64  << abort(FatalError);
65  }
66 
67  if
68  (
69  minLayerThickness_ < VSMALL
70  || maxLayerThickness_ < minLayerThickness_
71  )
72  {
74  (
75  "void Foam::layerAdditionRemoval::checkDefinition()"
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  {
86  FatalErrorIn("void Foam::layerAdditionRemoval::checkDefinition()")
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 Switch 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_(NULL),
155  facesPairingPtr_(NULL),
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, Switch(dict.lookup("active"))),
172  faceZoneID_(dict.lookup("faceZoneName"), ptc.mesh().faceZones()),
173  minLayerThickness_(readScalar(dict.lookup("minLayerThickness"))),
174  maxLayerThickness_(readScalar(dict.lookup("maxLayerThickness"))),
175  thicknessFromVolume_
176  (
177  dict.lookupOrDefault<Switch>("thicknessFromVolume", true)
178  ),
179  oldLayerThickness_(readOldThickness(dict)),
180  pointsPairingPtr_(NULL),
181  facesPairingPtr_(NULL),
182  triggerRemoval_(-1),
183  triggerAddition_(-1)
184 {
185  checkDefinition();
186 }
187 
188 
189 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
190 
192 {
193  clearAddressing();
194 }
195 
196 
197 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
198 
200 {
201  // Protect from multiple calculation in the same time-step
202  if (triggerRemoval_ > -1 || triggerAddition_ > -1)
203  {
204  return true;
205  }
206 
207  // Go through all the cells in the master layer and calculate
208  // approximate layer thickness as the ratio of the cell volume and
209  // face area in the face zone.
210  // Layer addition:
211  // When the max thickness exceeds the threshold, trigger refinement.
212  // Layer removal:
213  // When the min thickness falls below the threshold, trigger removal.
214 
215  const polyMesh& mesh = topoChanger().mesh();
216 
217  const faceZone& fz = mesh.faceZones()[faceZoneID_.index()];
218  const labelList& mc = fz.masterCells();
219 
220  const scalarField& V = mesh.cellVolumes();
221  const vectorField& S = mesh.faceAreas();
222 
223  if (min(V) < -VSMALL)
224  {
225  FatalErrorIn("bool layerAdditionRemoval::changeTopology() const")
226  << "negative cell volume. Error in mesh motion before "
227  << "topological change.\n V: " << V
228  << abort(FatalError);
229  }
230 
231  scalar avgDelta = 0;
232  scalar minDelta = GREAT;
233  scalar maxDelta = 0;
234  label nDelta = 0;
235 
236  if (thicknessFromVolume_)
237  {
238  // Thickness calculated from cell volume/face area
239  forAll(fz, faceI)
240  {
241  scalar curDelta = V[mc[faceI]]/mag(S[fz[faceI]]);
242  avgDelta += curDelta;
243  minDelta = min(minDelta, curDelta);
244  maxDelta = max(maxDelta, curDelta);
245  }
246 
247  nDelta = fz.size();
248  }
249  else
250  {
251  // Thickness calculated from edges on layer
252  const Map<label>& zoneMeshPointMap = fz().meshPointMap();
253 
254  // Edges with only one point on zone
255  forAll(mc, faceI)
256  {
257  const cell& cFaces = mesh.cells()[mc[faceI]];
258  const edgeList cellEdges(cFaces.edges(mesh.faces()));
259 
260  forAll(cellEdges, i)
261  {
262  const edge& e = cellEdges[i];
263 
264  if (zoneMeshPointMap.found(e[0]))
265  {
266  if (!zoneMeshPointMap.found(e[1]))
267  {
268  scalar curDelta = e.mag(mesh.points());
269  avgDelta += curDelta;
270  nDelta++;
271  minDelta = min(minDelta, curDelta);
272  maxDelta = max(maxDelta, curDelta);
273  }
274  }
275  else
276  {
277  if (zoneMeshPointMap.found(e[1]))
278  {
279  scalar curDelta = e.mag(mesh.points());
280  avgDelta += curDelta;
281  nDelta++;
282  minDelta = min(minDelta, curDelta);
283  maxDelta = max(maxDelta, curDelta);
284  }
285  }
286  }
287  }
288  }
289 
290  reduce(minDelta, minOp<scalar>());
291  reduce(maxDelta, maxOp<scalar>());
292  reduce(avgDelta, sumOp<scalar>());
293  reduce(nDelta, sumOp<label>());
294 
295  avgDelta /= nDelta;
296 
297  if (debug)
298  {
299  Pout<< "bool layerAdditionRemoval::changeTopology() const "
300  << " for object " << name() << " : " << nl
301  << "Layer thickness: min: " << minDelta
302  << " max: " << maxDelta << " avg: " << avgDelta
303  << " old thickness: " << oldLayerThickness_ << nl
304  << "Removal threshold: " << minLayerThickness_
305  << " addition threshold: " << maxLayerThickness_ << endl;
306  }
307 
308  bool topologicalChange = false;
309 
310  // If the thickness is decreasing and crosses the min thickness,
311  // trigger removal
312  if (oldLayerThickness_ < 0)
313  {
314  if (debug)
315  {
316  Pout<< "First step. No addition/removal" << endl;
317  }
318 
319  // No topological changes allowed before first mesh motion
320  oldLayerThickness_ = avgDelta;
321 
322  topologicalChange = false;
323  }
324  else if (avgDelta < oldLayerThickness_)
325  {
326  // Layers moving towards removal
327  if (minDelta < minLayerThickness_)
328  {
329  // Check layer pairing
330  if (setLayerPairing())
331  {
332  // A mesh layer detected. Check that collapse is valid
333  if (validCollapse())
334  {
335  // At this point, info about moving the old mesh
336  // in a way to collapse the cells in the removed
337  // layer is available. Not sure what to do with
338  // it.
339 
340  if (debug)
341  {
342  Pout<< "bool layerAdditionRemoval::changeTopology() "
343  << " const for object " << name() << " : "
344  << "Triggering layer removal" << endl;
345  }
346 
347  triggerRemoval_ = mesh.time().timeIndex();
348 
349  // Old thickness looses meaning.
350  // Set it up to indicate layer removal
351  oldLayerThickness_ = GREAT;
352 
353  topologicalChange = true;
354  }
355  else
356  {
357  // No removal, clear addressing
358  clearAddressing();
359  }
360  }
361  }
362  else
363  {
364  oldLayerThickness_ = avgDelta;
365  }
366  }
367  else
368  {
369  // Layers moving towards addition
370  if (maxDelta > maxLayerThickness_)
371  {
372  if (debug)
373  {
374  Pout<< "bool layerAdditionRemoval::changeTopology() const "
375  << " for object " << name() << " : "
376  << "Triggering layer addition" << endl;
377  }
378 
379  triggerAddition_ = mesh.time().timeIndex();
380 
381  // Old thickness looses meaning.
382  // Set it up to indicate layer removal
383  oldLayerThickness_ = 0;
384 
385  topologicalChange = true;
386  }
387  else
388  {
389  oldLayerThickness_ = avgDelta;
390  }
391  }
392 
393  return topologicalChange;
394 }
395 
396 
398 {
399  // Insert the layer addition/removal instructions
400  // into the topological change
401 
402  if (triggerRemoval_ == topoChanger().mesh().time().timeIndex())
403  {
404  removeCellLayer(ref);
405 
406  // Clear addressing. This also resets the addition/removal data
407  if (debug)
408  {
409  Pout<< "layerAdditionRemoval::setRefinement(polyTopoChange&) "
410  << "for object " << name() << " : "
411  << "Clearing addressing after layer removal" << endl;
412  }
413 
414  triggerRemoval_ = -1;
415  clearAddressing();
416  }
417 
418  if (triggerAddition_ == topoChanger().mesh().time().timeIndex())
419  {
420  addCellLayer(ref);
421 
422  // Clear addressing. This also resets the addition/removal data
423  if (debug)
424  {
425  Pout<< "layerAdditionRemoval::setRefinement(polyTopoChange&) "
426  << "for object " << name() << " : "
427  << "Clearing addressing after layer addition" << endl;
428  }
429 
430  triggerAddition_ = -1;
431  clearAddressing();
432  }
433 }
434 
435 
437 {
438  if (debug)
439  {
440  Pout<< "layerAdditionRemoval::updateMesh(const mapPolyMesh&) "
441  << "for object " << name() << " : "
442  << "Clearing addressing on external request";
443 
444  if (pointsPairingPtr_ || facesPairingPtr_)
445  {
446  Pout<< "Pointers set." << endl;
447  }
448  else
449  {
450  Pout<< "Pointers not set." << endl;
451  }
452  }
453 
454  // Mesh has changed topologically. Update local topological data
455  faceZoneID_.update(topoChanger().mesh().faceZones());
456 
457  clearAddressing();
458 }
459 
460 
462 {
463  if (t < VSMALL || maxLayerThickness_ < t)
464  {
466  (
467  "void layerAdditionRemoval::setMinLayerThickness"
468  "("
469  "const scalar"
470  ") const"
471  ) << "Incorrect layer thickness definition."
472  << abort(FatalError);
473  }
474 
475  minLayerThickness_ = t;
476 }
477 
478 
480 {
481  if (t < minLayerThickness_)
482  {
484  (
485  "void layerAdditionRemoval::setMaxLayerThickness"
486  "("
487  "const scalar"
488  ") const"
489  ) << "Incorrect layer thickness definition."
490  << abort(FatalError);
491  }
492 
493  maxLayerThickness_ = t;
494 }
495 
496 
498 {
499  os << nl << type() << nl
500  << name()<< nl
501  << faceZoneID_ << nl
502  << minLayerThickness_ << nl
503  << oldLayerThickness_ << nl
504  << maxLayerThickness_ << nl
505  << thicknessFromVolume_ << endl;
506 }
507 
508 
510 {
511  os << nl << name() << nl << token::BEGIN_BLOCK << nl
512  << " type " << type()
514  << " faceZoneName " << faceZoneID_.name()
516  << " minLayerThickness " << minLayerThickness_
518  << " maxLayerThickness " << maxLayerThickness_
520  << " thicknessFromVolume " << thicknessFromVolume_
522  << " oldLayerThickness " << oldLayerThickness_
524  << " active " << active()
526  << token::END_BLOCK << endl;
527 }
528 
529 
530 // ************************************************************************* //
const vectorField & faceAreas() const
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
edgeList edges(const faceUList &) const
Return cell edges.
Definition: cell.C:120
dimensioned< scalar > mag(const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Virtual base class for mesh modifiers.
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
const labelList & masterCells() const
Return labels of master cells (cells next to the master face.
Definition: faceZone.C:334
void deleteDemandDrivenData(DataPtr &dataPtr)
label timeIndex
Definition: getTimeIndex.H:4
A class for handling words, derived from string.
Definition: word.H:59
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:59
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
const cellList & cells() const
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
const polyMesh & mesh() const
Return the mesh reference.
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Namespace for OpenFOAM.
virtual ~layerAdditionRemoval()
Destructor.
scalar mag(const pointField &) const
Return scalar magnitude.
Definition: edgeI.H:163
static const char nl
Definition: Ostream.H:260
const double e
Elementary charge.
Definition: doubleFloat.H:78
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void writeDict(Ostream &) const
Write dictionary.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1035
label timeIndex() const
Return current time index.
Definition: TimeState.C:73
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define forAll(list, i)
Definition: UList.H:421
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Macros for easy insertion into run-time selection tables.
List of mesh modifiers defining the mesh dynamics.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
error FatalError
virtual void write(Ostream &) const
Write.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
Direct mesh changes based on v1.3 polyTopoChange syntax.
virtual void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
const Time & time() const
Return time.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1060
const scalarField & cellVolumes() const
void setMaxLayerThickness(const scalar t) const
Set max layer thickness which triggers removal.
void setMinLayerThickness(const scalar t) const
Set min layer thickness which triggers removal.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual void setRefinement(polyTopoChange &) const
Insert the layer addition/removal instructions.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
virtual bool changeTopology() const
Check for topology change.
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53