addCellLayer.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 \*---------------------------------------------------------------------------*/
27 
28 #include "layerAdditionRemoval.H"
29 #include "polyMesh.H"
30 #include "primitiveMesh.H"
31 #include "polyTopoChange.H"
32 #include "polyTopoChanger.H"
33 #include "polyAddPoint.H"
34 #include "polyAddCell.H"
35 #include "polyAddFace.H"
36 #include "polyModifyFace.H"
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
40 Foam::tmp<Foam::vectorField> Foam::layerAdditionRemoval::extrusionDir() const
41 {
42  const polyMesh& mesh = topoChanger().mesh();
43  const primitiveFacePatch& masterFaceLayer =
44  mesh.faceZones()[faceZoneID_.index()]();
45 
46  const pointField& points = mesh.points();
47  const labelList& mp = masterFaceLayer.meshPoints();
48 
49  tmp<vectorField> textrusionDir(new vectorField(mp.size()));
50  vectorField& extrusionDir = textrusionDir.ref();
51 
52  if (setLayerPairing())
53  {
54  if (debug)
55  {
56  Pout<< "void layerAdditionRemoval::extrusionDir() const "
57  << " for object " << name() << " : "
58  << "Using edges for point insertion" << endl;
59  }
60 
61  // Detected a valid layer. Grab the point and face collapse mapping
62  const labelList& ptc = pointsPairing();
63 
64  forAll(extrusionDir, mpI)
65  {
66  extrusionDir[mpI] = points[ptc[mpI]] - points[mp[mpI]];
67  }
68  }
69  else
70  {
71  if (debug)
72  {
73  Pout<< "void layerAdditionRemoval::extrusionDir() const "
74  << " for object " << name() << " : "
75  << "A valid layer could not be found in front of "
76  << "the addition face layer. Using face-based "
77  << "point normals for point addition"
78  << endl;
79  }
80 
81  extrusionDir = minLayerThickness_*masterFaceLayer.pointNormals();
82  }
83 
84  return textrusionDir;
85 }
86 
87 
88 void Foam::layerAdditionRemoval::addCellLayer
89 (
90  polyTopoChange& ref
91 ) const
92 {
93  // Insert the layer addition instructions into the topological change
94 
95  // Algorithm:
96  // 1. For every point in the master zone add a new point
97  // 2. For every face in the master zone add a cell
98  // 3. Go through all the edges of the master zone. For all internal edges,
99  // add a face with the correct orientation and make it internal.
100  // For all boundary edges, find the patch it belongs to and add the face
101  // to this patch
102  // 4. Create a set of new faces on the patch image and assign them to be
103  // between the old master cells and the newly created cells
104  // 5. Modify all the faces in the patch such that they are located
105  // between the old slave cells and newly created cells
106 
107  if (debug)
108  {
109  Pout<< "void layerAdditionRemoval::addCellLayer("
110  << "polyTopoChange& ref) const for object " << name() << " : "
111  << "Adding cell layer" << endl;
112  }
113 
114  // Create the points
115 
116  const polyMesh& mesh = topoChanger().mesh();
117  const primitiveFacePatch& masterFaceLayer =
118  mesh.faceZones()[faceZoneID_.index()]();
119 
120  const pointField& points = mesh.points();
121  const labelList& mp = masterFaceLayer.meshPoints();
122 
123  // Get the extrusion direction for the added points
124 
125  const vectorField pointOffsets(extrusionDir());
126 
127  // Add the new points
128  labelList addedPoints(mp.size());
129 
130  forAll(mp, pointi)
131  {
132  // Add the point nominal thickness away from the master zone point
133  // and grab the label
134  addedPoints[pointi] =
135  ref.setAction
136  (
137  polyAddPoint
138  (
139  points[mp[pointi]] // point
140  + addDelta_*pointOffsets[pointi],
141  mp[pointi], // master point
142  -1, // zone for point
143  true // supports a cell
144  )
145  );
146  }
147 
148  if (debug > 1)
149  {
150  Pout<< "mp: " << mp << " addedPoints: " << addedPoints << endl;
151  }
152 
153  // Create the cells
154 
155  const labelList& mc =
156  mesh.faceZones()[faceZoneID_.index()].masterCells();
157  const labelList& sc =
158  mesh.faceZones()[faceZoneID_.index()].slaveCells();
159 
160  const labelList& mf = mesh.faceZones()[faceZoneID_.index()];
161  const boolList& mfFlip = mesh.faceZones()[faceZoneID_.index()].flipMap();
162 
163  const labelList& own = mesh.faceOwner();
164  const labelList& nei = mesh.faceNeighbour();
165 
166  labelList addedCells(mf.size());
167 
168  forAll(mf, facei)
169  {
170  label celli = mc[facei];
171  label zoneI = mesh.cellZones().whichZone(celli);
172 
173  addedCells[facei] =
174  ref.setAction
175  (
176  polyAddCell
177  (
178  -1, // master point
179  -1, // master edge
180  mf[facei], // master face
181  -1, // master cell
182  zoneI // zone for cell
183  )
184  );
185  }
186 
187  // Create the new faces
188 
189  // Grab the local faces from the master zone
190  const faceList& zoneFaces = masterFaceLayer.localFaces();
191 
192  // Create the new faces by copying the master zone. All the added
193  // faces need to point into the newly created cells, which means
194  // all the faces from the master layer are flipped. The flux flip
195  // is determined from the original master layer face and the face
196  // owner: if the master cell is equal to the face owner the flux
197  // remains the same; otherwise it is flipped
198 
199  forAll(zoneFaces, facei)
200  {
201  const face oldFace = zoneFaces[facei].reverseFace();
202 
203  face newFace(oldFace.size());
204 
205  forAll(oldFace, pointi)
206  {
207  newFace[pointi] = addedPoints[oldFace[pointi]];
208  }
209 
210  bool flipFaceFlux = false;
211 
212  // Flip the face as necessary
213  if
214  (
215  !mesh.isInternalFace(mf[facei])
216  || mc[facei] == nei[mf[facei]]
217  )
218  {
219  flipFaceFlux = true;
220  }
221 
222  // Add the face
223  ref.setAction
224  (
225  polyAddFace
226  (
227  newFace, // face
228  mc[facei], // owner
229  addedCells[facei], // neighbour
230  -1, // master point
231  -1, // master edge
232  mf[facei], // master face for addition
233  flipFaceFlux, // flux flip
234  -1, // patch for face
235  -1, // zone for face
236  false // face zone flip
237  )
238  );
239 
240  if (debug > 1)
241  {
242  Pout<< "adding face: " << newFace
243  << " own: " << mc[facei]
244  << " nei: " << addedCells[facei]
245  << endl;
246  }
247  }
248 
249  // Modify the faces from the master zone for the new neighbour
250 
251  const faceList& faces = mesh.faces();
252 
253  // Pout<< "mfFlip: " << mfFlip << endl;
254 
255  forAll(mf, facei)
256  {
257  const label curfaceID = mf[facei];
258 
259  // If the face is internal, modify its owner to be the newly
260  // created cell. No flip is necessary
261  if (!mesh.isInternalFace(curfaceID))
262  {
263  ref.setAction
264  (
265  polyModifyFace
266  (
267  faces[curfaceID], // modified face
268  curfaceID, // label of face being modified
269  addedCells[facei], // owner
270  -1, // neighbour
271  false, // face flip
272  mesh.boundaryMesh().whichPatch(curfaceID),// patch for face
273  false, // remove from zone
274  faceZoneID_.index(), // zone for face
275  mfFlip[facei] // face flip in zone
276  )
277  );
278 
279  if (debug > 1)
280  {
281  Pout<< "Modifying a boundary face. Face: " << curfaceID
282  << " flip: " << mfFlip[facei]
283  << endl;
284  }
285  }
286 
287  // If slave cell is owner, the face remains the same (but with
288  // a new neighbour - the newly created cell). Otherwise, the
289  // face is flipped.
290  else if (sc[facei] == own[curfaceID])
291  {
292  // Orientation is good, just change neighbour
293  ref.setAction
294  (
295  polyModifyFace
296  (
297  faces[curfaceID], // modified face
298  curfaceID, // label of face being modified
299  own[curfaceID], // owner
300  addedCells[facei], // neighbour
301  false, // face flip
302  mesh.boundaryMesh().whichPatch(curfaceID),// patch for face
303  false, // remove from zone
304  faceZoneID_.index(), // zone for face
305  mfFlip[facei] // face flip in zone
306  )
307  );
308 
309  if (debug > 1)
310  {
311  Pout<< "modify face, no flip " << curfaceID
312  << " own: " << own[curfaceID]
313  << " nei: " << addedCells[facei]
314  << endl;
315  }
316  }
317  else
318  {
319  // Change in orientation; flip face
320  ref.setAction
321  (
322  polyModifyFace
323  (
324  faces[curfaceID].reverseFace(), // modified face
325  curfaceID, // label of face being modified
326  nei[curfaceID], // owner
327  addedCells[facei], // neighbour
328  true, // face flip
329  mesh.boundaryMesh().whichPatch(curfaceID), // patch for face
330  false, // remove from zone
331  faceZoneID_.index(), // zone for face
332  !mfFlip[facei] // face flip in zone
333  )
334  );
335 
336  if (debug > 1)
337  {
338  Pout<< "modify face, with flip " << curfaceID
339  << " own: " << own[curfaceID]
340  << " nei: " << addedCells[facei]
341  << endl;
342  }
343  }
344  }
345 
346  // Create new faces from edges
347 
348  const edgeList& zoneLocalEdges = masterFaceLayer.edges();
349 
350  const labelListList& edgeFaces = masterFaceLayer.edgeFaces();
351 
352  label nInternalEdges = masterFaceLayer.nInternalEdges();
353 
354  const labelList& meshEdges =
355  mesh.faceZones()[faceZoneID_.index()].meshEdges();
356 
357  // Do all internal edges
358  for (label curEdgeID = 0; curEdgeID < nInternalEdges; curEdgeID++)
359  {
360  face newFace(4);
361 
362  newFace[0] = mp[zoneLocalEdges[curEdgeID].start()];
363  newFace[1] = mp[zoneLocalEdges[curEdgeID].end()];
364  newFace[2] = addedPoints[zoneLocalEdges[curEdgeID].end()];
365  newFace[3] = addedPoints[zoneLocalEdges[curEdgeID].start()];
366 
367  ref.setAction
368  (
369  polyAddFace
370  (
371  newFace, // face
372  addedCells[edgeFaces[curEdgeID][0]], // owner
373  addedCells[edgeFaces[curEdgeID][1]], // neighbour
374  -1, // master point
375  meshEdges[curEdgeID], // master edge
376  -1, // master face
377  false, // flip flux
378  -1, // patch for face
379  -1, // zone for face
380  false // face zone flip
381  )
382  );
383 
384  if (debug > 1)
385  {
386  Pout<< "Add internal face off edge: " << newFace
387  << " own: " << addedCells[edgeFaces[curEdgeID][0]]
388  << " nei: " << addedCells[edgeFaces[curEdgeID][1]]
389  << endl;
390  }
391  }
392 
393  // Prepare creation of faces from boundary edges.
394  // Note:
395  // A tricky part of the algorithm is finding the patch into which the
396  // newly created face will be added. For this purpose, take the edge
397  // and grab all the faces coming from it. From the set of faces
398  // eliminate the internal faces and find the boundary face from the rest.
399  // If there is more than one boundary face (excluding the ones in
400  // the master zone), the patch with the lower label is selected.
401 
402  const labelListList& meshEdgeFaces = mesh.edgeFaces();
403 
404  const faceZoneMesh& zoneMesh = mesh.faceZones();
405 
406  // Do all boundary edges
407 
408  for
409  (
410  label curEdgeID = nInternalEdges;
411  curEdgeID < zoneLocalEdges.size();
412  curEdgeID++
413  )
414  {
415  face newFace(4);
416  newFace[0] = mp[zoneLocalEdges[curEdgeID].start()];
417  newFace[1] = mp[zoneLocalEdges[curEdgeID].end()];
418  newFace[2] = addedPoints[zoneLocalEdges[curEdgeID].end()];
419  newFace[3] = addedPoints[zoneLocalEdges[curEdgeID].start()];
420 
421  // Determine the patch for insertion
422  const labelList& curFaces = meshEdgeFaces[meshEdges[curEdgeID]];
423 
424  label patchID = -1;
425  label zoneID = -1;
426 
427  forAll(curFaces, facei)
428  {
429  const label cf = curFaces[facei];
430 
431  if (!mesh.isInternalFace(cf))
432  {
433  // Face not internal. Check if it is in the zone
434  if (zoneMesh.whichZone(cf) != faceZoneID_.index())
435  {
436  // Found the face in a boundary patch which is not in zone
438  zoneID = mesh.faceZones().whichZone(cf);
439 
440  break;
441  }
442  }
443  }
444 
445  if (patchID < 0)
446  {
448  << "Cannot find patch for edge " << meshEdges[curEdgeID]
449  << ". Edge: " << mesh.edges()[meshEdges[curEdgeID]]
450  << abort(FatalError);
451  }
452 
453  ref.setAction
454  (
455  polyAddFace
456  (
457  newFace, // face
458  addedCells[edgeFaces[curEdgeID][0]], // owner
459  -1, // neighbour
460  -1, // master point
461  meshEdges[curEdgeID], // master edge
462  -1, // master face
463  false, // flip flux
464  patchID, // patch for face
465  zoneID, // zone
466  false // zone face flip
467  )
468  );
469 
470  if (debug > 1)
471  {
472  Pout<< "add boundary face: " << newFace
473  << " into patch " << patchID
474  << " own: " << addedCells[edgeFaces[curEdgeID][0]]
475  << endl;
476  }
477  }
478 
479  // Modify the remaining faces of the master cells to reconnect to the new
480  // layer of faces.
481  // Algorithm: Go through all the cells of the master zone and make
482  // a map of faces to avoid duplicates. Do not insert the faces in
483  // the master patch (as they have already been dealt with). Make
484  // a master layer point renumbering map, which for every point in
485  // the master layer gives its new label. Loop through all faces in
486  // the map and attempt to renumber them using the master layer
487  // point renumbering map. Once the face is renumbered, compare it
488  // with the original face; if they are the same, the face has not
489  // changed; if not, modify the face but keep all of its old
490  // attributes (apart from the vertex numbers).
491 
492  // Create the map of faces in the master cell layer
493  labelHashSet masterCellFaceMap(primitiveMesh::facesPerCell_*mc.size());
494 
495  const cellList& cells = mesh.cells();
496 
497  forAll(mc, celli)
498  {
499  const labelList& curFaces = cells[mc[celli]];
500 
501  forAll(curFaces, facei)
502  {
503  // Check if the face belongs to the master zone; if not add it
504  if (zoneMesh.whichZone(curFaces[facei]) != faceZoneID_.index())
505  {
506  masterCellFaceMap.insert(curFaces[facei]);
507  }
508  }
509  }
510 
511  // Create the master layer point map
512  Map<label> masterLayerPointMap(2*mp.size());
513 
514  forAll(mp, pointi)
515  {
516  masterLayerPointMap.insert
517  (
518  mp[pointi],
519  addedPoints[pointi]
520  );
521  }
522 
523  // Grab the list of faces of the master layer
524  const labelList masterCellFaces = masterCellFaceMap.toc();
525 
526  forAll(masterCellFaces, facei)
527  {
528  // Attempt to renumber the face using the masterLayerPointMap.
529  // Missing point remain the same
530 
531  const label curFaceID = masterCellFaces[facei];
532 
533  const face& oldFace = faces[curFaceID];
534 
535  face newFace(oldFace.size());
536 
537  bool changed = false;
538 
539  forAll(oldFace, pointi)
540  {
541  if (masterLayerPointMap.found(oldFace[pointi]))
542  {
543  changed = true;
544 
545  newFace[pointi] = masterLayerPointMap.find(oldFace[pointi])();
546  }
547  else
548  {
549  newFace[pointi] = oldFace[pointi];
550  }
551  }
552 
553  // If the face has changed, create a modification entry
554  if (changed)
555  {
556  // Get face zone and its flip
557  label modifiedFaceZone = mesh.faceZones().whichZone(curFaceID);
558  bool modifiedFaceZoneFlip = false;
559 
560  if (modifiedFaceZone >= 0)
561  {
562  modifiedFaceZoneFlip =
563  mesh.faceZones()[modifiedFaceZone].flipMap()
564  [
565  mesh.faceZones()[modifiedFaceZone].whichFace(curFaceID)
566  ];
567  }
568 
569  if (mesh.isInternalFace(curFaceID))
570  {
571  ref.setAction
572  (
573  polyModifyFace
574  (
575  newFace, // modified face
576  curFaceID, // label of face being modified
577  own[curFaceID], // owner
578  nei[curFaceID], // neighbour
579  false, // face flip
580  -1, // patch for face
581  false, // remove from zone
582  modifiedFaceZone, // zone for face
583  modifiedFaceZoneFlip // face flip in zone
584  )
585  );
586 
587  if (debug > 1)
588  {
589  Pout<< "modifying stick-out face. Internal Old face: "
590  << oldFace
591  << " new face: " << newFace
592  << " own: " << own[curFaceID]
593  << " nei: " << nei[curFaceID]
594  << endl;
595  }
596  }
597  else
598  {
599  ref.setAction
600  (
601  polyModifyFace
602  (
603  newFace, // modified face
604  curFaceID, // label of face being modified
605  own[curFaceID], // owner
606  -1, // neighbour
607  false, // face flip
608  mesh.boundaryMesh().whichPatch(curFaceID),
609  // patch for face
610  false, // remove from zone
611  modifiedFaceZone, // zone for face
612  modifiedFaceZoneFlip // face flip in zone
613  )
614  );
615 
616  if (debug > 1)
617  {
618  Pout<< "modifying stick-out face. Boundary Old face: "
619  << oldFace
620  << " new face: " << newFace
621  << " own: " << own[curFaceID]
622  << " patch: "
623  << mesh.boundaryMesh().whichPatch(curFaceID)
624  << endl;
625  }
626  }
627  }
628  }
629 
630  if (debug)
631  {
632  Pout<< "void layerAdditionRemoval::addCellLayer(polyTopoChange&) const "
633  << " for object " << name() << ": "
634  << "Finished adding cell layer" << endl;
635  }
636 }
637 
638 
639 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
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::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::constant::atomic::mp
const dimensionedScalar mp
Proton mass.
Foam::edgeList
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
Foam::primitiveFacePatch
PrimitivePatch< List< face >, const pointField & > primitiveFacePatch
A PrimitivePatch with List storage for the faces, const reference for the point field.
Definition: primitiveFacePatch.H:51
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::polyMeshModifier::topoChanger
const polyTopoChanger & topoChanger() const
Return reference to morph engine.
Definition: polyMeshModifier.C:63
polyAddFace.H
Foam::primitiveMesh::edgeFaces
const labelListList & edgeFaces() const
Definition: primitiveMeshEdgeFaces.C:33
polyTopoChanger.H
polyTopoChange.H
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:492
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:69
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
primitiveMesh.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::primitiveMesh::edges
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Definition: primitiveMeshEdges.C:505
ref
rDeltaT ref()
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
polyMesh.H
layerAdditionRemoval.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
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:486
Foam::primitiveMesh::facesPerCell_
static const unsigned facesPerCell_
Estimated number of faces per cell.
Definition: primitiveMesh.H:404
Foam::faceZoneMesh
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
Definition: faceZoneMeshFwd.H:44
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
patchID
label patchID
Definition: boundaryProcessorFaPatchPoints.H:5
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
Foam::cellList
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:47
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:810
polyAddPoint.H
Foam::ZoneMesh::whichZone
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:315
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
zoneID
const labelIOList & zoneID
Definition: interpolatedFaces.H:22
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
polyModifyFace.H
Foam::polyMeshModifier::name
const word & name() const
Return name of this modifier.
Definition: polyMeshModifier.H:150
polyAddCell.H
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::polyTopoChanger::mesh
const polyMesh & mesh() const
Return the mesh reference.
Definition: polyTopoChanger.H:111
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:102
points
const pointField & points
Definition: gmvOutputHeader.H:1
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys and label hasher.
Definition: HashSet.H:409
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113
Foam::DynamicID::index
label index() const
Return index of first matching zone.
Definition: DynamicID.H:110