attachInterface.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 "attachDetach.H"
27 #include "polyMesh.H"
28 #include "primitiveMesh.H"
29 #include "polyTopoChange.H"
30 #include "polyTopoChanger.H"
31 #include "polyRemovePoint.H"
32 #include "polyRemoveFace.H"
33 #include "polyModifyFace.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 const Foam::scalar Foam::attachDetach::positionDifference_ = 1e-8;
38 
39 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40 
41 void Foam::attachDetach::attachInterface
42 (
43  polyTopoChange& ref
44 ) const
45 {
46  // Algorithm:
47  // 1. Create the reverse patch out of the slave faces.
48  // 2. Go through all the mesh points from the master and slave patch.
49  // If the point labels are different, insert them into the point
50  // renumbering list and remove them from the mesh.
51  // 3. Remove all faces from the slave patch
52  // 4. Modify all the faces from the master patch by making them internal
53  // between the faceCell cells for the two patches. If the master owner
54  // is higher than the slave owner, turn the face around
55  // 5. Get all the faces attached to the slave patch points.
56  // If they have not been removed, renumber them using the
57  // point renumbering list.
58 
59  if (debug)
60  {
61  Pout<< "void attachDetach::attachInterface("
62  << "polyTopoChange& ref) const "
63  << " for object " << name() << " : "
64  << "Attaching interface" << endl;
65  }
66 
67  const polyMesh& mesh = topoChanger().mesh();
68  const faceList& faces = mesh.faces();
69  const labelList& own = mesh.faceOwner();
70  const labelList& nei = mesh.faceNeighbour();
71 
72  const polyPatch& masterPatch = mesh.boundaryMesh()[masterPatchID_.index()];
73  const polyPatch& slavePatch = mesh.boundaryMesh()[slavePatchID_.index()];
74 
75  const label masterPatchStart = masterPatch.start();
76  const label slavePatchStart = slavePatch.start();
77 
78  const labelList& slaveMeshPoints = slavePatch.meshPoints();
79 
80  const Map<label>& removedPointMap = pointMatchMap();
81 
82  const labelList removedPoints = removedPointMap.toc();
83 
84  forAll(removedPoints, pointi)
85  {
86  // Pout<< "Removing point:" << removedPoints[pointi]
87  // << " currently at:" << ref.points()[removedPoints[pointi]]
88  // << endl;
89 
90  ref.setAction(polyRemovePoint(removedPoints[pointi]));
91  }
92 
93 // Pout<< "Points to be mapped: " << removedPoints << endl;
94  // Remove all faces from the slave patch
95  forAll(slavePatch, i)
96  {
97  // Pout<< "Removing face " << i + slavePatchStart
98  // << " with verts:" << ref.faces()[i + slavePatchStart]
99  // << " at:"
100  // << UIndirectList<point>
101  // (
102  // ref.points(),
103  // ref.faces()[i + slavePatchStart]
104  // )
105  // << endl;
106  ref.setAction(polyRemoveFace(i + slavePatchStart));
107  }
108 
109  // Modify the faces from the master patch
110  const labelList& masterFaceCells = masterPatch.faceCells();
111  const labelList& slaveFaceCells = slavePatch.faceCells();
112 
113  const boolList& mfFlip = mesh.faceZones()[faceZoneID_.index()].flipMap();
114 
115  forAll(masterFaceCells, facei)
116  {
117  // If slave neighbour is greater than master, face does not need
118  // turning. Modify it to become internal
119  if (masterFaceCells[facei] < slaveFaceCells[facei])
120  {
121  ref.setAction
122  (
123  polyModifyFace
124  (
125  faces[masterPatchStart + facei], // modified face
126  masterPatchStart + facei, // label of face being modified
127  masterFaceCells[facei], // owner
128  slaveFaceCells[facei], // neighbour
129  false, // face flip
130  -1, // patch for face
131  false, // remove from zone
132  faceZoneID_.index(), // zone for face
133  mfFlip[facei] // face flip in zone
134  )
135  );
136  }
137  else
138  {
139  // Flip required
140  ref.setAction
141  (
142  polyModifyFace
143  (
144  faces[masterPatchStart + facei].reverseFace(), // mod face
145  masterPatchStart + facei, // label of face being modified
146  slaveFaceCells[facei], // owner
147  masterFaceCells[facei], // neighbour
148  true, // face flip
149  -1, // patch for face
150  false, // remove from zone
151  faceZoneID_.index(), // zone for face
152  !mfFlip[facei] // face flip in zone
153  )
154  );
155  }
156  }
157 
158  // Renumber faces affected by point removal
159 // Pout<< "slaveMeshPoints: " << slaveMeshPoints << endl;
160  // Make a map of faces that need to be renumbered
161  labelHashSet facesToModifyMap
162  (
163  slaveMeshPoints.size()*primitiveMesh::facesPerPoint_
164  );
165 
166  const labelListList& pf = mesh.pointFaces();
167 
168  // Grab all the faces off the points in the slave patch. If the face has
169  // not been removed, add it to the map of faces to renumber
170  forAll(slaveMeshPoints, pointi)
171  {
172  const labelList& curFaces = pf[slaveMeshPoints[pointi]];
173 
174  forAll(curFaces, facei)
175  {
176  if (!ref.faceRemoved(curFaces[facei]))
177  {
178  facesToModifyMap.insert(curFaces[facei]);
179  }
180  }
181  }
182 
183  // Grab the faces to be renumbered
184  const labelList ftm = facesToModifyMap.toc();
185 
186  forAll(ftm, facei)
187  {
188  // For every face to modify, copy the face and re-map the vertices.
189  // It is known all the faces will be changed since they hang off
190  // re-mapped vertices
191  label curFaceID = ftm[facei];
192 
193  face newFace(faces[curFaceID]);
194 
195  forAll(newFace, pointi)
196  {
198  removedPointMap.find(newFace[pointi]);
199 
200  if (rpmIter != removedPointMap.end())
201  {
202  // Point mapped. Replace it
203  newFace[pointi] = rpmIter();
204  }
205  }
206 
207  // Pout<< "face label: " << curFaceID
208  // << " old face: " << faces[curFaceID]
209  // << " new face: " << newFace
210  // << endl;
211 
212  // Get face zone and its flip
213  label modifiedFaceZone = mesh.faceZones().whichZone(curFaceID);
214  bool modifiedFaceZoneFlip = false;
215 
216  if (modifiedFaceZone >= 0)
217  {
218  modifiedFaceZoneFlip =
219  mesh.faceZones()[modifiedFaceZone].flipMap()
220  [
221  mesh.faceZones()[modifiedFaceZone].whichFace(curFaceID)
222  ];
223  }
224 
225 
226  label patchID = mesh.boundaryMesh().whichPatch(curFaceID);
227  label neiCell;
228  if (patchID == -1)
229  {
230  neiCell = nei[curFaceID];
231  }
232  else
233  {
234  neiCell = -1;
235  }
236 
237 
238  // Modify the face
239  ref.setAction
240  (
241  polyModifyFace
242  (
243  newFace, // modified face
244  curFaceID, // label of face being modified
245  own[curFaceID], // owner
246  neiCell, // neighbour
247  false, // face flip
248  patchID, // patch for face
249  false, // remove from zone
250  modifiedFaceZone, // zone for face
251  modifiedFaceZoneFlip // face flip in zone
252  )
253  );
254  }
255 
256  if (debug)
257  {
258  Pout<< "void attachDetach::attachInterface("
259  << "polyTopoChange& ref) const "
260  << " for object " << name() << " : "
261  << "Finished attaching interface" << endl;
262  }
263 }
264 
265 
267 (
268  pointField& motionPoints
269 ) const
270 {
271  const Map<label>& removedPointMap = pointMatchMap();
272 
273  const labelList removedPoints = removedPointMap.toc();
274 
275  if (debug)
276  {
277  Pout<< "void attachDetach::modifyMotionPoints("
278  << "pointField& motionPoints) const "
279  << " for object " << name() << " : "
280  << "Adjusting motion points." << endl;
281 
282  // Calculate the difference in motion point positions
283  scalar pointDiff = 0;
284 
285  forAll(removedPoints, pointi)
286  {
287  pointDiff +=
288  mag
289  (
290  motionPoints[removedPoints[pointi]]
291  - motionPoints[removedPointMap.find(removedPoints[pointi])()]
292  );
293  }
294 
295  if (pointDiff > removedPoints.size()*positionDifference_)
296  {
297  Pout<< "Point motion difference = " << pointDiff << endl;
298  }
299  }
300 
301  // Put the slave point on top of the master point
302  forAll(removedPoints, pointi)
303  {
304  motionPoints[removedPoints[pointi]] =
305  motionPoints[removedPointMap.find(removedPoints[pointi])()];
306  }
307 
308 }
309 
310 
311 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:163
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
HashTable< label, label, Hash< label > >::const_iterator const_iterator
Definition: Map.H:59
static const unsigned facesPerPoint_
Estimated number of faces per point.
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:208
volScalarField & e
Elementary charge.
Definition: createFields.H:11
List< label > labelList
A List of labels.
Definition: labelList.H:56
label index() const
Return index of first matching zone.
Definition: DynamicID.H:107
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.
virtual void modifyMotionPoints(pointField &motionPoints) const
Modify motion points to comply with the topological change.
dimensioned< scalar > mag(const dimensioned< Type > &)
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:202
const polyMesh & mesh() const
Return the mesh reference.