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  Copyright (C) 2020-2023 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 "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  {
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  << "Incorrect layer thickness definition."
75  << abort(FatalError);
76  }
77 
78  if
79  (
81  (
82  topoChanger().mesh().faceZones()[faceZoneID_.index()].empty()
83  )
84  )
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 void Foam::layerAdditionRemoval::clearAddressing() const
101 {
102  pointsPairingPtr_.reset(nullptr);
103  facesPairingPtr_.reset(nullptr);
104 }
105 
106 
107 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
108 
109 Foam::layerAdditionRemoval::layerAdditionRemoval
110 (
111  const word& name,
112  const label index,
113  const polyTopoChanger& ptc,
114  const word& zoneName,
115  const scalar minThickness,
116  const scalar maxThickness,
117  const bool thicknessFromVolume
118 )
119 :
120  polyMeshModifier(name, index, ptc, true),
121  faceZoneID_(zoneName, ptc.mesh().faceZones()),
122  minLayerThickness_(minThickness),
123  maxLayerThickness_(maxThickness),
124  thicknessFromVolume_(thicknessFromVolume),
125  oldLayerThickness_(-1.0),
126  pointsPairingPtr_(nullptr),
127  facesPairingPtr_(nullptr),
128  triggerRemoval_(-1),
129  triggerAddition_(-1)
130 {
131  checkDefinition();
132 }
133 
134 
135 Foam::layerAdditionRemoval::layerAdditionRemoval
136 (
137  const word& name,
138  const dictionary& dict,
139  const label index,
140  const polyTopoChanger& ptc
141 )
142 :
143  polyMeshModifier(name, index, ptc, dict.get<bool>("active")),
144  faceZoneID_(dict.lookup("faceZoneName"), ptc.mesh().faceZones()),
145  minLayerThickness_(dict.get<scalar>("minLayerThickness")),
146  maxLayerThickness_(dict.get<scalar>("maxLayerThickness")),
147  thicknessFromVolume_(dict.getOrDefault("thicknessFromVolume", true)),
148  oldLayerThickness_(dict.getOrDefault<scalar>("oldLayerThickness", -1)),
149  pointsPairingPtr_(nullptr),
150  facesPairingPtr_(nullptr),
151  triggerRemoval_(-1),
152  triggerAddition_(-1)
153 {
154  checkDefinition();
155 }
156 
157 
158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
159 
161 {
162  // Protect from multiple calculation in the same time-step
163  if (triggerRemoval_ > -1 || triggerAddition_ > -1)
164  {
165  return true;
166  }
167 
168  // Go through all the cells in the master layer and calculate
169  // approximate layer thickness as the ratio of the cell volume and
170  // face area in the face zone.
171  // Layer addition:
172  // When the max thickness exceeds the threshold, trigger refinement.
173  // Layer removal:
174  // When the min thickness falls below the threshold, trigger removal.
175 
176  const polyMesh& mesh = topoChanger().mesh();
177 
178  const faceZone& fz = mesh.faceZones()[faceZoneID_.index()];
179  const labelList& mc = fz.masterCells();
180 
181  const scalarField& V = mesh.cellVolumes();
182  const vectorField& S = mesh.faceAreas();
183 
184  if (min(V) < -VSMALL)
185  {
187  << "negative cell volume. Error in mesh motion before "
188  << "topological change.\n V: " << V
189  << abort(FatalError);
190  }
191 
192  scalar avgDelta = 0;
193  scalar minDelta = GREAT;
194  scalar maxDelta = 0;
195  label nDelta = 0;
196 
197  if (thicknessFromVolume_)
198  {
199  // Thickness calculated from cell volume/face area
200  forAll(fz, facei)
201  {
202  scalar curDelta = V[mc[facei]]/mag(S[fz[facei]]);
203  avgDelta += curDelta;
204  minDelta = min(minDelta, curDelta);
205  maxDelta = max(maxDelta, curDelta);
206  }
207 
208  nDelta = fz.size();
209  }
210  else
211  {
212  // Thickness calculated from edges on layer
213  const Map<label>& zoneMeshPointMap = fz().meshPointMap();
214 
215  // Edges with only one point on zone
216  forAll(mc, facei)
217  {
218  const cell& cFaces = mesh.cells()[mc[facei]];
219  const edgeList cellEdges(cFaces.edges(mesh.faces()));
220 
221  forAll(cellEdges, i)
222  {
223  const edge& e = cellEdges[i];
224 
225  if (zoneMeshPointMap.found(e[0]))
226  {
227  if (!zoneMeshPointMap.found(e[1]))
228  {
229  scalar curDelta = e.mag(mesh.points());
230  avgDelta += curDelta;
231  nDelta++;
232  minDelta = min(minDelta, curDelta);
233  maxDelta = max(maxDelta, curDelta);
234  }
235  }
236  else
237  {
238  if (zoneMeshPointMap.found(e[1]))
239  {
240  scalar curDelta = e.mag(mesh.points());
241  avgDelta += curDelta;
242  nDelta++;
243  minDelta = min(minDelta, curDelta);
244  maxDelta = max(maxDelta, curDelta);
245  }
246  }
247  }
248  }
249  }
250 
251  reduce(minDelta, minOp<scalar>());
252  reduce(maxDelta, maxOp<scalar>());
253  reduce(avgDelta, sumOp<scalar>());
254  reduce(nDelta, sumOp<label>());
255 
256  avgDelta /= nDelta;
257 
258  if (debug)
259  {
260  Pout<< "bool layerAdditionRemoval::changeTopology() const "
261  << " for object " << name() << " : " << nl
262  << "Layer thickness: min: " << minDelta
263  << " max: " << maxDelta << " avg: " << avgDelta
264  << " old thickness: " << oldLayerThickness_ << nl
265  << "Removal threshold: " << minLayerThickness_
266  << " addition threshold: " << maxLayerThickness_ << endl;
267  }
268 
269  bool topologicalChange = false;
270 
271  // If the thickness is decreasing and crosses the min thickness,
272  // trigger removal
273  if (oldLayerThickness_ < 0)
274  {
275  if (debug)
276  {
277  Pout<< "First step. No addition/removal" << endl;
278  }
279 
280  // No topological changes allowed before first mesh motion
281  oldLayerThickness_ = avgDelta;
282 
283  topologicalChange = false;
284  }
285  else if (avgDelta < oldLayerThickness_)
286  {
287  // Layers moving towards removal
288  if (minDelta < minLayerThickness_)
289  {
290  // Check layer pairing
291  if (setLayerPairing())
292  {
293  // A mesh layer detected. Check that collapse is valid
294  if (validCollapse())
295  {
296  // At this point, info about moving the old mesh
297  // in a way to collapse the cells in the removed
298  // layer is available. Not sure what to do with
299  // it.
300 
301  if (debug)
302  {
303  Pout<< "bool layerAdditionRemoval::changeTopology() "
304  << " const for object " << name() << " : "
305  << "Triggering layer removal" << endl;
306  }
307 
308  triggerRemoval_ = mesh.time().timeIndex();
309 
310  // Old thickness looses meaning.
311  // Set it up to indicate layer removal
312  oldLayerThickness_ = GREAT;
313 
314  topologicalChange = true;
315  }
316  else
317  {
318  // No removal, clear addressing
319  clearAddressing();
320  }
321  }
322  }
323  else
324  {
325  oldLayerThickness_ = avgDelta;
326  }
327  }
328  else
329  {
330  // Layers moving towards addition
331  if (maxDelta > maxLayerThickness_)
332  {
333  if (debug)
334  {
335  Pout<< "bool layerAdditionRemoval::changeTopology() const "
336  << " for object " << name() << " : "
337  << "Triggering layer addition" << endl;
338  }
339 
340  triggerAddition_ = mesh.time().timeIndex();
341 
342  // Old thickness looses meaning.
343  // Set it up to indicate layer removal
344  oldLayerThickness_ = 0;
345 
346  topologicalChange = true;
347  }
348  else
349  {
350  oldLayerThickness_ = avgDelta;
351  }
352  }
353 
354  return topologicalChange;
355 }
356 
357 
358 void Foam::layerAdditionRemoval::setRefinement(polyTopoChange& ref) const
359 {
360  // Insert the layer addition/removal instructions
361  // into the topological change
362 
363  if (triggerRemoval_ == topoChanger().mesh().time().timeIndex())
364  {
365  removeCellLayer(ref);
366 
367  // Clear addressing. This also resets the addition/removal data
368  if (debug)
369  {
370  Pout<< "layerAdditionRemoval::setRefinement(polyTopoChange&) "
371  << "for object " << name() << " : "
372  << "Clearing addressing after layer removal" << endl;
373  }
374 
375  triggerRemoval_ = -1;
376  clearAddressing();
377  }
378 
379  if (triggerAddition_ == topoChanger().mesh().time().timeIndex())
380  {
381  addCellLayer(ref);
382 
383  // Clear addressing. This also resets the addition/removal data
384  if (debug)
385  {
386  Pout<< "layerAdditionRemoval::setRefinement(polyTopoChange&) "
387  << "for object " << name() << " : "
388  << "Clearing addressing after layer addition" << endl;
389  }
391  triggerAddition_ = -1;
392  clearAddressing();
393  }
394 }
395 
396 
397 void Foam::layerAdditionRemoval::updateMesh(const mapPolyMesh&)
398 {
399  if (debug)
400  {
401  Pout<< "layerAdditionRemoval::updateMesh(const mapPolyMesh&) "
402  << "for object " << name() << " : "
403  << "Clearing addressing on external request";
404 
405  if (pointsPairingPtr_ || facesPairingPtr_)
406  {
407  Pout<< "Pointers set." << endl;
408  }
409  else
410  {
411  Pout<< "Pointers not set." << endl;
412  }
413  }
414 
415  // Mesh has changed topologically. Update local topological data
416  faceZoneID_.update(topoChanger().mesh().faceZones());
417 
418  clearAddressing();
419 }
420 
421 
422 void Foam::layerAdditionRemoval::setMinLayerThickness(const scalar t) const
423 {
424  if (t < VSMALL || maxLayerThickness_ < t)
425  {
427  << "Incorrect layer thickness definition."
429  }
430 
431  minLayerThickness_ = t;
432 }
433 
434 
435 void Foam::layerAdditionRemoval::setMaxLayerThickness(const scalar t) const
436 {
437  if (t < minLayerThickness_)
438  {
440  << "Incorrect layer thickness definition."
442  }
443 
444  maxLayerThickness_ = t;
445 }
446 
447 
448 void Foam::layerAdditionRemoval::write(Ostream& os) const
449 {
450  os << nl << type() << nl
451  << name()<< nl
452  << faceZoneID_ << nl
453  << minLayerThickness_ << nl
454  << oldLayerThickness_ << nl
455  << maxLayerThickness_ << nl
456  << thicknessFromVolume_ << endl;
457 }
458 
459 
461 {
462  os << nl;
463 
464  os.beginBlock(name());
465 
466  os.writeEntry("type", type());
467  os.writeEntry("faceZoneName", faceZoneID_.name());
468  os.writeEntry("minLayerThickness", minLayerThickness_);
469  os.writeEntry("maxLayerThickness", maxLayerThickness_);
470  os.writeEntry("thicknessFromVolume", thicknessFromVolume_);
471  os.writeEntry("oldLayerThickness", oldLayerThickness_);
472  os.writeEntry("active", active());
473 
474  os.endBlock();
475 }
476 
477 
478 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
const wordRe & name() const noexcept
The selector name.
Definition: DynamicID.H:123
virtual void write(Ostream &) const
Write.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
List< edge > edgeList
List of edge.
Definition: edgeList.H:32
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual bool changeTopology() const
Check for topology change.
virtual void setRefinement(polyTopoChange &) const
Insert the layer addition/removal instructions.
virtual void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
const cellList & cells() const
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Lookup type of boundary radiation properties.
Definition: lookup.H:57
Macros for easy insertion into run-time selection tables.
void setMaxLayerThickness(const scalar t) const
Set max layer thickness which triggers removal.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
bool active() const noexcept
Has the zone been found.
Definition: DynamicID.H:147
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
void setMinLayerThickness(const scalar t) const
Set min layer thickness which triggers removal.
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
List of mesh modifiers defining the mesh dynamics.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:108
label timeIndex() const noexcept
Return the current time index.
Definition: TimeStateI.H:43
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
Virtual base class for mesh modifiers.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
label index() const
The index of the first matching items, -1 if no matches.
Definition: DynamicID.H:139
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
rDeltaT ref()
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:670
defineTypeNameAndDebug(combustionModel, 0)
const polyTopoChanger & topoChanger() const
Return reference to morph engine.
virtual void writeDict(Ostream &) const
Write dictionary.
const vectorField & faceAreas() const
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
const word & name() const
Return name of this modifier.
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:90
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Namespace for OpenFOAM.
label timeIndex
Definition: getTimeIndex.H:24
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
const scalarField & cellVolumes() const