oldCyclicPolyPatch.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 "oldCyclicPolyPatch.H"
28 #include "polyBoundaryMesh.H"
29 #include "polyMesh.H"
30 #include "demandDrivenData.H"
31 #include "OFstream.H"
32 #include "patchZones.H"
33 #include "matchPoints.H"
34 #include "Time.H"
35 #include "transformList.H"
36 #include "cyclicPolyPatch.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(oldCyclicPolyPatch, 0);
43 
44  addToRunTimeSelectionTable(polyPatch, oldCyclicPolyPatch, word);
45  addToRunTimeSelectionTable(polyPatch, oldCyclicPolyPatch, dictionary);
46 }
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 Foam::pointField Foam::oldCyclicPolyPatch::calcFaceCentres
51 (
52  const UList<face>& faces,
53  const pointField& points
54 )
55 {
56  pointField ctrs(faces.size());
57 
58  forAll(faces, facei)
59  {
60  ctrs[facei] = faces[facei].centre(points);
61  }
62 
63  return ctrs;
64 }
65 
66 
67 Foam::pointField Foam::oldCyclicPolyPatch::getAnchorPoints
68 (
69  const UList<face>& faces,
70  const pointField& points
71 )
72 {
73  pointField anchors(faces.size());
74 
75  forAll(faces, facei)
76  {
77  anchors[facei] = points[faces[facei][0]];
78  }
79 
80  return anchors;
81 }
82 
83 
84 Foam::label Foam::oldCyclicPolyPatch::findMaxArea
85 (
86  const pointField& points,
87  const faceList& faces
88 )
89 {
90  label maxI = -1;
91  scalar maxAreaSqr = -GREAT;
92 
93  forAll(faces, facei)
94  {
95  scalar areaSqr = magSqr(faces[facei].normal(points));
96 
97  if (areaSqr > maxAreaSqr)
98  {
99  maxAreaSqr = areaSqr;
100  maxI = facei;
101  }
102  }
103  return maxI;
104 }
105 
106 
107 bool Foam::oldCyclicPolyPatch::getGeometricHalves
108 (
109  const primitivePatch& pp,
110  labelList& half0ToPatch,
111  labelList& half1ToPatch
112 ) const
113 {
114  // Get geometric zones of patch by looking at normals.
115  // Method 1: any edge with sharpish angle is edge between two halves.
116  // (this will handle e.g. wedge geometries).
117  // Also two fully disconnected regions will be handled this way.
118  // Method 2: sort faces into two halves based on face normal.
119 
120  // Calculate normals
121  const vectorField& faceNormals = pp.faceNormals();
122 
123  // Find edges with sharp angles.
124  boolList regionEdge(pp.nEdges(), false);
125 
126  const labelListList& edgeFaces = pp.edgeFaces();
127 
128  label nRegionEdges = 0;
129 
130  forAll(edgeFaces, edgeI)
131  {
132  const labelList& eFaces = edgeFaces[edgeI];
133 
134  // Check manifold edges for sharp angle.
135  // (Non-manifold already handled by patchZones)
136  if (eFaces.size() == 2)
137  {
138  if ((faceNormals[eFaces[0]] & faceNormals[eFaces[1]])< featureCos_)
139  {
140  regionEdge[edgeI] = true;
141 
142  nRegionEdges++;
143  }
144  }
145  }
146 
147 
148  // For every face determine zone it is connected to (without crossing
149  // any regionEdge)
150  patchZones ppZones(pp, regionEdge);
151 
152  if (debug)
153  {
154  Pout<< "oldCyclicPolyPatch::getGeometricHalves : "
155  << "Found " << nRegionEdges << " edges on patch " << name()
156  << " where the cos of the angle between two connected faces"
157  << " was less than " << featureCos_ << nl
158  << "Patch divided by these and by single sides edges into "
159  << ppZones.nZones() << " parts." << endl;
160 
161 
162  // Dumping zones to obj files.
163 
164  labelList nZoneFaces(ppZones.nZones());
165 
166  for (label zoneI = 0; zoneI < ppZones.nZones(); zoneI++)
167  {
168  OFstream stream
169  (
170  boundaryMesh().mesh().time().path()
171  /name()+"_zone_"+Foam::name(zoneI)+".obj"
172  );
173  Pout<< "oldCyclicPolyPatch::getGeometricHalves : Writing zone "
174  << zoneI << " face centres to OBJ file " << stream.name()
175  << endl;
176 
177  labelList zoneFaces(findIndices(ppZones, zoneI));
178 
179  forAll(zoneFaces, i)
180  {
181  writeOBJ(stream, pp[zoneFaces[i]].centre(pp.points()));
182  }
183 
184  nZoneFaces[zoneI] = zoneFaces.size();
185  }
186  }
187 
188 
189  if (ppZones.nZones() == 2)
190  {
191  half0ToPatch = findIndices(ppZones, 0);
192  half1ToPatch = findIndices(ppZones, 1);
193  }
194  else
195  {
196  if (debug)
197  {
198  Pout<< "oldCyclicPolyPatch::getGeometricHalves :"
199  << " falling back to face-normal comparison" << endl;
200  }
201  label n0Faces = 0;
202  half0ToPatch.setSize(pp.size());
203 
204  label n1Faces = 0;
205  half1ToPatch.setSize(pp.size());
206 
207  // Compare to face 0 normal.
208  forAll(faceNormals, facei)
209  {
210  if ((faceNormals[facei] & faceNormals[0]) > 0)
211  {
212  half0ToPatch[n0Faces++] = facei;
213  }
214  else
215  {
216  half1ToPatch[n1Faces++] = facei;
217  }
218  }
219  half0ToPatch.setSize(n0Faces);
220  half1ToPatch.setSize(n1Faces);
221 
222  if (debug)
223  {
224  Pout<< "oldCyclicPolyPatch::getGeometricHalves :"
225  << " Number of faces per zone:("
226  << n0Faces << ' ' << n1Faces << ')' << endl;
227  }
228  }
229 
230  if (half0ToPatch.size() != half1ToPatch.size())
231  {
232  fileName casePath(boundaryMesh().mesh().time().path());
233 
234  // Dump halves
235  {
236  fileName nm0(casePath/name()+"_half0_faces.obj");
237  Pout<< "oldCyclicPolyPatch::getGeometricHalves : Writing half0"
238  << " faces to OBJ file " << nm0 << endl;
239  writeOBJ(nm0, UIndirectList<face>(pp, half0ToPatch)(), pp.points());
240 
241  fileName nm1(casePath/name()+"_half1_faces.obj");
242  Pout<< "oldCyclicPolyPatch::getGeometricHalves : Writing half1"
243  << " faces to OBJ file " << nm1 << endl;
244  writeOBJ(nm1, UIndirectList<face>(pp, half1ToPatch)(), pp.points());
245  }
246 
247  // Dump face centres
248  {
249  OFstream str0(casePath/name()+"_half0.obj");
250  Pout<< "oldCyclicPolyPatch::getGeometricHalves : Writing half0"
251  << " face centres to OBJ file " << str0.name() << endl;
252 
253  forAll(half0ToPatch, i)
254  {
255  writeOBJ(str0, pp[half0ToPatch[i]].centre(pp.points()));
256  }
257 
258  OFstream str1(casePath/name()+"_half1.obj");
259  Pout<< "oldCyclicPolyPatch::getGeometricHalves : Writing half1"
260  << " face centres to OBJ file " << str1.name() << endl;
261  forAll(half1ToPatch, i)
262  {
263  writeOBJ(str1, pp[half1ToPatch[i]].centre(pp.points()));
264  }
265  }
266 
268  << "Patch " << name() << " gets decomposed in two zones of"
269  << "inequal size: " << half0ToPatch.size()
270  << " and " << half1ToPatch.size() << endl
271  << "This means that the patch is either not two separate regions"
272  << " or one region where the angle between the different regions"
273  << " is not sufficiently sharp." << endl
274  << "Please adapt the featureCos setting." << endl
275  << "Continuing with incorrect face ordering from now on!" << endl;
276 
277  return false;
278  }
279  else
280  {
281  return true;
282  }
283 }
284 
285 
286 void Foam::oldCyclicPolyPatch::getCentresAndAnchors
287 (
288  const primitivePatch& pp,
289  const faceList& half0Faces,
290  const faceList& half1Faces,
291 
292  pointField& ppPoints,
293  pointField& half0Ctrs,
294  pointField& half1Ctrs,
295  pointField& anchors0,
296  scalarField& tols
297 ) const
298 {
299  // Get geometric data on both halves.
300  half0Ctrs = calcFaceCentres(half0Faces, pp.points());
301  anchors0 = getAnchorPoints(half0Faces, pp.points());
302  half1Ctrs = calcFaceCentres(half1Faces, pp.points());
303 
304  switch (transform())
305  {
306  case ROTATIONAL:
307  {
308  label face0 = getConsistentRotationFace(half0Ctrs);
309  label face1 = getConsistentRotationFace(half1Ctrs);
310 
311  vector n0 = ((half0Ctrs[face0] - rotationCentre_) ^ rotationAxis_);
312  vector n1 = ((half1Ctrs[face1] - rotationCentre_) ^ -rotationAxis_);
313  n0 /= mag(n0) + VSMALL;
314  n1 /= mag(n1) + VSMALL;
315 
316  if (debug)
317  {
318  Pout<< "oldCyclicPolyPatch::getCentresAndAnchors :"
319  << " Specified rotation :"
320  << " n0:" << n0 << " n1:" << n1 << endl;
321  }
322 
323  // Rotation (around origin)
324  const tensor reverseT(rotationTensor(n0, -n1));
325 
326  // Rotation
327  forAll(half0Ctrs, facei)
328  {
329  half0Ctrs[facei] = Foam::transform(reverseT, half0Ctrs[facei]);
330  anchors0[facei] = Foam::transform(reverseT, anchors0[facei]);
331  }
332 
333  ppPoints = Foam::transform(reverseT, pp.points());
334 
335  break;
336  }
337  //- Problem: usually specified translation is not accurate enough
338  //- To get proper match so keep automatic determination over here.
339  //case TRANSLATIONAL:
340  //{
341  // // Transform 0 points.
342  //
343  // if (debug)
344  // {
345  // Pout<< "oldCyclicPolyPatch::getCentresAndAnchors :"
346  // << "Specified translation : " << separationVector_
347  // << endl;
348  // }
349  //
350  // half0Ctrs += separationVector_;
351  // anchors0 += separationVector_;
352  // break;
353  //}
354  default:
355  {
356  // Assumes that cyclic is planar. This is also the initial
357  // condition for patches without faces.
358 
359  // Determine the face with max area on both halves. These
360  // two faces are used to determine the transformation tensors
361  label max0I = findMaxArea(pp.points(), half0Faces);
362  vector n0 = half0Faces[max0I].normal(pp.points());
363  n0 /= mag(n0) + VSMALL;
364 
365  label max1I = findMaxArea(pp.points(), half1Faces);
366  vector n1 = half1Faces[max1I].normal(pp.points());
367  n1 /= mag(n1) + VSMALL;
368 
369  if (mag(n0 & n1) < 1-matchTolerance())
370  {
371  if (debug)
372  {
373  Pout<< "oldCyclicPolyPatch::getCentresAndAnchors :"
374  << " Detected rotation :"
375  << " n0:" << n0 << " n1:" << n1 << endl;
376  }
377 
378  // Rotation (around origin)
379  const tensor reverseT(rotationTensor(n0, -n1));
380 
381  // Rotation
382  forAll(half0Ctrs, facei)
383  {
384  half0Ctrs[facei] = Foam::transform
385  (
386  reverseT,
387  half0Ctrs[facei]
388  );
389  anchors0[facei] = Foam::transform
390  (
391  reverseT,
392  anchors0[facei]
393  );
394  }
395  ppPoints = Foam::transform(reverseT, pp.points());
396  }
397  else
398  {
399  // Parallel translation. Get average of all used points.
400 
401  primitiveFacePatch half0(half0Faces, pp.points());
402  const pointField& half0Pts = half0.localPoints();
403  const point ctr0(sum(half0Pts)/half0Pts.size());
404 
405  primitiveFacePatch half1(half1Faces, pp.points());
406  const pointField& half1Pts = half1.localPoints();
407  const point ctr1(sum(half1Pts)/half1Pts.size());
408 
409  if (debug)
410  {
411  Pout<< "oldCyclicPolyPatch::getCentresAndAnchors :"
412  << " Detected translation :"
413  << " n0:" << n0 << " n1:" << n1
414  << " ctr0:" << ctr0 << " ctr1:" << ctr1 << endl;
415  }
416 
417  half0Ctrs += ctr1 - ctr0;
418  anchors0 += ctr1 - ctr0;
419  ppPoints = pp.points() + ctr1 - ctr0;
420  }
421  break;
422  }
423  }
424 
425 
426  // Calculate typical distance per face
427  tols = matchTolerance()*calcFaceTol(half1Faces, pp.points(), half1Ctrs);
428 }
429 
430 
431 bool Foam::oldCyclicPolyPatch::matchAnchors
432 (
433  const bool report,
434  const primitivePatch& pp,
435  const labelList& half0ToPatch,
436  const pointField& anchors0,
437 
438  const labelList& half1ToPatch,
439  const faceList& half1Faces,
440  const labelList& from1To0,
441 
442  const scalarField& tols,
443 
444  labelList& faceMap,
445  labelList& rotation
446 ) const
447 {
448  // Set faceMap such that half0 faces get first and corresponding half1
449  // faces last.
450 
451  forAll(half0ToPatch, half0Facei)
452  {
453  // Label in original patch
454  label patchFacei = half0ToPatch[half0Facei];
455 
456  faceMap[patchFacei] = half0Facei;
457 
458  // No rotation
459  rotation[patchFacei] = 0;
460  }
461 
462  bool fullMatch = true;
463 
464  forAll(from1To0, half1Facei)
465  {
466  label patchFacei = half1ToPatch[half1Facei];
467 
468  // This face has to match the corresponding one on half0.
469  label half0Facei = from1To0[half1Facei];
470 
471  label newFacei = half0Facei + pp.size()/2;
472 
473  faceMap[patchFacei] = newFacei;
474 
475  // Rotate patchFacei such that its f[0] aligns with that of
476  // the corresponding face
477  // (which after shuffling will be at position half0Facei)
478 
479  const point& wantedAnchor = anchors0[half0Facei];
480 
481  rotation[newFacei] = getRotation
482  (
483  pp.points(),
484  half1Faces[half1Facei],
485  wantedAnchor,
486  tols[half1Facei]
487  );
488 
489  if (rotation[newFacei] == -1)
490  {
491  fullMatch = false;
492 
493  if (report)
494  {
495  const face& f = half1Faces[half1Facei];
497  << "Patch:" << name() << " : "
498  << "Cannot find point on face " << f
499  << " with vertices:"
500  << UIndirectList<point>(pp.points(), f)()
501  << " that matches point " << wantedAnchor
502  << " when matching the halves of cyclic patch " << name()
503  << endl
504  << "Continuing with incorrect face ordering from now on!"
505  << endl;
506  }
507  }
508  }
509  return fullMatch;
510 }
511 
512 
513 Foam::label Foam::oldCyclicPolyPatch::getConsistentRotationFace
514 (
515  const pointField& faceCentres
516 ) const
517 {
518  const scalarField magRadSqr
519  (
520  magSqr((faceCentres - rotationCentre_) ^ rotationAxis_)
521  );
522  scalarField axisLen
523  (
524  (faceCentres - rotationCentre_) & rotationAxis_
525  );
526  axisLen = axisLen - min(axisLen);
527  const scalarField magLenSqr
528  (
529  magRadSqr + axisLen*axisLen
530  );
531 
532  label rotFace = -1;
533  scalar maxMagLenSqr = -GREAT;
534  scalar maxMagRadSqr = -GREAT;
535  forAll(faceCentres, i)
536  {
537  if (magLenSqr[i] >= maxMagLenSqr)
538  {
539  if (magRadSqr[i] > maxMagRadSqr)
540  {
541  rotFace = i;
542  maxMagLenSqr = magLenSqr[i];
543  maxMagRadSqr = magRadSqr[i];
544  }
545  }
546  }
547 
548  if (debug)
549  {
550  Info<< "getConsistentRotationFace(const pointField&)" << nl
551  << " rotFace = " << rotFace << nl
552  << " point = " << faceCentres[rotFace] << endl;
553  }
554 
555  return rotFace;
556 }
557 
558 
559 // * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * * * //
560 
562 (
563  const word& name,
564  const label size,
565  const label start,
566  const label index,
567  const polyBoundaryMesh& bm,
568  const word& patchType,
569  const transformType transform
570 )
571 :
572  coupledPolyPatch(name, size, start, index, bm, patchType, transform),
573  featureCos_(0.9),
574  rotationAxis_(Zero),
575  rotationCentre_(Zero),
576  separationVector_(Zero)
577 {}
578 
579 
581 (
582  const word& name,
583  const dictionary& dict,
584  const label index,
585  const polyBoundaryMesh& bm,
586  const word& patchType
587 )
588 :
589  coupledPolyPatch(name, dict, index, bm, patchType),
590  featureCos_(0.9),
591  rotationAxis_(Zero),
592  rotationCentre_(Zero),
593  separationVector_(Zero)
594 {
595  if (dict.found("neighbourPatch"))
596  {
598  (
599  dict
600  ) << "Found \"neighbourPatch\" entry when reading cyclic patch "
601  << name << endl
602  << "Is this mesh already with split cyclics?" << endl
603  << "If so run a newer version that supports it"
604  << ", if not comment out the \"neighbourPatch\" entry and re-run"
605  << exit(FatalIOError);
606  }
607 
608  dict.readIfPresent("featureCos", featureCos_);
609 
610  switch (transform())
611  {
612  case ROTATIONAL:
613  {
614  dict.lookup("rotationAxis") >> rotationAxis_;
615  dict.lookup("rotationCentre") >> rotationCentre_;
616  break;
617  }
618  case TRANSLATIONAL:
619  {
620  dict.lookup("separationVector") >> separationVector_;
621  break;
622  }
623  default:
624  {
625  // no additional info required
626  }
627  }
628 }
629 
630 
632 (
633  const oldCyclicPolyPatch& pp,
634  const polyBoundaryMesh& bm
635 )
636 :
637  coupledPolyPatch(pp, bm),
638  featureCos_(pp.featureCos_),
639  rotationAxis_(pp.rotationAxis_),
640  rotationCentre_(pp.rotationCentre_),
641  separationVector_(pp.separationVector_)
642 {}
643 
644 
646 (
647  const oldCyclicPolyPatch& pp,
648  const polyBoundaryMesh& bm,
649  const label index,
650  const label newSize,
651  const label newStart
652 )
653 :
654  coupledPolyPatch(pp, bm, index, newSize, newStart),
655  featureCos_(pp.featureCos_),
656  rotationAxis_(pp.rotationAxis_),
657  rotationCentre_(pp.rotationCentre_),
658  separationVector_(pp.separationVector_)
659 {}
660 
661 
662 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
663 
665 {}
666 
667 
668 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
669 
671 {
673 }
674 
675 
677 (
678  const primitivePatch& referPatch,
679  const pointField& thisCtrs,
680  const vectorField& thisAreas,
681  const pointField& thisCc,
682  const pointField& nbrCtrs,
683  const vectorField& nbrAreas,
684  const pointField& nbrCc
685 )
686 {}
687 
688 
690 {}
691 
692 
694 (
695  PstreamBuffers& pBufs,
696  const pointField& p
697 )
698 {
699  polyPatch::initMovePoints(pBufs, p);
700 }
701 
702 
704 (
705  PstreamBuffers& pBufs,
706  const pointField& p
707 )
708 {
709  polyPatch::movePoints(pBufs, p);
710 }
711 
712 
714 {
716 }
717 
718 
720 {
721  polyPatch::updateMesh(pBufs);
722 }
723 
724 
726 (
728  const primitivePatch& pp
729 ) const
730 {}
731 
732 
733 // Return new ordering. Ordering is -faceMap: for every face index
734 // the new face -rotation:for every new face the clockwise shift
735 // of the original face. Return false if nothing changes (faceMap
736 // is identity, rotation is 0)
738 (
740  const primitivePatch& pp,
741  labelList& faceMap,
742  labelList& rotation
743 ) const
744 {
745  faceMap.setSize(pp.size());
746  faceMap = -1;
747 
748  rotation.setSize(pp.size());
749  rotation = 0;
750 
751  if (pp.empty())
752  {
753  // No faces, nothing to change.
754  return false;
755  }
756 
757  if (pp.size()&1)
758  {
760  << "Size of cyclic " << name() << " should be a multiple of 2"
761  << ". It is " << pp.size() << abort(FatalError);
762  }
763 
764  label halfSize = pp.size()/2;
765 
766  // Supplied primitivePatch already with new points.
767  // Cyclics are limited to one transformation tensor
768  // currently anyway (i.e. straight plane) so should not be too big a
769  // problem.
770 
771 
772  // Indices of faces on half0
773  labelList half0ToPatch;
774  // Indices of faces on half1
775  labelList half1ToPatch;
776 
777 
778  // 1. Test if already correctly oriented by starting from trivial ordering.
779  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
780 
781  half0ToPatch = identity(halfSize);
782  half1ToPatch = half0ToPatch + halfSize;
783 
784  // Get faces
785  faceList half0Faces(UIndirectList<face>(pp, half0ToPatch));
786  faceList half1Faces(UIndirectList<face>(pp, half1ToPatch));
787 
788  // Get geometric quantities
789  pointField half0Ctrs, half1Ctrs, anchors0, ppPoints;
790  scalarField tols;
791  getCentresAndAnchors
792  (
793  pp,
794  half0Faces,
795  half1Faces,
796 
797  ppPoints,
798  half0Ctrs,
799  half1Ctrs,
800  anchors0,
801  tols
802  );
803 
804  // Geometric match of face centre vectors
805  labelList from1To0;
806  bool matchedAll = matchPoints
807  (
808  half1Ctrs,
809  half0Ctrs,
810  tols,
811  false,
812  from1To0
813  );
814 
815  if (debug)
816  {
817  Pout<< "oldCyclicPolyPatch::order : test if already ordered:"
818  << matchedAll << endl;
819 
820  // Dump halves
821  fileName nm0("match1_"+name()+"_half0_faces.obj");
822  Pout<< "oldCyclicPolyPatch::order : Writing half0"
823  << " faces to OBJ file " << nm0 << endl;
824  writeOBJ(nm0, half0Faces, ppPoints);
825 
826  fileName nm1("match1_"+name()+"_half1_faces.obj");
827  Pout<< "oldCyclicPolyPatch::order : Writing half1"
828  << " faces to OBJ file " << nm1 << endl;
829  writeOBJ(nm1, half1Faces, pp.points());
830 
831  OFstream ccStr
832  (
833  boundaryMesh().mesh().time().path()
834  /"match1_"+ name() + "_faceCentres.obj"
835  );
836  Pout<< "oldCyclicPolyPatch::order : "
837  << "Dumping currently found cyclic match as lines between"
838  << " corresponding face centres to file " << ccStr.name()
839  << endl;
840 
841  // Recalculate untransformed face centres
842  //pointField rawHalf0Ctrs = calcFaceCentres(half0Faces, pp.points());
843  label vertI = 0;
844 
845  forAll(half1Ctrs, i)
846  {
847  //if (from1To0[i] != -1)
848  {
849  // Write edge between c1 and c0
850  //const point& c0 = rawHalf0Ctrs[from1To0[i]];
851  //const point& c0 = half0Ctrs[from1To0[i]];
852  const point& c0 = half0Ctrs[i];
853  const point& c1 = half1Ctrs[i];
854  writeOBJ(ccStr, c0, c1, vertI);
855  }
856  }
857  }
858 
859 
860  // 2. Ordered in pairs (so 0,1 coupled and 2,3 etc.)
861  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
862 
863  if (!matchedAll)
864  {
865  label facei = 0;
866  for (label i = 0; i < halfSize; i++)
867  {
868  half0ToPatch[i] = facei++;
869  half1ToPatch[i] = facei++;
870  }
871 
872  // And redo all matching
873  half0Faces = UIndirectList<face>(pp, half0ToPatch);
874  half1Faces = UIndirectList<face>(pp, half1ToPatch);
875 
876  getCentresAndAnchors
877  (
878  pp,
879  half0Faces,
880  half1Faces,
881 
882  ppPoints,
883  half0Ctrs,
884  half1Ctrs,
885  anchors0,
886  tols
887  );
888 
889  // Geometric match of face centre vectors
890  matchedAll = matchPoints
891  (
892  half1Ctrs,
893  half0Ctrs,
894  tols,
895  false,
896  from1To0
897  );
898 
899  if (debug)
900  {
901  Pout<< "oldCyclicPolyPatch::order : test if pairwise ordered:"
902  << matchedAll << endl;
903  // Dump halves
904  fileName nm0("match2_"+name()+"_half0_faces.obj");
905  Pout<< "oldCyclicPolyPatch::order : Writing half0"
906  << " faces to OBJ file " << nm0 << endl;
907  writeOBJ(nm0, half0Faces, ppPoints);
908 
909  fileName nm1("match2_"+name()+"_half1_faces.obj");
910  Pout<< "oldCyclicPolyPatch::order : Writing half1"
911  << " faces to OBJ file " << nm1 << endl;
912  writeOBJ(nm1, half1Faces, pp.points());
913 
914  OFstream ccStr
915  (
916  boundaryMesh().mesh().time().path()
917  /"match2_"+name()+"_faceCentres.obj"
918  );
919  Pout<< "oldCyclicPolyPatch::order : "
920  << "Dumping currently found cyclic match as lines between"
921  << " corresponding face centres to file " << ccStr.name()
922  << endl;
923 
924  // Recalculate untransformed face centres
925  label vertI = 0;
926 
927  forAll(half1Ctrs, i)
928  {
929  if (from1To0[i] != -1)
930  {
931  // Write edge between c1 and c0
932  const point& c0 = half0Ctrs[from1To0[i]];
933  const point& c1 = half1Ctrs[i];
934  writeOBJ(ccStr, c0, c1, vertI);
935  }
936  }
937  }
938  }
939 
940 
941  // 3. Baffles(coincident faces) converted into cyclics (e.g. jump)
942  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
943 
944  if (!matchedAll)
945  {
946  label baffleI = 0;
947 
948  forAll(pp, facei)
949  {
950  const face& f = pp.localFaces()[facei];
951  const labelList& pFaces = pp.pointFaces()[f[0]];
952 
953  label matchedFacei = -1;
954 
955  forAll(pFaces, i)
956  {
957  label otherFacei = pFaces[i];
958 
959  if (otherFacei > facei)
960  {
961  const face& otherF = pp.localFaces()[otherFacei];
962 
963  // Note: might pick up two similar oriented faces
964  // (but that is illegal anyway)
965  if (f == otherF)
966  {
967  matchedFacei = otherFacei;
968  break;
969  }
970  }
971  }
972 
973  if (matchedFacei != -1)
974  {
975  half0ToPatch[baffleI] = facei;
976  half1ToPatch[baffleI] = matchedFacei;
977  baffleI++;
978  }
979  }
980 
981  if (baffleI == halfSize)
982  {
983  // And redo all matching
984  half0Faces = UIndirectList<face>(pp, half0ToPatch);
985  half1Faces = UIndirectList<face>(pp, half1ToPatch);
986 
987  getCentresAndAnchors
988  (
989  pp,
990  half0Faces,
991  half1Faces,
992 
993  ppPoints,
994  half0Ctrs,
995  half1Ctrs,
996  anchors0,
997  tols
998  );
999 
1000  // Geometric match of face centre vectors
1001  matchedAll = matchPoints
1002  (
1003  half1Ctrs,
1004  half0Ctrs,
1005  tols,
1006  false,
1007  from1To0
1008  );
1009 
1010  if (debug)
1011  {
1012  Pout<< "oldCyclicPolyPatch::order : test if baffles:"
1013  << matchedAll << endl;
1014  // Dump halves
1015  fileName nm0("match3_"+name()+"_half0_faces.obj");
1016  Pout<< "oldCyclicPolyPatch::order : Writing half0"
1017  << " faces to OBJ file " << nm0 << endl;
1018  writeOBJ(nm0, half0Faces, ppPoints);
1019 
1020  fileName nm1("match3_"+name()+"_half1_faces.obj");
1021  Pout<< "oldCyclicPolyPatch::order : Writing half1"
1022  << " faces to OBJ file " << nm1 << endl;
1023  writeOBJ(nm1, half1Faces, pp.points());
1024 
1025  OFstream ccStr
1026  (
1027  boundaryMesh().mesh().time().path()
1028  /"match3_"+ name() + "_faceCentres.obj"
1029  );
1030  Pout<< "oldCyclicPolyPatch::order : "
1031  << "Dumping currently found cyclic match as lines between"
1032  << " corresponding face centres to file " << ccStr.name()
1033  << endl;
1034 
1035  // Recalculate untransformed face centres
1036  label vertI = 0;
1037 
1038  forAll(half1Ctrs, i)
1039  {
1040  if (from1To0[i] != -1)
1041  {
1042  // Write edge between c1 and c0
1043  const point& c0 = half0Ctrs[from1To0[i]];
1044  const point& c1 = half1Ctrs[i];
1045  writeOBJ(ccStr, c0, c1, vertI);
1046  }
1047  }
1048  }
1049  }
1050  }
1051 
1052 
1053  // 4. Automatic geometric ordering
1054  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1055 
1056  if (!matchedAll)
1057  {
1058  // Split faces according to feature angle or topology
1059  bool okSplit = getGeometricHalves(pp, half0ToPatch, half1ToPatch);
1060 
1061  if (!okSplit)
1062  {
1063  // Did not split into two equal parts.
1064  return false;
1065  }
1066 
1067  // And redo all matching
1068  half0Faces = UIndirectList<face>(pp, half0ToPatch);
1069  half1Faces = UIndirectList<face>(pp, half1ToPatch);
1070 
1071  getCentresAndAnchors
1072  (
1073  pp,
1074  half0Faces,
1075  half1Faces,
1076 
1077  ppPoints,
1078  half0Ctrs,
1079  half1Ctrs,
1080  anchors0,
1081  tols
1082  );
1083 
1084  // Geometric match of face centre vectors
1085  matchedAll = matchPoints
1086  (
1087  half1Ctrs,
1088  half0Ctrs,
1089  tols,
1090  false,
1091  from1To0
1092  );
1093 
1094  if (debug)
1095  {
1096  Pout<< "oldCyclicPolyPatch::order : automatic ordering result:"
1097  << matchedAll << endl;
1098  // Dump halves
1099  fileName nm0("match4_"+name()+"_half0_faces.obj");
1100  Pout<< "oldCyclicPolyPatch::order : Writing half0"
1101  << " faces to OBJ file " << nm0 << endl;
1102  writeOBJ(nm0, half0Faces, ppPoints);
1103 
1104  fileName nm1("match4_"+name()+"_half1_faces.obj");
1105  Pout<< "oldCyclicPolyPatch::order : Writing half1"
1106  << " faces to OBJ file " << nm1 << endl;
1107  writeOBJ(nm1, half1Faces, pp.points());
1108 
1109  OFstream ccStr
1110  (
1111  boundaryMesh().mesh().time().path()
1112  /"match4_"+ name() + "_faceCentres.obj"
1113  );
1114  Pout<< "oldCyclicPolyPatch::order : "
1115  << "Dumping currently found cyclic match as lines between"
1116  << " corresponding face centres to file " << ccStr.name()
1117  << endl;
1118 
1119  // Recalculate untransformed face centres
1120  label vertI = 0;
1121 
1122  forAll(half1Ctrs, i)
1123  {
1124  if (from1To0[i] != -1)
1125  {
1126  // Write edge between c1 and c0
1127  const point& c0 = half0Ctrs[from1To0[i]];
1128  const point& c1 = half1Ctrs[i];
1129  writeOBJ(ccStr, c0, c1, vertI);
1130  }
1131  }
1132  }
1133  }
1134 
1135 
1136  if (!matchedAll || debug)
1137  {
1138  // Dump halves
1139  fileName nm0(name()+"_half0_faces.obj");
1140  Pout<< "oldCyclicPolyPatch::order : Writing half0"
1141  << " faces to OBJ file " << nm0 << endl;
1142  writeOBJ(nm0, half0Faces, pp.points());
1143 
1144  fileName nm1(name()+"_half1_faces.obj");
1145  Pout<< "oldCyclicPolyPatch::order : Writing half1"
1146  << " faces to OBJ file " << nm1 << endl;
1147  writeOBJ(nm1, half1Faces, pp.points());
1148 
1149  OFstream ccStr
1150  (
1151  boundaryMesh().mesh().time().path()
1152  /name() + "_faceCentres.obj"
1153  );
1154  Pout<< "oldCyclicPolyPatch::order : "
1155  << "Dumping currently found cyclic match as lines between"
1156  << " corresponding face centres to file " << ccStr.name()
1157  << endl;
1158 
1159  // Recalculate untransformed face centres
1160  //pointField rawHalf0Ctrs = calcFaceCentres(half0Faces, pp.points());
1161  label vertI = 0;
1162 
1163  forAll(half1Ctrs, i)
1164  {
1165  if (from1To0[i] != -1)
1166  {
1167  // Write edge between c1 and c0
1168  //const point& c0 = rawHalf0Ctrs[from1To0[i]];
1169  const point& c0 = half0Ctrs[from1To0[i]];
1170  const point& c1 = half1Ctrs[i];
1171  writeOBJ(ccStr, c0, c1, vertI);
1172  }
1173  }
1174  }
1175 
1176 
1177  if (!matchedAll)
1178  {
1180  << "Patch:" << name() << " : "
1181  << "Cannot match vectors to faces on both sides of patch" << endl
1182  << " Perhaps your faces do not match?"
1183  << " The obj files written contain the current match." << endl
1184  << " Continuing with incorrect face ordering from now on!"
1185  << endl;
1186 
1187  return false;
1188  }
1189 
1190 
1191  // Set faceMap such that half0 faces get first and corresponding half1
1192  // faces last.
1193  matchAnchors
1194  (
1195  true, // report if anchor matching error
1196  pp,
1197  half0ToPatch,
1198  anchors0,
1199  half1ToPatch,
1200  half1Faces,
1201  from1To0,
1202  tols,
1203  faceMap,
1204  rotation
1205  );
1206 
1207  // Return false if no change neccesary, true otherwise.
1208 
1209  forAll(faceMap, facei)
1210  {
1211  if (faceMap[facei] != facei || rotation[facei] != 0)
1212  {
1213  return true;
1214  }
1215  }
1216 
1217  return false;
1218 }
1219 
1220 
1222 {
1223  // Replacement of polyPatch::write to write 'cyclic' instead of type():
1224  os.writeKeyword("type") << cyclicPolyPatch::typeName
1225  << token::END_STATEMENT << nl;
1227  os.writeKeyword("nFaces") << size() << token::END_STATEMENT << nl;
1228  os.writeKeyword("startFace") << start() << token::END_STATEMENT << nl;
1229 
1230 
1231  os.writeKeyword("featureCos") << featureCos_ << token::END_STATEMENT << nl;
1232  switch (transform())
1233  {
1234  case ROTATIONAL:
1235  {
1236  os.writeKeyword("rotationAxis") << rotationAxis_
1237  << token::END_STATEMENT << nl;
1238  os.writeKeyword("rotationCentre") << rotationCentre_
1239  << token::END_STATEMENT << nl;
1240  break;
1241  }
1242  case TRANSLATIONAL:
1243  {
1244  os.writeKeyword("separationVector") << separationVector_
1245  << token::END_STATEMENT << nl;
1246  break;
1247  }
1248  default:
1249  {
1250  // no additional info to write
1251  }
1252  }
1253 
1255  << "Please run foamUpgradeCyclics to convert these old-style"
1256  << " cyclics into two separate cyclics patches."
1257  << endl;
1258 }
1259 
1260 
1261 // ************************************************************************* //
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
#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
A class for handling file names.
Definition: fileName.H:69
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: polyPatch.H:115
Output to file stream.
Definition: OFstream.H:82
List< face > faceList
Definition: faceListFwd.H:43
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
&#39;old&#39; style cyclic polyPatch with all faces in single patch. Does ordering but cannot be used to run...
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:47
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
Definition: polyPatch.H:108
labelList findIndices(const ListType &, typename ListType::const_reference, const label start=0)
Find all occurences of given element. Linear search.
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
Determine correspondence between points. See below.
Macros for easy insertion into run-time selection tables.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
virtual void movePoints(PstreamBuffers &, const pointField &p)
Correct patches after moving points.
Definition: polyPatch.C:57
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
virtual void movePoints(PstreamBuffers &, const pointField &)
Correct patches after moving points.
A list of faces which address into the list of points.
bool matchPoints(const UList< point > &pts0, const UList< point > &pts1, const UList< scalar > &matchDistance, const bool verbose, labelList &from0To1, const point &origin=point::zero)
Determine correspondence between pointFields. Gets passed.
Definition: matchPoints.C:33
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
const Field< PointType > & points() const
Return reference to global points.
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:91
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: polyPatch.C:62
errorManip< error > abort(error &err)
Definition: errorManip.H:131
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
Foam::polyBoundaryMesh.
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
const labelListList & pointFaces() const
Return point-face addressing.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:262
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
defineTypeNameAndDebug(combustionModel, 0)
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
labelList f(nPoints)
Buffers for inter-processor communications streams (UOPstream, UIPstream).
PrimitivePatch< face, List, const pointField & > primitiveFacePatch
Foam::primitiveFacePatch.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
PrimitivePatch< face, SubList, const pointField & > primitivePatch
Addressing for a faceList slice.
Spatial transformation functions for primitive fields.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: polyPatch.H:100
void write(Ostream &) const
Write patchIdentifier as a dictionary.
void setSize(const label)
Reset size of List.
Definition: List.C:281
Template functions to aid in the implementation of demand driven data.
vector point
Point is a vector.
Definition: point.H:41
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
A List with indirect addressing.
Definition: fvMatrix.H:106
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
oldCyclicPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform=UNKNOWN)
Construct from components.
IOerror FatalIOError