oldCyclicPolyPatch.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 "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].area(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  const label max0I = findMaxArea(pp.points(), half0Faces);
362  const vector n0 = half0Faces[max0I].normal(pp.points());
363 
364  const label max1I = findMaxArea(pp.points(), half1Faces);
365  const vector n1 = half1Faces[max1I].normal(pp.points());
366 
367  if (mag(n0 & n1) < 1-matchTolerance())
368  {
369  if (debug)
370  {
371  Pout<< "oldCyclicPolyPatch::getCentresAndAnchors :"
372  << " Detected rotation :"
373  << " n0:" << n0 << " n1:" << n1 << endl;
374  }
375 
376  // Rotation (around origin)
377  const tensor reverseT(rotationTensor(n0, -n1));
378 
379  // Rotation
380  forAll(half0Ctrs, facei)
381  {
382  half0Ctrs[facei] = Foam::transform
383  (
384  reverseT,
385  half0Ctrs[facei]
386  );
387  anchors0[facei] = Foam::transform
388  (
389  reverseT,
390  anchors0[facei]
391  );
392  }
393  ppPoints = Foam::transform(reverseT, pp.points());
394  }
395  else
396  {
397  // Parallel translation. Get average of all used points.
398 
399  primitiveFacePatch half0(half0Faces, pp.points());
400  const pointField& half0Pts = half0.localPoints();
401  const point ctr0(sum(half0Pts)/half0Pts.size());
402 
403  primitiveFacePatch half1(half1Faces, pp.points());
404  const pointField& half1Pts = half1.localPoints();
405  const point ctr1(sum(half1Pts)/half1Pts.size());
406 
407  if (debug)
408  {
409  Pout<< "oldCyclicPolyPatch::getCentresAndAnchors :"
410  << " Detected translation :"
411  << " n0:" << n0 << " n1:" << n1
412  << " ctr0:" << ctr0 << " ctr1:" << ctr1 << endl;
413  }
414 
415  half0Ctrs += ctr1 - ctr0;
416  anchors0 += ctr1 - ctr0;
417  ppPoints = pp.points() + ctr1 - ctr0;
418  }
419  break;
420  }
421  }
422 
423 
424  // Calculate typical distance per face
425  tols = matchTolerance()*calcFaceTol(half1Faces, pp.points(), half1Ctrs);
426 }
427 
428 
429 bool Foam::oldCyclicPolyPatch::matchAnchors
430 (
431  const bool report,
432  const primitivePatch& pp,
433  const labelList& half0ToPatch,
434  const pointField& anchors0,
435 
436  const labelList& half1ToPatch,
437  const faceList& half1Faces,
438  const labelList& from1To0,
439 
440  const scalarField& tols,
441 
442  labelList& faceMap,
443  labelList& rotation
444 ) const
445 {
446  // Set faceMap such that half0 faces get first and corresponding half1
447  // faces last.
448 
449  forAll(half0ToPatch, half0Facei)
450  {
451  // Label in original patch
452  label patchFacei = half0ToPatch[half0Facei];
453 
454  faceMap[patchFacei] = half0Facei;
455 
456  // No rotation
457  rotation[patchFacei] = 0;
458  }
459 
460  bool fullMatch = true;
461 
462  forAll(from1To0, half1Facei)
463  {
464  label patchFacei = half1ToPatch[half1Facei];
465 
466  // This face has to match the corresponding one on half0.
467  label half0Facei = from1To0[half1Facei];
468 
469  label newFacei = half0Facei + pp.size()/2;
470 
471  faceMap[patchFacei] = newFacei;
472 
473  // Rotate patchFacei such that its f[0] aligns with that of
474  // the corresponding face
475  // (which after shuffling will be at position half0Facei)
476 
477  const point& wantedAnchor = anchors0[half0Facei];
478 
479  rotation[newFacei] = getRotation
480  (
481  pp.points(),
482  half1Faces[half1Facei],
483  wantedAnchor,
484  tols[half1Facei]
485  );
486 
487  if (rotation[newFacei] == -1)
488  {
489  fullMatch = false;
490 
491  if (report)
492  {
493  const face& f = half1Faces[half1Facei];
495  << "Patch:" << name() << " : "
496  << "Cannot find point on face " << f
497  << " with vertices:"
498  << UIndirectList<point>(pp.points(), f)()
499  << " that matches point " << wantedAnchor
500  << " when matching the halves of cyclic patch " << name()
501  << endl
502  << "Continuing with incorrect face ordering from now on!"
503  << endl;
504  }
505  }
506  }
507  return fullMatch;
508 }
509 
510 
511 Foam::label Foam::oldCyclicPolyPatch::getConsistentRotationFace
512 (
513  const pointField& faceCentres
514 ) const
515 {
516  const scalarField magRadSqr
517  (
518  magSqr((faceCentres - rotationCentre_) ^ rotationAxis_)
519  );
520  scalarField axisLen
521  (
522  (faceCentres - rotationCentre_) & rotationAxis_
523  );
524  axisLen = axisLen - min(axisLen);
525  const scalarField magLenSqr
526  (
527  magRadSqr + axisLen*axisLen
528  );
529 
530  label rotFace = -1;
531  scalar maxMagLenSqr = -great;
532  scalar maxMagRadSqr = -great;
533  forAll(faceCentres, i)
534  {
535  if (magLenSqr[i] >= maxMagLenSqr)
536  {
537  if (magRadSqr[i] > maxMagRadSqr)
538  {
539  rotFace = i;
540  maxMagLenSqr = magLenSqr[i];
541  maxMagRadSqr = magRadSqr[i];
542  }
543  }
544  }
545 
546  if (debug)
547  {
548  Info<< "getConsistentRotationFace(const pointField&)" << nl
549  << " rotFace = " << rotFace << nl
550  << " point = " << faceCentres[rotFace] << endl;
551  }
552 
553  return rotFace;
554 }
555 
556 
557 // * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * * * //
558 
560 (
561  const word& name,
562  const label size,
563  const label start,
564  const label index,
565  const polyBoundaryMesh& bm,
566  const word& patchType,
567  const transformType transform
568 )
569 :
570  coupledPolyPatch(name, size, start, index, bm, patchType, transform),
571  featureCos_(0.9),
572  rotationAxis_(Zero),
573  rotationCentre_(Zero),
574  separationVector_(Zero)
575 {}
576 
577 
579 (
580  const word& name,
581  const dictionary& dict,
582  const label index,
583  const polyBoundaryMesh& bm,
584  const word& patchType
585 )
586 :
587  coupledPolyPatch(name, dict, index, bm, patchType),
588  featureCos_(0.9),
589  rotationAxis_(Zero),
590  rotationCentre_(Zero),
591  separationVector_(Zero)
592 {
593  if (dict.found("neighbourPatch"))
594  {
596  (
597  dict
598  ) << "Found \"neighbourPatch\" entry when reading cyclic patch "
599  << name << endl
600  << "Is this mesh already with split cyclics?" << endl
601  << "If so run a newer version that supports it"
602  << ", if not comment out the \"neighbourPatch\" entry and re-run"
603  << exit(FatalIOError);
604  }
605 
606  dict.readIfPresent("featureCos", featureCos_);
607 
608  switch (transform())
609  {
610  case ROTATIONAL:
611  {
612  dict.lookup("rotationAxis") >> rotationAxis_;
613  dict.lookup("rotationCentre") >> rotationCentre_;
614  break;
615  }
616  case TRANSLATIONAL:
617  {
618  dict.lookup("separationVector") >> separationVector_;
619  break;
620  }
621  default:
622  {
623  // no additional info required
624  }
625  }
626 }
627 
628 
630 (
631  const oldCyclicPolyPatch& pp,
632  const polyBoundaryMesh& bm
633 )
634 :
635  coupledPolyPatch(pp, bm),
636  featureCos_(pp.featureCos_),
637  rotationAxis_(pp.rotationAxis_),
638  rotationCentre_(pp.rotationCentre_),
639  separationVector_(pp.separationVector_)
640 {}
641 
642 
644 (
645  const oldCyclicPolyPatch& pp,
646  const polyBoundaryMesh& bm,
647  const label index,
648  const label newSize,
649  const label newStart
650 )
651 :
652  coupledPolyPatch(pp, bm, index, newSize, newStart),
653  featureCos_(pp.featureCos_),
654  rotationAxis_(pp.rotationAxis_),
655  rotationCentre_(pp.rotationCentre_),
656  separationVector_(pp.separationVector_)
657 {}
658 
659 
660 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
661 
663 {}
664 
665 
666 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
667 
669 {
671 }
672 
673 
675 (
676  const primitivePatch& referPatch,
677  const pointField& thisCtrs,
678  const vectorField& thisAreas,
679  const pointField& thisCc,
680  const pointField& nbrCtrs,
681  const vectorField& nbrAreas,
682  const pointField& nbrCc
683 )
684 {}
685 
686 
688 {}
689 
690 
692 (
693  PstreamBuffers& pBufs,
694  const pointField& p
695 )
696 {
697  polyPatch::initMovePoints(pBufs, p);
698 }
699 
700 
702 (
703  PstreamBuffers& pBufs,
704  const pointField& p
705 )
706 {
707  polyPatch::movePoints(pBufs, p);
708 }
709 
710 
712 {
714 }
715 
716 
718 {
719  polyPatch::updateMesh(pBufs);
720 }
721 
722 
724 (
726  const primitivePatch& pp
727 ) const
728 {}
729 
730 
731 // Return new ordering. Ordering is -faceMap: for every face index
732 // the new face -rotation:for every new face the clockwise shift
733 // of the original face. Return false if nothing changes (faceMap
734 // is identity, rotation is 0)
736 (
738  const primitivePatch& pp,
739  labelList& faceMap,
740  labelList& rotation
741 ) const
742 {
743  faceMap.setSize(pp.size());
744  faceMap = -1;
745 
746  rotation.setSize(pp.size());
747  rotation = 0;
748 
749  if (pp.empty())
750  {
751  // No faces, nothing to change.
752  return false;
753  }
754 
755  if (pp.size()&1)
756  {
758  << "Size of cyclic " << name() << " should be a multiple of 2"
759  << ". It is " << pp.size() << abort(FatalError);
760  }
761 
762  label halfSize = pp.size()/2;
763 
764  // Supplied primitivePatch already with new points.
765  // Cyclics are limited to one transformation tensor
766  // currently anyway (i.e. straight plane) so should not be too big a
767  // problem.
768 
769 
770  // Indices of faces on half0
771  labelList half0ToPatch;
772  // Indices of faces on half1
773  labelList half1ToPatch;
774 
775 
776  // 1. Test if already correctly oriented by starting from trivial ordering.
777  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
778 
779  half0ToPatch = identity(halfSize);
780  half1ToPatch = half0ToPatch + halfSize;
781 
782  // Get faces
783  faceList half0Faces(UIndirectList<face>(pp, half0ToPatch));
784  faceList half1Faces(UIndirectList<face>(pp, half1ToPatch));
785 
786  // Get geometric quantities
787  pointField half0Ctrs, half1Ctrs, anchors0, ppPoints;
788  scalarField tols;
789  getCentresAndAnchors
790  (
791  pp,
792  half0Faces,
793  half1Faces,
794 
795  ppPoints,
796  half0Ctrs,
797  half1Ctrs,
798  anchors0,
799  tols
800  );
801 
802  // Geometric match of face centre vectors
803  labelList from1To0;
804  bool matchedAll = matchPoints
805  (
806  half1Ctrs,
807  half0Ctrs,
808  tols,
809  false,
810  from1To0
811  );
812 
813  if (debug)
814  {
815  Pout<< "oldCyclicPolyPatch::order : test if already ordered:"
816  << matchedAll << endl;
817 
818  // Dump halves
819  fileName nm0("match1_"+name()+"_half0_faces.obj");
820  Pout<< "oldCyclicPolyPatch::order : Writing half0"
821  << " faces to OBJ file " << nm0 << endl;
822  writeOBJ(nm0, half0Faces, ppPoints);
823 
824  fileName nm1("match1_"+name()+"_half1_faces.obj");
825  Pout<< "oldCyclicPolyPatch::order : Writing half1"
826  << " faces to OBJ file " << nm1 << endl;
827  writeOBJ(nm1, half1Faces, pp.points());
828 
829  OFstream ccStr
830  (
831  boundaryMesh().mesh().time().path()
832  /"match1_"+ name() + "_faceCentres.obj"
833  );
834  Pout<< "oldCyclicPolyPatch::order : "
835  << "Dumping currently found cyclic match as lines between"
836  << " corresponding face centres to file " << ccStr.name()
837  << endl;
838 
839  // Recalculate untransformed face centres
840  // pointField rawHalf0Ctrs = calcFaceCentres(half0Faces, pp.points());
841  label vertI = 0;
842 
843  forAll(half1Ctrs, i)
844  {
845  // if (from1To0[i] != -1)
846  {
847  // Write edge between c1 and c0
848  // const point& c0 = rawHalf0Ctrs[from1To0[i]];
849  // const point& c0 = half0Ctrs[from1To0[i]];
850  const point& c0 = half0Ctrs[i];
851  const point& c1 = half1Ctrs[i];
852  writeOBJ(ccStr, c0, c1, vertI);
853  }
854  }
855  }
856 
857 
858  // 2. Ordered in pairs (so 0,1 coupled and 2,3 etc.)
859  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
860 
861  if (!matchedAll)
862  {
863  label facei = 0;
864  for (label i = 0; i < halfSize; i++)
865  {
866  half0ToPatch[i] = facei++;
867  half1ToPatch[i] = facei++;
868  }
869 
870  // And redo all matching
871  half0Faces = UIndirectList<face>(pp, half0ToPatch);
872  half1Faces = UIndirectList<face>(pp, half1ToPatch);
873 
874  getCentresAndAnchors
875  (
876  pp,
877  half0Faces,
878  half1Faces,
879 
880  ppPoints,
881  half0Ctrs,
882  half1Ctrs,
883  anchors0,
884  tols
885  );
886 
887  // Geometric match of face centre vectors
888  matchedAll = matchPoints
889  (
890  half1Ctrs,
891  half0Ctrs,
892  tols,
893  false,
894  from1To0
895  );
896 
897  if (debug)
898  {
899  Pout<< "oldCyclicPolyPatch::order : test if pairwise ordered:"
900  << matchedAll << endl;
901  // Dump halves
902  fileName nm0("match2_"+name()+"_half0_faces.obj");
903  Pout<< "oldCyclicPolyPatch::order : Writing half0"
904  << " faces to OBJ file " << nm0 << endl;
905  writeOBJ(nm0, half0Faces, ppPoints);
906 
907  fileName nm1("match2_"+name()+"_half1_faces.obj");
908  Pout<< "oldCyclicPolyPatch::order : Writing half1"
909  << " faces to OBJ file " << nm1 << endl;
910  writeOBJ(nm1, half1Faces, pp.points());
911 
912  OFstream ccStr
913  (
914  boundaryMesh().mesh().time().path()
915  /"match2_"+name()+"_faceCentres.obj"
916  );
917  Pout<< "oldCyclicPolyPatch::order : "
918  << "Dumping currently found cyclic match as lines between"
919  << " corresponding face centres to file " << ccStr.name()
920  << endl;
921 
922  // Recalculate untransformed face centres
923  label vertI = 0;
924 
925  forAll(half1Ctrs, i)
926  {
927  if (from1To0[i] != -1)
928  {
929  // Write edge between c1 and c0
930  const point& c0 = half0Ctrs[from1To0[i]];
931  const point& c1 = half1Ctrs[i];
932  writeOBJ(ccStr, c0, c1, vertI);
933  }
934  }
935  }
936  }
937 
938 
939  // 3. Baffles(coincident faces) converted into cyclics (e.g. jump)
940  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
941 
942  if (!matchedAll)
943  {
944  label baffleI = 0;
945 
946  forAll(pp, facei)
947  {
948  const face& f = pp.localFaces()[facei];
949  const labelList& pFaces = pp.pointFaces()[f[0]];
950 
951  label matchedFacei = -1;
952 
953  forAll(pFaces, i)
954  {
955  label otherFacei = pFaces[i];
956 
957  if (otherFacei > facei)
958  {
959  const face& otherF = pp.localFaces()[otherFacei];
960 
961  // Note: might pick up two similar oriented faces
962  // (but that is illegal anyway)
963  if (f == otherF)
964  {
965  matchedFacei = otherFacei;
966  break;
967  }
968  }
969  }
970 
971  if (matchedFacei != -1)
972  {
973  half0ToPatch[baffleI] = facei;
974  half1ToPatch[baffleI] = matchedFacei;
975  baffleI++;
976  }
977  }
978 
979  if (baffleI == halfSize)
980  {
981  // And redo all matching
982  half0Faces = UIndirectList<face>(pp, half0ToPatch);
983  half1Faces = UIndirectList<face>(pp, half1ToPatch);
984 
985  getCentresAndAnchors
986  (
987  pp,
988  half0Faces,
989  half1Faces,
990 
991  ppPoints,
992  half0Ctrs,
993  half1Ctrs,
994  anchors0,
995  tols
996  );
997 
998  // Geometric match of face centre vectors
999  matchedAll = matchPoints
1000  (
1001  half1Ctrs,
1002  half0Ctrs,
1003  tols,
1004  false,
1005  from1To0
1006  );
1007 
1008  if (debug)
1009  {
1010  Pout<< "oldCyclicPolyPatch::order : test if baffles:"
1011  << matchedAll << endl;
1012  // Dump halves
1013  fileName nm0("match3_"+name()+"_half0_faces.obj");
1014  Pout<< "oldCyclicPolyPatch::order : Writing half0"
1015  << " faces to OBJ file " << nm0 << endl;
1016  writeOBJ(nm0, half0Faces, ppPoints);
1017 
1018  fileName nm1("match3_"+name()+"_half1_faces.obj");
1019  Pout<< "oldCyclicPolyPatch::order : Writing half1"
1020  << " faces to OBJ file " << nm1 << endl;
1021  writeOBJ(nm1, half1Faces, pp.points());
1022 
1023  OFstream ccStr
1024  (
1025  boundaryMesh().mesh().time().path()
1026  /"match3_"+ name() + "_faceCentres.obj"
1027  );
1028  Pout<< "oldCyclicPolyPatch::order : "
1029  << "Dumping currently found cyclic match as lines between"
1030  << " corresponding face centres to file " << ccStr.name()
1031  << endl;
1032 
1033  // Recalculate untransformed face centres
1034  label vertI = 0;
1035 
1036  forAll(half1Ctrs, i)
1037  {
1038  if (from1To0[i] != -1)
1039  {
1040  // Write edge between c1 and c0
1041  const point& c0 = half0Ctrs[from1To0[i]];
1042  const point& c1 = half1Ctrs[i];
1043  writeOBJ(ccStr, c0, c1, vertI);
1044  }
1045  }
1046  }
1047  }
1048  }
1049 
1050 
1051  // 4. Automatic geometric ordering
1052  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1053 
1054  if (!matchedAll)
1055  {
1056  // Split faces according to feature angle or topology
1057  bool okSplit = getGeometricHalves(pp, half0ToPatch, half1ToPatch);
1058 
1059  if (!okSplit)
1060  {
1061  // Did not split into two equal parts.
1062  return false;
1063  }
1064 
1065  // And redo all matching
1066  half0Faces = UIndirectList<face>(pp, half0ToPatch);
1067  half1Faces = UIndirectList<face>(pp, half1ToPatch);
1068 
1069  getCentresAndAnchors
1070  (
1071  pp,
1072  half0Faces,
1073  half1Faces,
1074 
1075  ppPoints,
1076  half0Ctrs,
1077  half1Ctrs,
1078  anchors0,
1079  tols
1080  );
1081 
1082  // Geometric match of face centre vectors
1083  matchedAll = matchPoints
1084  (
1085  half1Ctrs,
1086  half0Ctrs,
1087  tols,
1088  false,
1089  from1To0
1090  );
1091 
1092  if (debug)
1093  {
1094  Pout<< "oldCyclicPolyPatch::order : automatic ordering result:"
1095  << matchedAll << endl;
1096  // Dump halves
1097  fileName nm0("match4_"+name()+"_half0_faces.obj");
1098  Pout<< "oldCyclicPolyPatch::order : Writing half0"
1099  << " faces to OBJ file " << nm0 << endl;
1100  writeOBJ(nm0, half0Faces, ppPoints);
1101 
1102  fileName nm1("match4_"+name()+"_half1_faces.obj");
1103  Pout<< "oldCyclicPolyPatch::order : Writing half1"
1104  << " faces to OBJ file " << nm1 << endl;
1105  writeOBJ(nm1, half1Faces, pp.points());
1106 
1107  OFstream ccStr
1108  (
1109  boundaryMesh().mesh().time().path()
1110  /"match4_"+ name() + "_faceCentres.obj"
1111  );
1112  Pout<< "oldCyclicPolyPatch::order : "
1113  << "Dumping currently found cyclic match as lines between"
1114  << " corresponding face centres to file " << ccStr.name()
1115  << endl;
1116 
1117  // Recalculate untransformed face centres
1118  label vertI = 0;
1119 
1120  forAll(half1Ctrs, i)
1121  {
1122  if (from1To0[i] != -1)
1123  {
1124  // Write edge between c1 and c0
1125  const point& c0 = half0Ctrs[from1To0[i]];
1126  const point& c1 = half1Ctrs[i];
1127  writeOBJ(ccStr, c0, c1, vertI);
1128  }
1129  }
1130  }
1131  }
1132 
1133 
1134  if (!matchedAll || debug)
1135  {
1136  // Dump halves
1137  fileName nm0(name()+"_half0_faces.obj");
1138  Pout<< "oldCyclicPolyPatch::order : Writing half0"
1139  << " faces to OBJ file " << nm0 << endl;
1140  writeOBJ(nm0, half0Faces, pp.points());
1141 
1142  fileName nm1(name()+"_half1_faces.obj");
1143  Pout<< "oldCyclicPolyPatch::order : Writing half1"
1144  << " faces to OBJ file " << nm1 << endl;
1145  writeOBJ(nm1, half1Faces, pp.points());
1146 
1147  OFstream ccStr
1148  (
1149  boundaryMesh().mesh().time().path()
1150  /name() + "_faceCentres.obj"
1151  );
1152  Pout<< "oldCyclicPolyPatch::order : "
1153  << "Dumping currently found cyclic match as lines between"
1154  << " corresponding face centres to file " << ccStr.name()
1155  << endl;
1156 
1157  // Recalculate untransformed face centres
1158  // pointField rawHalf0Ctrs = calcFaceCentres(half0Faces, pp.points());
1159  label vertI = 0;
1160 
1161  forAll(half1Ctrs, i)
1162  {
1163  if (from1To0[i] != -1)
1164  {
1165  // Write edge between c1 and c0
1166  // const point& c0 = rawHalf0Ctrs[from1To0[i]];
1167  const point& c0 = half0Ctrs[from1To0[i]];
1168  const point& c1 = half1Ctrs[i];
1169  writeOBJ(ccStr, c0, c1, vertI);
1170  }
1171  }
1172  }
1173 
1174 
1175  if (!matchedAll)
1176  {
1178  << "Patch:" << name() << " : "
1179  << "Cannot match vectors to faces on both sides of patch" << endl
1180  << " Perhaps your faces do not match?"
1181  << " The obj files written contain the current match." << endl
1182  << " Continuing with incorrect face ordering from now on!"
1183  << endl;
1184 
1185  return false;
1186  }
1187 
1188 
1189  // Set faceMap such that half0 faces get first and corresponding half1
1190  // faces last.
1191  matchAnchors
1192  (
1193  true, // report if anchor matching error
1194  pp,
1195  half0ToPatch,
1196  anchors0,
1197  half1ToPatch,
1198  half1Faces,
1199  from1To0,
1200  tols,
1201  faceMap,
1202  rotation
1203  );
1204 
1205  // Return false if no change necessary, true otherwise.
1206 
1207  forAll(faceMap, facei)
1208  {
1209  if (faceMap[facei] != facei || rotation[facei] != 0)
1210  {
1211  return true;
1212  }
1213  }
1214 
1215  return false;
1216 }
1217 
1218 
1220 {
1221  // Replacement of polyPatch::write to write 'cyclic' instead of type():
1222  os.writeKeyword("type") << cyclicPolyPatch::typeName
1223  << token::END_STATEMENT << nl;
1225  os.writeKeyword("nFaces") << size() << token::END_STATEMENT << nl;
1226  os.writeKeyword("startFace") << start() << token::END_STATEMENT << nl;
1227 
1228 
1229  os.writeKeyword("featureCos") << featureCos_ << token::END_STATEMENT << nl;
1230  switch (transform())
1231  {
1232  case ROTATIONAL:
1233  {
1234  os.writeKeyword("rotationAxis") << rotationAxis_
1235  << token::END_STATEMENT << nl;
1236  os.writeKeyword("rotationCentre") << rotationCentre_
1237  << token::END_STATEMENT << nl;
1238  break;
1239  }
1240  case TRANSLATIONAL:
1241  {
1242  os.writeKeyword("separationVector") << separationVector_
1243  << token::END_STATEMENT << nl;
1244  break;
1245  }
1246  default:
1247  {
1248  // no additional info to write
1249  }
1250  }
1251 
1253  << "Please run foamUpgradeCyclics to convert these old-style"
1254  << " cyclics into two separate cyclics patches."
1255  << endl;
1256 }
1257 
1258 
1259 // ************************************************************************* //
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:256
&#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 occurrences 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.
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:97
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:265
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.
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
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