cellSplitter.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 "cellSplitter.H"
27 #include "polyMesh.H"
28 #include "polyTopoChange.H"
29 #include "polyAddCell.H"
30 #include "polyAddFace.H"
31 #include "polyAddPoint.H"
32 #include "polyModifyFace.H"
33 #include "mapPolyMesh.H"
34 #include "meshTools.H"
35 
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 defineTypeNameAndDebug(cellSplitter, 0);
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 void Foam::cellSplitter::getFaceInfo
48 (
49  const label facei,
50  label& patchID,
51  label& zoneID,
52  label& zoneFlip
53 ) const
54 {
55  patchID = -1;
56 
57  if (!mesh_.isInternalFace(facei))
58  {
59  patchID = mesh_.boundaryMesh().whichPatch(facei);
60  }
61 
62  zoneID = mesh_.faceZones().whichZone(facei);
63 
64  zoneFlip = false;
65 
66  if (zoneID >= 0)
67  {
68  const faceZone& fZone = mesh_.faceZones()[zoneID];
69 
70  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
71  }
72 }
73 
74 
75 // Find the new owner of facei (since the original cell has been split into
76 // newCells
77 Foam::label Foam::cellSplitter::newOwner
78 (
79  const label facei,
80  const Map<labelList>& cellToCells
81 ) const
82 {
83  label oldOwn = mesh_.faceOwner()[facei];
84 
85  Map<labelList>::const_iterator fnd = cellToCells.find(oldOwn);
86 
87  if (fnd == cellToCells.end())
88  {
89  // Unsplit cell
90  return oldOwn;
91  }
92  else
93  {
94  // Look up index of face in the cells' faces.
95 
96  const labelList& newCells = fnd();
97 
98  const cell& cFaces = mesh_.cells()[oldOwn];
99 
100  return newCells[findIndex(cFaces, facei)];
101  }
102 }
103 
104 
105 Foam::label Foam::cellSplitter::newNeighbour
106 (
107  const label facei,
108  const Map<labelList>& cellToCells
109 ) const
110 {
111  label oldNbr = mesh_.faceNeighbour()[facei];
112 
113  Map<labelList>::const_iterator fnd = cellToCells.find(oldNbr);
114 
115  if (fnd == cellToCells.end())
116  {
117  // Unsplit cell
118  return oldNbr;
119  }
120  else
121  {
122  // Look up index of face in the cells' faces.
123 
124  const labelList& newCells = fnd();
125 
126  const cell& cFaces = mesh_.cells()[oldNbr];
127 
128  return newCells[findIndex(cFaces, facei)];
129  }
130 }
131 
132 
133 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
134 
135 // Construct from components
136 Foam::cellSplitter::cellSplitter(const polyMesh& mesh)
137 :
138  mesh_(mesh),
139  addedPoints_()
140 {}
141 
142 
143 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
144 
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
152 (
153  const Map<point>& cellToMidPoint,
154  polyTopoChange& meshMod
155 )
156 {
157  addedPoints_.clear();
158  addedPoints_.resize(cellToMidPoint.size());
159 
160 
161  //
162  // Introduce cellToMidPoints.
163  //
164 
165  forAllConstIter(Map<point>, cellToMidPoint, iter)
166  {
167  label celli = iter.key();
168 
169  label anchorPoint = mesh_.cellPoints()[celli][0];
170 
171  label addedPointi =
172  meshMod.setAction
173  (
174  polyAddPoint
175  (
176  iter(), // point
177  anchorPoint, // master point
178  -1, // zone for point
179  true // supports a cell
180  )
181  );
182  addedPoints_.insert(celli, addedPointi);
183 
184  // Pout<< "Added point " << addedPointi
185  // << iter() << " in cell " << celli << " with centre "
186  // << mesh_.cellCentres()[celli] << endl;
187  }
188 
189 
190  //
191  // Add cells (first one is modified original cell)
192  //
193 
194  Map<labelList> cellToCells(cellToMidPoint.size());
195 
196  forAllConstIter(Map<point>, cellToMidPoint, iter)
197  {
198  label celli = iter.key();
199 
200  const cell& cFaces = mesh_.cells()[celli];
201 
202  // Cells created for this cell.
203  labelList newCells(cFaces.size());
204 
205  // First pyramid is the original cell
206  newCells[0] = celli;
207 
208  // Add other pyramids
209  for (label i = 1; i < cFaces.size(); i++)
210  {
211  label addedCelli =
212  meshMod.setAction
213  (
214  polyAddCell
215  (
216  -1, // master point
217  -1, // master edge
218  -1, // master face
219  celli, // master cell
220  -1 // zone
221  )
222  );
223 
224  newCells[i] = addedCelli;
225  }
226 
227  cellToCells.insert(celli, newCells);
228 
229  // Pout<< "Split cell " << celli
230  // << " with centre " << mesh_.cellCentres()[celli] << nl
231  // << " faces:" << cFaces << nl
232  // << " into :" << newCells << endl;
233  }
234 
235 
236  //
237  // Introduce internal faces. These go from edges of the cell to the mid
238  // point.
239  //
240 
241  forAllConstIter(Map<point>, cellToMidPoint, iter)
242  {
243  label celli = iter.key();
244 
245  label midPointi = addedPoints_[celli];
246 
247  const cell& cFaces = mesh_.cells()[celli];
248 
249  const labelList& cEdges = mesh_.cellEdges()[celli];
250 
251  forAll(cEdges, i)
252  {
253  label edgeI = cEdges[i];
254  const edge& e = mesh_.edges()[edgeI];
255 
256  // Get the faces on the cell using the edge
257  label face0, face1;
258  meshTools::getEdgeFaces(mesh_, celli, edgeI, face0, face1);
259 
260  // Get the cells on both sides of the face by indexing into cFaces.
261  // (since newly created cells are stored in cFaces order)
262  const labelList& newCells = cellToCells[celli];
263 
264  label cell0 = newCells[findIndex(cFaces, face0)];
265  label cell1 = newCells[findIndex(cFaces, face1)];
266 
267  if (cell0 < cell1)
268  {
269  // Construct face to midpoint that is pointing away from
270  // (pyramid split off from) celli
271 
272  const face& f0 = mesh_.faces()[face0];
273 
274  label index = findIndex(f0, e[0]);
275 
276  bool edgeInFaceOrder = (f0[f0.fcIndex(index)] == e[1]);
277 
278  // Check if celli is the face owner
279 
280  face newF(3);
281  if (edgeInFaceOrder == (mesh_.faceOwner()[face0] == celli))
282  {
283  // edge used in face order.
284  newF[0] = e[1];
285  newF[1] = e[0];
286  newF[2] = midPointi;
287  }
288  else
289  {
290  newF[0] = e[0];
291  newF[1] = e[1];
292  newF[2] = midPointi;
293  }
294 
295  // Now newF points away from cell0
296  meshMod.setAction
297  (
298  polyAddFace
299  (
300  newF, // face
301  cell0, // owner
302  cell1, // neighbour
303  -1, // master point
304  -1, // master edge
305  face0, // master face for addition
306  false, // flux flip
307  -1, // patch for face
308  -1, // zone for face
309  false // face zone flip
310  )
311  );
312  }
313  else
314  {
315  // Construct face to midpoint that is pointing away from
316  // (pyramid split off from) celli
317 
318  const face& f1 = mesh_.faces()[face1];
319 
320  label index = findIndex(f1, e[0]);
321 
322  bool edgeInFaceOrder = (f1[f1.fcIndex(index)] == e[1]);
323 
324  // Check if celli is the face owner
325 
326  face newF(3);
327  if (edgeInFaceOrder == (mesh_.faceOwner()[face1] == celli))
328  {
329  // edge used in face order.
330  newF[0] = e[1];
331  newF[1] = e[0];
332  newF[2] = midPointi;
333  }
334  else
335  {
336  newF[0] = e[0];
337  newF[1] = e[1];
338  newF[2] = midPointi;
339  }
340 
341  // Now newF points away from cell1
342  meshMod.setAction
343  (
344  polyAddFace
345  (
346  newF, // face
347  cell1, // owner
348  cell0, // neighbour
349  -1, // master point
350  -1, // master edge
351  face0, // master face for addition
352  false, // flux flip
353  -1, // patch for face
354  -1, // zone for face
355  false // face zone flip
356  )
357  );
358  }
359  }
360  }
361 
362 
363  //
364  // Update all existing faces for split owner or neighbour.
365  //
366 
367 
368  // Mark off affected face.
369  boolList faceUpToDate(mesh_.nFaces(), true);
370 
371  forAllConstIter(Map<point>, cellToMidPoint, iter)
372  {
373  label celli = iter.key();
374 
375  const cell& cFaces = mesh_.cells()[celli];
376 
377  forAll(cFaces, i)
378  {
379  label facei = cFaces[i];
380 
381  faceUpToDate[facei] = false;
382  }
383  }
384 
385  forAll(faceUpToDate, facei)
386  {
387  if (!faceUpToDate[facei])
388  {
389  const face& f = mesh_.faces()[facei];
390 
391  if (mesh_.isInternalFace(facei))
392  {
393  label newOwn = newOwner(facei, cellToCells);
394  label newNbr = newNeighbour(facei, cellToCells);
395 
396  if (newOwn < newNbr)
397  {
398  meshMod.setAction
399  (
400  polyModifyFace
401  (
402  f,
403  facei,
404  newOwn, // owner
405  newNbr, // neighbour
406  false, // flux flip
407  -1, // patch for face
408  false, // remove from zone
409  -1, // zone for face
410  false // face zone flip
411  )
412  );
413  }
414  else
415  {
416  meshMod.setAction
417  (
418  polyModifyFace
419  (
420  f.reverseFace(),
421  facei,
422  newNbr, // owner
423  newOwn, // neighbour
424  false, // flux flip
425  -1, // patch for face
426  false, // remove from zone
427  -1, // zone for face
428  false // face zone flip
429  )
430  );
431  }
432 
433  }
434  else
435  {
436  label newOwn = newOwner(facei, cellToCells);
437 
438  label patchID, zoneID, zoneFlip;
439  getFaceInfo(facei, patchID, zoneID, zoneFlip);
440 
441  meshMod.setAction
442  (
443  polyModifyFace
444  (
445  mesh_.faces()[facei],
446  facei,
447  newOwn, // owner
448  -1, // neighbour
449  false, // flux flip
450  patchID, // patch for face
451  false, // remove from zone
452  zoneID, // zone for face
453  zoneFlip // face zone flip
454  )
455  );
456  }
457 
458  faceUpToDate[facei] = true;
459  }
460  }
461 }
462 
463 
464 void Foam::cellSplitter::updateMesh(const mapPolyMesh& morphMap)
465 {
466  // Create copy since we're deleting entries. Only if both cell and added
467  // point get mapped do they get inserted.
468  Map<label> newAddedPoints(addedPoints_.size());
469 
470  forAllConstIter(Map<label>, addedPoints_, iter)
471  {
472  label oldCelli = iter.key();
473 
474  label newCelli = morphMap.reverseCellMap()[oldCelli];
475 
476  label oldPointi = iter();
477 
478  label newPointi = morphMap.reversePointMap()[oldPointi];
479 
480  if (newCelli >= 0 && newPointi >= 0)
481  {
482  newAddedPoints.insert(newCelli, newPointi);
483  }
484  }
485 
486  // Copy
487  addedPoints_.transfer(newAddedPoints);
488 }
489 
490 
491 // ************************************************************************* //
~cellSplitter()
Destructor.
#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 setRefinement(const Map< point > &cellToMidPoint, polyTopoChange &meshMod)
Insert mesh changes into meshMod.
HashTable< T, label, Hash< label > >::const_iterator const_iterator
Definition: Map.H:59
cellSplitter(const polyMesh &mesh)
Construct from mesh.
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
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
dynamicFvMesh & mesh
List< label > labelList
A List of labels.
Definition: labelList.H:56
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
defineTypeNameAndDebug(combustionModel, 0)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
void getEdgeFaces(const primitiveMesh &, const label celli, const label edgeI, label &face0, label &face1)
Get faces on cell using edgeI. Throws error if no two found.
Definition: meshTools.C:458
Namespace for OpenFOAM.