mirrorFvMesh.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 "mirrorFvMesh.H"
27 #include "Time.H"
28 #include "plane.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 Foam::mirrorFvMesh::mirrorFvMesh(const IOobject& io, const word& dictName)
33 :
34  fvMesh(io),
35  mirrorMeshDict_
36  (
37  IOobject
38  (
39  dictName,
40  time().system(),
41  *this,
42  IOobject::MUST_READ,
43  IOobject::NO_WRITE
44  )
45  )
46 {
47  plane mirrorPlane(mirrorMeshDict_);
48 
49  scalar planeTolerance
50  (
51  readScalar(mirrorMeshDict_.lookup("planeTolerance"))
52  );
53 
54  const pointField& oldPoints = points();
55  const faceList& oldFaces = faces();
56  const cellList& oldCells = cells();
57  const label nOldInternalFaces = nInternalFaces();
58  const polyPatchList& oldPatches = boundaryMesh();
59 
60  // Mirror the points
61  Info<< "Mirroring points. Old points: " << oldPoints.size();
62 
63  pointField newPoints(2*oldPoints.size());
64  label nNewPoints = 0;
65 
66  labelList mirrorPointLookup(oldPoints.size(), -1);
67 
68  // Grab the old points
69  forAll(oldPoints, pointi)
70  {
71  newPoints[nNewPoints] = oldPoints[pointi];
72  nNewPoints++;
73  }
74 
75  forAll(oldPoints, pointi)
76  {
77  scalar alpha =
78  mirrorPlane.normalIntersect
79  (
80  oldPoints[pointi],
81  mirrorPlane.normal()
82  );
83 
84  // Check plane on tolerance
85  if (mag(alpha) > planeTolerance)
86  {
87  // The point gets mirrored
88  newPoints[nNewPoints] =
89  oldPoints[pointi] + 2.0*alpha*mirrorPlane.normal();
90 
91  // remember the point correspondence
92  mirrorPointLookup[pointi] = nNewPoints;
93  nNewPoints++;
94  }
95  else
96  {
97  // The point is on the plane and does not get mirrored
98  // Adjust plane location
99  newPoints[nNewPoints] =
100  oldPoints[pointi] + alpha*mirrorPlane.normal();
101 
102  mirrorPointLookup[pointi] = pointi;
103  }
104  }
105 
106  // Reset the size of the point list
107  Info<< " New points: " << nNewPoints << endl;
108  newPoints.setSize(nNewPoints);
109 
110  // Construct new to old map
111  pointMapPtr_.reset(new labelList(newPoints.size()));
112  labelList& pointMap = pointMapPtr_();
113  // Insert old points
114  forAll(oldPoints, oldPointi)
115  {
116  pointMap[oldPointi] = oldPointi;
117  }
118  forAll(mirrorPointLookup, oldPointi)
119  {
120  pointMap[mirrorPointLookup[oldPointi]] = oldPointi;
121  }
122 
123 
124 
125  Info<< "Mirroring faces. Old faces: " << oldFaces.size();
126 
127  // Algorithm:
128  // During mirroring, the faces that were previously boundary faces
129  // in the mirror plane may become ineternal faces. In order to
130  // deal with the ordering of the faces, the algorithm is split
131  // into two parts. For original faces, the internal faces are
132  // distributed to their owner cells. Once all internal faces are
133  // distributed, the boundary faces are visited and if they are in
134  // the mirror plane they are added to the master cells (the future
135  // boundary faces are not touched). After the first phase, the
136  // internal faces are collected in the cell order and numbering
137  // information is added. Then, the internal faces are mirrored
138  // and the face numbering data is stored for the mirrored section.
139  // Once all the internal faces are mirrored, the boundary faces
140  // are added by mirroring the faces patch by patch.
141 
142  // Distribute internal faces
143  labelListList newCellFaces(oldCells.size());
144 
145  const labelUList& oldOwnerStart = lduAddr().ownerStartAddr();
146 
147  forAll(newCellFaces, celli)
148  {
149  labelList& curFaces = newCellFaces[celli];
150 
151  const label s = oldOwnerStart[celli];
152  const label e = oldOwnerStart[celli + 1];
153 
154  curFaces.setSize(e - s);
155 
156  forAll(curFaces, i)
157  {
158  curFaces[i] = s + i;
159  }
160  }
161 
162  // Distribute boundary faces. Remember the faces that have been inserted
163  // as internal
164  boolListList insertedBouFace(oldPatches.size());
165 
166  forAll(oldPatches, patchi)
167  {
168  const polyPatch& curPatch = oldPatches[patchi];
169 
170  if (curPatch.coupled())
171  {
173  << "Found coupled patch " << curPatch.name() << endl
174  << " Mirroring faces on coupled patches destroys"
175  << " the ordering. This might be fixed by running a dummy"
176  << " createPatch afterwards." << endl;
177  }
178 
179  boolList& curInsBouFace = insertedBouFace[patchi];
180 
181  curInsBouFace.setSize(curPatch.size());
182  curInsBouFace = false;
183 
184  // Get faceCells for face insertion
185  const labelUList& curFaceCells = curPatch.faceCells();
186  const label curStart = curPatch.start();
187 
188  forAll(curPatch, facei)
189  {
190  // Find out if the mirrored face is identical to the
191  // original. If so, the face needs to become internal and
192  // added to its owner cell
193  const face& origFace = curPatch[facei];
194 
195  face mirrorFace(origFace.size());
196  forAll(mirrorFace, pointi)
197  {
198  mirrorFace[pointi] = mirrorPointLookup[origFace[pointi]];
199  }
200 
201  if (origFace == mirrorFace)
202  {
203  // The mirror is identical to current face. This will
204  // become an internal face
205  const label oldSize = newCellFaces[curFaceCells[facei]].size();
206 
207  newCellFaces[curFaceCells[facei]].setSize(oldSize + 1);
208  newCellFaces[curFaceCells[facei]][oldSize] = curStart + facei;
209 
210  curInsBouFace[facei] = true;
211  }
212  }
213  }
214 
215  // Construct the new list of faces. Boundary faces are added
216  // last, cush that each patch is mirrored separately. The
217  // addressing is stored in two separate arrays: first for the
218  // original cells (face order has changed) and then for the
219  // mirrored cells.
220  labelList masterFaceLookup(oldFaces.size(), -1);
221  labelList mirrorFaceLookup(oldFaces.size(), -1);
222 
223  faceList newFaces(2*oldFaces.size());
224  label nNewFaces = 0;
225 
226  // Insert original (internal) faces
227  forAll(newCellFaces, celli)
228  {
229  const labelList& curCellFaces = newCellFaces[celli];
230 
231  forAll(curCellFaces, cfI)
232  {
233  newFaces[nNewFaces] = oldFaces[curCellFaces[cfI]];
234  masterFaceLookup[curCellFaces[cfI]] = nNewFaces;
235 
236  nNewFaces++;
237  }
238  }
239 
240  // Mirror internal faces
241  for (label facei = 0; facei < nOldInternalFaces; facei++)
242  {
243  const face& oldFace = oldFaces[facei];
244  face& nf = newFaces[nNewFaces];
245  nf.setSize(oldFace.size());
246 
247  nf[0] = mirrorPointLookup[oldFace[0]];
248 
249  for (label i = 1; i < oldFace.size(); i++)
250  {
251  nf[i] = mirrorPointLookup[oldFace[oldFace.size() - i]];
252  }
253 
254  mirrorFaceLookup[facei] = nNewFaces;
255  nNewFaces++;
256  }
257 
258  // Mirror boundary faces patch by patch
259 
260 
261  labelList newToOldPatch(boundary().size(), -1);
262  labelList newPatchSizes(boundary().size(), -1);
263  labelList newPatchStarts(boundary().size(), -1);
264  label nNewPatches = 0;
265 
267  {
268  const label curPatchSize = boundaryMesh()[patchi].size();
269  const label curPatchStart = boundaryMesh()[patchi].start();
270  const boolList& curInserted = insertedBouFace[patchi];
271 
272  newPatchStarts[nNewPatches] = nNewFaces;
273 
274  // Master side
275  for (label facei = 0; facei < curPatchSize; facei++)
276  {
277  // Check if the face has already been added. If not, add it and
278  // insert the numbering details.
279  if (!curInserted[facei])
280  {
281  newFaces[nNewFaces] = oldFaces[curPatchStart + facei];
282 
283  masterFaceLookup[curPatchStart + facei] = nNewFaces;
284  nNewFaces++;
285  }
286  }
287 
288  // Mirror side
289  for (label facei = 0; facei < curPatchSize; facei++)
290  {
291  // Check if the face has already been added. If not, add it and
292  // insert the numbering details.
293  if (!curInserted[facei])
294  {
295  const face& oldFace = oldFaces[curPatchStart + facei];
296  face& nf = newFaces[nNewFaces];
297  nf.setSize(oldFace.size());
298 
299  nf[0] = mirrorPointLookup[oldFace[0]];
300 
301  for (label i = 1; i < oldFace.size(); i++)
302  {
303  nf[i] = mirrorPointLookup[oldFace[oldFace.size() - i]];
304  }
305 
306  mirrorFaceLookup[curPatchStart + facei] = nNewFaces;
307  nNewFaces++;
308  }
309  else
310  {
311  // Grab the index of the master face for the mirror side
312  mirrorFaceLookup[curPatchStart + facei] =
313  masterFaceLookup[curPatchStart + facei];
314  }
315  }
316 
317  // If patch exists, grab the name and type of the original patch
318  if (nNewFaces > newPatchStarts[nNewPatches])
319  {
320  newToOldPatch[nNewPatches] = patchi;
321 
322  newPatchSizes[nNewPatches] =
323  nNewFaces - newPatchStarts[nNewPatches];
324 
325  nNewPatches++;
326  }
327  }
328 
329  // Tidy up the lists
330  newFaces.setSize(nNewFaces);
331  Info<< " New faces: " << nNewFaces << endl;
332 
333  newToOldPatch.setSize(nNewPatches);
334  newPatchSizes.setSize(nNewPatches);
335  newPatchStarts.setSize(nNewPatches);
336 
337  Info<< "Mirroring patches. Old patches: " << boundary().size()
338  << " New patches: " << nNewPatches << endl;
339 
340  Info<< "Mirroring cells. Old cells: " << oldCells.size()
341  << " New cells: " << 2*oldCells.size() << endl;
342 
343  cellList newCells(2*oldCells.size());
344  label nNewCells = 0;
345 
346  // Construct new to old cell map
347  cellMapPtr_.reset(new labelList(newCells.size()));
348  labelList& cellMap = cellMapPtr_();
349 
350  // Grab the original cells. Take care of face renumbering.
351  forAll(oldCells, celli)
352  {
353  const cell& oc = oldCells[celli];
354 
355  cell& nc = newCells[nNewCells];
356  nc.setSize(oc.size());
357 
358  forAll(oc, i)
359  {
360  nc[i] = masterFaceLookup[oc[i]];
361  }
362 
363  cellMap[nNewCells] = celli;
364 
365  nNewCells++;
366  }
367 
368  // Mirror the cells
369  forAll(oldCells, celli)
370  {
371  const cell& oc = oldCells[celli];
372 
373  cell& nc = newCells[nNewCells];
374  nc.setSize(oc.size());
375 
376  forAll(oc, i)
377  {
378  nc[i] = mirrorFaceLookup[oc[i]];
379  }
380 
381  cellMap[nNewCells] = celli;
382 
383  nNewCells++;
384  }
385 
386  // Mirror the cell shapes
387  Info<< "Mirroring cell shapes." << endl;
388 
389  Info<< nl << "Creating new mesh" << endl;
390  mirrorMeshPtr_.reset
391  (
392  new fvMesh
393  (
394  io,
395  xferMove(newPoints),
396  xferMove(newFaces),
397  xferMove(newCells)
398  )
399  );
400 
401  fvMesh& pMesh = mirrorMeshPtr_();
402 
403  // Add the boundary patches
404  List<polyPatch*> p(newPatchSizes.size());
405 
406  forAll(p, patchi)
407  {
408  p[patchi] = boundaryMesh()[newToOldPatch[patchi]].clone
409  (
410  pMesh.boundaryMesh(),
411  patchi,
412  newPatchSizes[patchi],
413  newPatchStarts[patchi]
414  ).ptr();
415  }
416 
417  pMesh.addPatches(p);
418 }
419 
420 
421 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
422 
424 {}
425 
426 
427 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
PtrList< polyPatch > polyPatchList
container classes for polyPatch
Definition: polyPatchList.H:45
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:270
#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
List< List< bool > > boolListList
Definition: boolList.H:51
label nInternalFaces() const
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
const cellList & cells() const
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
UList< label > labelUList
Definition: UList.H:64
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1003
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Xfer< T > xferMove(T &)
Construct by transferring the contents of the arg.
List< label > labelList
A List of labels.
Definition: labelList.H:56
faceListList boundary(nPatches)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1028
static const char nl
Definition: Ostream.H:265
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void setSize(const label)
Reset size of List.
Definition: List.C:281
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField & p
List< cell > cellList
list of cells
Definition: cellList.H:42
~mirrorFvMesh()
Destructor.
int system(const std::string &command)
Execute the specified command.
Definition: POSIX.C:1217