decoupleSlidingInterface.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-2016 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 \*---------------------------------------------------------------------------*/
25 
26 #include "slidingInterface.H"
27 #include "polyMesh.H"
28 #include "primitiveMesh.H"
29 #include "polyTopoChange.H"
30 #include "polyTopoChanger.H"
31 #include "polyModifyFace.H"
32 #include "polyModifyPoint.H"
33 
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 
36 void Foam::slidingInterface::decoupleInterface
37 (
38  polyTopoChange& ref
39 ) const
40 {
41  if (debug)
42  {
43  Pout<< "void slidingInterface::decoupleInterface("
44  << "polyTopoChange& ref) const : "
45  << "Decoupling sliding interface " << name() << endl;
46  }
47 
48  if (!attached_)
49  {
50  if (debug)
51  {
52  Pout<< "void slidingInterface::decoupleInterface("
53  << "polyTopoChange& ref) const : "
54  << "Interface already decoupled." << endl;
55  }
56 
57  return;
58  }
59 
60  // Clear previous couple
61  clearCouple(ref);
62 
63  const polyMesh& mesh = topoChanger().mesh();
64  const faceList& faces = mesh.faces();
65  const cellList& cells = mesh.cells();
66 
67  const labelList& own = mesh.faceOwner();
68  const labelList& nei = mesh.faceNeighbour();
69 
70  // Master side
71 
72  const primitiveFacePatch& masterPatch =
73  mesh.faceZones()[masterFaceZoneID_.index()]();
74 
75  const labelList& masterPatchAddr =
76  mesh.faceZones()[masterFaceZoneID_.index()];
77 
78  const boolList& masterPatchFlip =
79  mesh.faceZones()[masterFaceZoneID_.index()].flipMap();
80 
81  const labelList& masterFc = masterFaceCells();
82 
83  // Recover faces in master patch
84 
85  forAll(masterPatchAddr, facei)
86  {
87  // Make a copy of the face and turn it if necessary
88  face newFace = faces[masterPatchAddr[facei]];
89 
90  if (masterPatchFlip[facei])
91  {
92  newFace.flip();
93  }
94 
95  ref.setAction
96  (
97  polyModifyFace
98  (
99  newFace, // new face
100  masterPatchAddr[facei], // master face index
101  masterFc[facei], // owner
102  -1, // neighbour
103  false, // flux flip
104  masterPatchID_.index(), // patch ID
105  false, // remove from zone
106  masterFaceZoneID_.index(), // zone ID
107  false // zone flip. Face corrected
108  )
109  );
110 
111  // Pout<< "Modifying master patch face no "
112  // << masterPatchAddr[facei]
113  // << " face: " << faces[masterPatchAddr[facei]]
114  // << " old owner: " << own[masterPatchAddr[facei]]
115  // << " new owner: " << masterFc[facei]
116  // << endl;
117  }
118 
119  // Slave side
120 
121  const primitiveFacePatch& slavePatch =
122  mesh.faceZones()[slaveFaceZoneID_.index()]();
123 
124  const labelList& slavePatchAddr =
125  mesh.faceZones()[slaveFaceZoneID_.index()];
126 
127  const boolList& slavePatchFlip =
128  mesh.faceZones()[slaveFaceZoneID_.index()].flipMap();
129 
130  const labelList& slaveFc = slaveFaceCells();
131 
132  // Grab retired point mapping
133  const Map<label>& rpm = retiredPointMap();
134 
135  // Recover faces in slave patch
136 
137  forAll(slavePatchAddr, facei)
138  {
139  // Make a copy of face and turn it if necessary
140  face newFace = faces[slavePatchAddr[facei]];
141 
142  if (slavePatchFlip[facei])
143  {
144  newFace.flip();
145  }
146 
147  // Recover retired points on the slave side
148  forAll(newFace, pointi)
149  {
150  Map<label>::const_iterator rpmIter = rpm.find(newFace[pointi]);
151  if (rpmIter != rpm.end())
152  {
153  // Master of retired point; grab its original
154  // Pout<< "Reinstating retired point: " << newFace[pointi]
155  // << " with old: " << rpm.find(newFace[pointi])()
156  // << endl;
157 
158  newFace[pointi] = rpmIter();
159  }
160  }
161 
162  ref.setAction
163  (
164  polyModifyFace
165  (
166  newFace, // new face
167  slavePatchAddr[facei], // master face index
168  slaveFc[facei], // owner
169  -1, // neighbour
170  false, // flux flip
171  slavePatchID_.index(), // patch ID
172  false, // remove from zone
173  slaveFaceZoneID_.index(), // zone ID
174  false // zone flip. Face corrected
175  )
176  );
177  }
178 
179  // Re-create the master stick-out faces
180 
181  // Grab the list of faces in the layer
182  const labelList& masterStickOuts = masterStickOutFaces();
183 
184  forAll(masterStickOuts, facei)
185  {
186  // Renumber the face and remove additional points
187 
188  const label curFaceID = masterStickOuts[facei];
189 
190  const face& oldFace = faces[curFaceID];
191 
192  DynamicList<label> newFaceLabels(oldFace.size());
193 
194  bool changed = false;
195 
196  forAll(oldFace, pointi)
197  {
198  // Check if the point is removed
199  if (ref.pointRemoved(oldFace[pointi]))
200  {
201  // Point removed; skip it
202  changed = true;
203  }
204  else
205  {
206  newFaceLabels.append(oldFace[pointi]);
207  }
208  }
209 
210  if (changed)
211  {
212  if (newFaceLabels.size() < 3)
213  {
215  << "Face " << curFaceID << " reduced to less than "
216  << "3 points. Topological/cutting error." << nl
217  << "Old face: " << oldFace << " new face: " << newFaceLabels
218  << abort(FatalError);
219  }
220 
221  // Get face zone and its flip
222  label modifiedFaceZone = mesh.faceZones().whichZone(curFaceID);
223  bool modifiedFaceZoneFlip = false;
224 
225  if (modifiedFaceZone >= 0)
226  {
227  modifiedFaceZoneFlip =
228  mesh.faceZones()[modifiedFaceZone].flipMap()
229  [
230  mesh.faceZones()[modifiedFaceZone].whichFace(curFaceID)
231  ];
232  }
233 
234  face newFace;
235  newFace.transfer(newFaceLabels);
236 
237  // Pout<< "Modifying master stick-out face " << curFaceID
238  // << " old face: " << oldFace
239  // << " new face: " << newFace
240  // << endl;
241 
242  // Modify the face
243  ref.setAction
244  (
245  polyModifyFace
246  (
247  newFace, // modified face
248  curFaceID, // label of face being modified
249  own[curFaceID], // owner
250  nei[curFaceID], // neighbour
251  false, // face flip
252  mesh.boundaryMesh().whichPatch(curFaceID), // patch for face
253  false, // remove from zone
254  modifiedFaceZone, // zone for face
255  modifiedFaceZoneFlip // face flip in zone
256  )
257  );
258  }
259  }
260 
261  // Re-create the slave stick-out faces
262 
263  labelHashSet slaveLayerCellFaceMap
264  (
265  primitiveMesh::facesPerCell_*(masterPatch.size() + slavePatch.size())
266  );
267 
268  forAll(slaveFc, facei)
269  {
270  const labelList& curFaces = cells[slaveFc[facei]];
271 
272  forAll(curFaces, facei)
273  {
274  // Check if the face belongs to the slave face zone; and
275  // if it has been removed; if not add it
276  if
277  (
278  mesh.faceZones().whichZone(curFaces[facei])
279  != slaveFaceZoneID_.index()
280  && !ref.faceRemoved(curFaces[facei])
281 
282  )
283  {
284  slaveLayerCellFaceMap.insert(curFaces[facei]);
285  }
286  }
287  }
288 
289  // Grab the list of faces in the layer
290  const labelList& slaveStickOuts = slaveStickOutFaces();
291 
292  // Grab master point mapping
293  const Map<label>& masterPm = masterPatch.meshPointMap();
294 
295  forAll(slaveStickOuts, facei)
296  {
297  // Renumber the face and remove additional points
298 
299  const label curFaceID = slaveStickOuts[facei];
300 
301  const face& oldFace = faces[curFaceID];
302 
303  DynamicList<label> newFaceLabels(oldFace.size());
304 
305  bool changed = false;
306 
307  forAll(oldFace, pointi)
308  {
309  // Check if the point is removed or retired
310  if (rpm.found(oldFace[pointi]))
311  {
312  // Master of retired point; grab its original
313  changed = true;
314 
315  // Pout<< "Reinstating retired point: " << oldFace[pointi]
316  // << " with old: " << rpm.find(oldFace[pointi])()
317  // << endl;
318 
319  newFaceLabels.append(rpm.find(oldFace[pointi])());
320  }
321  else if (ref.pointRemoved(oldFace[pointi]))
322  {
323  // Point removed; skip it
324  changed = true;
325  }
326  else if (masterPm.found(oldFace[pointi]))
327  {
328  // Point from master patch only; skip it
329  changed = true;
330  }
331  else
332  {
333  newFaceLabels.append(oldFace[pointi]);
334  }
335  }
336 
337  if (changed)
338  {
339  if (newFaceLabels.size() < 3)
340  {
342  << "Face " << curFaceID << " reduced to less than "
343  << "3 points. Topological/cutting error." << nl
344  << "Old face: " << oldFace << " new face: " << newFaceLabels
345  << abort(FatalError);
346  }
347 
348  // Get face zone and its flip
349  label modifiedFaceZone = mesh.faceZones().whichZone(curFaceID);
350  bool modifiedFaceZoneFlip = false;
351 
352  if (modifiedFaceZone >= 0)
353  {
354  modifiedFaceZoneFlip =
355  mesh.faceZones()[modifiedFaceZone].flipMap()
356  [
357  mesh.faceZones()[modifiedFaceZone].whichFace(curFaceID)
358  ];
359  }
360 
361  face newFace;
362  newFace.transfer(newFaceLabels);
363 
364  // Pout<< "Modifying slave stick-out face " << curFaceID
365  // << " old face: " << oldFace
366  // << " new face: " << newFace
367  // << endl;
368 
369  // Modify the face
370  ref.setAction
371  (
372  polyModifyFace
373  (
374  newFace, // modified face
375  curFaceID, // label of face being modified
376  own[curFaceID], // owner
377  nei[curFaceID], // neighbour
378  false, // face flip
379  mesh.boundaryMesh().whichPatch(curFaceID), // patch for face
380  false, // remove from zone
381  modifiedFaceZone, // zone for face
382  modifiedFaceZoneFlip // face flip in zone
383  )
384  );
385  }
386  }
387 
388  // Bring all slave patch points back to life
389  const pointField& points = mesh.points();
390 
391  const labelList& slaveMeshPoints =
392  mesh.faceZones()[slaveFaceZoneID_.index()]().meshPoints();
393 
394  forAll(slaveMeshPoints, pointi)
395  {
396  ref.setAction
397  (
398  polyModifyPoint
399  (
400  slaveMeshPoints[pointi], // point ID
401  points[slaveMeshPoints[pointi]], // point
402  false, // remove from zone
403  mesh.pointZones().whichZone(slaveMeshPoints[pointi]), // zone
404  true // in a cell
405  )
406  );
407  }
408 
409  // Clear the retired point numbering
410  retiredPointMapPtr_->clear();
411 
412  // Finished decoupling
413  attached_ = false;
414 
415  if (debug)
416  {
417  Pout<< "void slidingInterface::coupleInterface("
418  << "polyTopoChange& ref) const : "
419  << "Finished decoupling sliding interface " << name() << endl;
420  }
421 }
422 
423 
424 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
HashTable< label, label, Hash< label > >::const_iterator const_iterator
Definition: Map.H:59
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
static const unsigned facesPerCell_
Estimated number of faces per cell.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
void clear()
Clear all entries from table.
Definition: HashTable.C:468
List< label > labelList
A List of labels.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:131
label index() const
Return index of first matching zone.
Definition: DynamicID.H:107
static const char nl
Definition: Ostream.H:262
PrimitivePatch< face, List, const pointField & > primitiveFacePatch
Foam::primitiveFacePatch.
const polyTopoChanger & topoChanger() const
Return reference to morph engine.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
const word & name() const
Return name of this modifier.
const polyMesh & mesh() const
Return the mesh reference.
List< cell > cellList
list of cells
Definition: cellList.H:42