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