processorPolyPatch.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 "processorPolyPatch.H"
28 #include "dictionary.H"
29 #include "SubField.H"
30 #include "demandDrivenData.H"
31 #include "matchPoints.H"
32 #include "OFstream.H"
33 #include "polyMesh.H"
34 #include "Time.H"
35 #include "transformList.H"
36 #include "PstreamBuffers.H"
37 #include "ConstCirculator.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(processorPolyPatch, 0);
44  addToRunTimeSelectionTable(polyPatch, processorPolyPatch, dictionary);
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
51 (
52  const word& name,
53  const label size,
54  const label start,
55  const label index,
56  const polyBoundaryMesh& bm,
57  const int myProcNo,
58  const int neighbProcNo,
59  const transformType transform,
60  const word& patchType
61 )
62 :
63  coupledPolyPatch(name, size, start, index, bm, patchType, transform),
64  myProcNo_(myProcNo),
65  neighbProcNo_(neighbProcNo),
66  neighbFaceCentres_(),
67  neighbFaceAreas_(),
68  neighbFaceCellCentres_()
69 {}
70 
71 
73 (
74  const label size,
75  const label start,
76  const label index,
77  const polyBoundaryMesh& bm,
78  const int myProcNo,
79  const int neighbProcNo,
80  const transformType transform,
81  const word& patchType
82 )
83 :
85  (
86  newName(myProcNo, neighbProcNo),
87  size,
88  start,
89  index,
90  bm,
91  patchType,
92  transform
93  ),
94  myProcNo_(myProcNo),
95  neighbProcNo_(neighbProcNo),
96  neighbFaceCentres_(),
97  neighbFaceAreas_(),
98  neighbFaceCellCentres_()
99 {}
100 
101 
103 (
104  const word& name,
105  const dictionary& dict,
106  const label index,
107  const polyBoundaryMesh& bm,
108  const word& patchType
109 )
110 :
111  coupledPolyPatch(name, dict, index, bm, patchType),
112  myProcNo_(readLabel(dict.lookup("myProcNo"))),
113  neighbProcNo_(readLabel(dict.lookup("neighbProcNo"))),
114  neighbFaceCentres_(),
115  neighbFaceAreas_(),
116  neighbFaceCellCentres_()
117 {}
118 
119 
121 (
122  const processorPolyPatch& pp,
123  const polyBoundaryMesh& bm
124 )
125 :
126  coupledPolyPatch(pp, bm),
127  myProcNo_(pp.myProcNo_),
128  neighbProcNo_(pp.neighbProcNo_),
129  neighbFaceCentres_(),
130  neighbFaceAreas_(),
131  neighbFaceCellCentres_()
132 {}
133 
134 
136 (
137  const processorPolyPatch& pp,
138  const polyBoundaryMesh& bm,
139  const label index,
140  const label newSize,
141  const label newStart
142 )
143 :
144  coupledPolyPatch(pp, bm, index, newSize, newStart),
145  myProcNo_(pp.myProcNo_),
146  neighbProcNo_(pp.neighbProcNo_),
147  neighbFaceCentres_(),
148  neighbFaceAreas_(),
149  neighbFaceCellCentres_()
150 {}
151 
152 
154 (
155  const processorPolyPatch& pp,
156  const polyBoundaryMesh& bm,
157  const label index,
158  const labelUList& mapAddressing,
159  const label newStart
160 )
161 :
162  coupledPolyPatch(pp, bm, index, mapAddressing, newStart),
163  myProcNo_(pp.myProcNo_),
164  neighbProcNo_(pp.neighbProcNo_),
165  neighbFaceCentres_(),
166  neighbFaceAreas_(),
167  neighbFaceCellCentres_()
168 {}
169 
170 
171 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
172 
174 {
175  neighbPointsPtr_.clear();
176  neighbEdgesPtr_.clear();
177 }
178 
179 
180 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
181 
183 (
184  const label myProcNo,
185  const label neighbProcNo
186 )
187 {
188  return
189  "procBoundary"
190  + Foam::name(myProcNo)
191  + "to"
192  + Foam::name(neighbProcNo);
193 }
194 
195 
197 {
198  if (Pstream::parRun())
199  {
200  UOPstream toNeighbProc(neighbProcNo(), pBufs);
201 
202  toNeighbProc
203  << faceCentres()
204  << faceAreas()
205  << faceCellCentres();
206  }
207 }
208 
209 
211 {
212  if (Pstream::parRun())
213  {
214  {
215  UIPstream fromNeighbProc(neighbProcNo(), pBufs);
216 
217  fromNeighbProc
218  >> neighbFaceCentres_
219  >> neighbFaceAreas_
220  >> neighbFaceCellCentres_;
221  }
222 
223  // My normals
224  vectorField faceNormals(size());
225 
226  // Neighbour normals
227  vectorField nbrFaceNormals(neighbFaceAreas_.size());
228 
229  // Face match tolerances
230  scalarField tols = calcFaceTol(*this, points(), faceCentres());
231 
232  // Calculate normals from areas and check
233  forAll(faceNormals, facei)
234  {
235  scalar magSf = mag(faceAreas()[facei]);
236  scalar nbrMagSf = mag(neighbFaceAreas_[facei]);
237  scalar avSf = (magSf + nbrMagSf)/2.0;
238 
239  // For small face area calculation the results of the area
240  // calculation have been found to only be accurate to ~1e-20
241  if (magSf < SMALL || nbrMagSf < SMALL)
242  {
243  // Undetermined normal. Use dummy normal to force separation
244  // check.
245  faceNormals[facei] = point(1, 0, 0);
246  nbrFaceNormals[facei] = -faceNormals[facei];
247  tols[facei] = GREAT;
248  }
249  else if (mag(magSf - nbrMagSf) > matchTolerance()*sqr(tols[facei]))
250  {
251  fileName nm
252  (
253  boundaryMesh().mesh().time().path()
254  /name()+"_faces.obj"
255  );
256 
257  Pout<< "processorPolyPatch::calcGeometry : Writing my "
258  << size()
259  << " faces to OBJ file " << nm << endl;
260 
261  writeOBJ(nm, *this, points());
262 
263  OFstream ccStr
264  (
265  boundaryMesh().mesh().time().path()
266  /name() + "_faceCentresConnections.obj"
267  );
268 
269  Pout<< "processorPolyPatch::calcGeometry :"
270  << " Dumping cell centre lines between"
271  << " corresponding face centres to OBJ file" << ccStr.name()
272  << endl;
273 
274  label vertI = 0;
275 
276  forAll(faceCentres(), facej)
277  {
278  const point& c0 = neighbFaceCentres_[facej];
279  const point& c1 = faceCentres()[facej];
280 
281  writeOBJ(ccStr, c0, c1, vertI);
282  }
283 
285  << "face " << facei << " area does not match neighbour by "
286  << 100*mag(magSf - nbrMagSf)/avSf
287  << "% -- possible face ordering problem." << endl
288  << "patch:" << name()
289  << " my area:" << magSf
290  << " neighbour area:" << nbrMagSf
291  << " matching tolerance:"
292  << matchTolerance()*sqr(tols[facei])
293  << endl
294  << "Mesh face:" << start()+facei
295  << " vertices:"
296  << UIndirectList<point>(points(), operator[](facei))()
297  << endl
298  << "If you are certain your matching is correct"
299  << " you can increase the 'matchTolerance' setting"
300  << " in the patch dictionary in the boundary file."
301  << endl
302  << "Rerun with processor debug flag set for"
303  << " more information." << exit(FatalError);
304  }
305  else
306  {
307  faceNormals[facei] = faceAreas()[facei]/magSf;
308  nbrFaceNormals[facei] = neighbFaceAreas_[facei]/nbrMagSf;
309  }
310  }
311 
312  calcTransformTensors
313  (
314  faceCentres(),
315  neighbFaceCentres_,
316  faceNormals,
317  nbrFaceNormals,
318  matchTolerance()*tols,
319  matchTolerance(),
320  transform()
321  );
322  }
323 }
324 
325 
327 (
328  PstreamBuffers& pBufs,
329  const pointField& p
330 )
331 {
332  polyPatch::movePoints(pBufs, p);
334 }
335 
336 
338 (
339  PstreamBuffers& pBufs,
340  const pointField&
341 )
342 {
344 }
345 
346 
348 {
350 
351  if (Pstream::parRun())
352  {
353  // Express all points as patch face and index in face.
354  labelList pointFace(nPoints());
355  labelList pointIndex(nPoints());
356 
357  for (label patchPointi = 0; patchPointi < nPoints(); patchPointi++)
358  {
359  label facei = pointFaces()[patchPointi][0];
360 
361  pointFace[patchPointi] = facei;
362 
363  const face& f = localFaces()[facei];
364 
365  pointIndex[patchPointi] = findIndex(f, patchPointi);
366  }
367 
368  // Express all edges as patch face and index in face.
369  labelList edgeFace(nEdges());
370  labelList edgeIndex(nEdges());
371 
372  for (label patchEdgeI = 0; patchEdgeI < nEdges(); patchEdgeI++)
373  {
374  label facei = edgeFaces()[patchEdgeI][0];
375 
376  edgeFace[patchEdgeI] = facei;
377 
378  const labelList& fEdges = faceEdges()[facei];
379 
380  edgeIndex[patchEdgeI] = findIndex(fEdges, patchEdgeI);
381  }
382 
383  UOPstream toNeighbProc(neighbProcNo(), pBufs);
384 
385  toNeighbProc
386  << pointFace
387  << pointIndex
388  << edgeFace
389  << edgeIndex;
390  }
391 }
392 
393 
395 {
396  // For completeness
397  polyPatch::updateMesh(pBufs);
398 
399  neighbPointsPtr_.clear();
400  neighbEdgesPtr_.clear();
401 
402  if (Pstream::parRun())
403  {
404  labelList nbrPointFace;
405  labelList nbrPointIndex;
406  labelList nbrEdgeFace;
407  labelList nbrEdgeIndex;
408 
409  {
410  // !Note: there is one situation where the opposite points and
411  // do not exactly match and that is while we are distributing
412  // meshes and multiple parts come together from different
413  // processors. This can temporarily create the situation that
414  // we have points which will be merged out later but we still
415  // need the face connectivity to be correct.
416  // So: cannot check here on same points and edges.
417 
418  UIPstream fromNeighbProc(neighbProcNo(), pBufs);
419 
420  fromNeighbProc
421  >> nbrPointFace
422  >> nbrPointIndex
423  >> nbrEdgeFace
424  >> nbrEdgeIndex;
425  }
426 
427  // Convert neighbour faces and indices into face back into
428  // my edges and points.
429 
430  // Convert points.
431  // ~~~~~~~~~~~~~~~
432 
433  neighbPointsPtr_.reset(new labelList(nPoints(), -1));
434  labelList& neighbPoints = neighbPointsPtr_();
435 
436  forAll(nbrPointFace, nbrPointi)
437  {
438  // Find face and index in face on this side.
439  const face& f = localFaces()[nbrPointFace[nbrPointi]];
440 
441  label index = (f.size() - nbrPointIndex[nbrPointi]) % f.size();
442  label patchPointi = f[index];
443 
444  if (neighbPoints[patchPointi] == -1)
445  {
446  // First reference of point
447  neighbPoints[patchPointi] = nbrPointi;
448  }
449  else if (neighbPoints[patchPointi] >= 0)
450  {
451  // Point already visited. Mark as duplicate.
452  neighbPoints[patchPointi] = -2;
453  }
454  }
455 
456  // Reset all duplicate entries to -1.
457  forAll(neighbPoints, patchPointi)
458  {
459  if (neighbPoints[patchPointi] == -2)
460  {
461  neighbPoints[patchPointi] = -1;
462  }
463  }
464 
465  // Convert edges.
466  // ~~~~~~~~~~~~~~
467 
468  neighbEdgesPtr_.reset(new labelList(nEdges(), -1));
469  labelList& neighbEdges = neighbEdgesPtr_();
470 
471  forAll(nbrEdgeFace, nbrEdgeI)
472  {
473  // Find face and index in face on this side.
474  const labelList& f = faceEdges()[nbrEdgeFace[nbrEdgeI]];
475  label index = (f.size() - nbrEdgeIndex[nbrEdgeI] - 1) % f.size();
476  label patchEdgeI = f[index];
477 
478  if (neighbEdges[patchEdgeI] == -1)
479  {
480  // First reference of edge
481  neighbEdges[patchEdgeI] = nbrEdgeI;
482  }
483  else if (neighbEdges[patchEdgeI] >= 0)
484  {
485  // Edge already visited. Mark as duplicate.
486  neighbEdges[patchEdgeI] = -2;
487  }
488  }
489 
490  // Reset all duplicate entries to -1.
491  forAll(neighbEdges, patchEdgeI)
492  {
493  if (neighbEdges[patchEdgeI] == -2)
494  {
495  neighbEdges[patchEdgeI] = -1;
496  }
497  }
498 
499  // Remove any addressing used for shared points/edges calculation
500  // since mostly not needed.
502  }
503 }
504 
505 
507 {
508  if (!neighbPointsPtr_.valid())
509  {
511  << "No extended addressing calculated for patch " << name()
512  << abort(FatalError);
513  }
514  return neighbPointsPtr_();
515 }
516 
517 
519 {
520  if (!neighbEdgesPtr_.valid())
521  {
523  << "No extended addressing calculated for patch " << name()
524  << abort(FatalError);
525  }
526  return neighbEdgesPtr_();
527 }
528 
529 
531 (
532  PstreamBuffers& pBufs,
533  const primitivePatch& pp
534 ) const
535 {
536  if
537  (
538  !Pstream::parRun()
539  || transform() == NOORDERING
540  )
541  {
542  return;
543  }
544 
545  if (debug)
546  {
547  fileName nm
548  (
549  boundaryMesh().mesh().time().path()
550  /name()+"_faces.obj"
551  );
552  Pout<< "processorPolyPatch::order : Writing my " << pp.size()
553  << " faces to OBJ file " << nm << endl;
554  writeOBJ(nm, pp, pp.points());
555 
556  // Calculate my face centres
557  const pointField& fc = pp.faceCentres();
558 
559  OFstream localStr
560  (
561  boundaryMesh().mesh().time().path()
562  /name() + "_localFaceCentres.obj"
563  );
564  Pout<< "processorPolyPatch::order : "
565  << "Dumping " << fc.size()
566  << " local faceCentres to " << localStr.name() << endl;
567 
568  forAll(fc, facei)
569  {
570  writeOBJ(localStr, fc[facei]);
571  }
572  }
573 
574  if (owner())
575  {
576  if (transform() == COINCIDENTFULLMATCH)
577  {
578  // Pass the patch points and faces across
579  UOPstream toNeighbour(neighbProcNo(), pBufs);
580  toNeighbour << pp.localPoints()
581  << pp.localFaces();
582  }
583  else
584  {
585  const pointField& ppPoints = pp.points();
586 
587  pointField anchors(getAnchorPoints(pp, ppPoints, transform()));
588 
589  // Get the average of the points of each face. This is needed in
590  // case the face centroid calculation is incorrect due to the face
591  // having a very high aspect ratio.
592  pointField facePointAverages(pp.size(), Zero);
593  forAll(pp, fI)
594  {
595  const labelList& facePoints = pp[fI];
596 
597  forAll(facePoints, pI)
598  {
599  facePointAverages[fI] += ppPoints[facePoints[pI]];
600  }
601 
602  facePointAverages[fI] /= facePoints.size();
603  }
604 
605  // Now send all info over to the neighbour
606  UOPstream toNeighbour(neighbProcNo(), pBufs);
607  toNeighbour << pp.faceCentres() << pp.faceNormals()
608  << anchors << facePointAverages;
609  }
610  }
611 }
612 
613 
615 (
616  const face& a,
617  const pointField& aPts,
618  const face& b,
619  const pointField& bPts,
620  const bool sameOrientation,
621  const scalar absTolSqr,
622  scalar& matchDistSqr
623 )
624 {
625  if (a.size() != b.size())
626  {
627  return -1;
628  }
629 
630  enum CirculatorBase::direction circulateDirection
632 
633  if (!sameOrientation)
634  {
635  circulateDirection = CirculatorBase::ANTICLOCKWISE;
636  }
637 
638  label matchFp = -1;
639 
640  scalar closestMatchDistSqr = sqr(GREAT);
641 
642  ConstCirculator<face> aCirc(a);
643  ConstCirculator<face> bCirc(b);
644 
645  do
646  {
647  const scalar diffSqr = magSqr(aPts[aCirc()] - bPts[bCirc()]);
648 
649  if (diffSqr < absTolSqr)
650  {
651  // Found a matching point. Set the fulcrum of b to the iterator
652  ConstCirculator<face> bCirc2 = bCirc;
653  ++aCirc;
654 
655  bCirc2.setFulcrumToIterator();
656 
657  if (!sameOrientation)
658  {
659  --bCirc2;
660  }
661  else
662  {
663  ++bCirc2;
664  }
665 
666  matchDistSqr = diffSqr;
667 
668  do
669  {
670  const scalar diffSqr2 = magSqr(aPts[aCirc()] - bPts[bCirc2()]);
671 
672  if (diffSqr2 > absTolSqr)
673  {
674  // No match.
675  break;
676  }
677 
678  matchDistSqr += diffSqr2;
679  }
680  while
681  (
683  bCirc2.circulate(circulateDirection)
684  );
685 
686  if (!aCirc.circulate())
687  {
688  if (matchDistSqr < closestMatchDistSqr)
689  {
690  closestMatchDistSqr = matchDistSqr;
691 
692  if (!sameOrientation)
693  {
694  matchFp = a.size() - bCirc.nRotations();
695  }
696  else
697  {
698  matchFp = bCirc.nRotations();
699  }
700 
701  if (closestMatchDistSqr == 0)
702  {
703  break;
704  }
705  }
706  }
707 
708  // Reset aCirc
709  aCirc.setIteratorToFulcrum();
710  }
711 
712  } while (bCirc.circulate(circulateDirection));
713 
714  matchDistSqr = closestMatchDistSqr;
715 
716  return matchFp;
717 }
718 
719 
721 (
722  PstreamBuffers& pBufs,
723  const primitivePatch& pp,
724  labelList& faceMap,
725  labelList& rotation
726 ) const
727 {
728  // Note: we only get the faces that originate from internal faces.
729 
730  if
731  (
732  !Pstream::parRun()
733  || transform() == NOORDERING
734  )
735  {
736  return false;
737  }
738 
739  faceMap.setSize(pp.size());
740  faceMap = -1;
741 
742  rotation.setSize(pp.size());
743  rotation = 0;
744 
745  bool change = false;
746 
747  if (owner())
748  {
749  // Do nothing (i.e. identical mapping, zero rotation).
750  // See comment at top.
751  forAll(faceMap, patchFacei)
752  {
753  faceMap[patchFacei] = patchFacei;
754  }
755 
756  if (transform() != COINCIDENTFULLMATCH)
757  {
758  const pointField& ppPoints = pp.points();
759 
760  pointField anchors(getAnchorPoints(pp, ppPoints, transform()));
761 
762  // Calculate typical distance from face centre
763  scalarField tols
764  (
765  matchTolerance()*calcFaceTol(pp, pp.points(), pp.faceCentres())
766  );
767 
768  forAll(faceMap, patchFacei)
769  {
770  const point& wantedAnchor = anchors[patchFacei];
771 
772  rotation[patchFacei] = getRotation
773  (
774  ppPoints,
775  pp[patchFacei],
776  wantedAnchor,
777  tols[patchFacei]
778  );
779 
780  if (rotation[patchFacei] > 0)
781  {
782  change = true;
783  }
784  }
785  }
786 
787  return change;
788  }
789  else
790  {
791  if (transform() == COINCIDENTFULLMATCH)
792  {
793  vectorField masterPts;
794  faceList masterFaces;
795 
796  {
797  UIPstream fromNeighbour(neighbProcNo(), pBufs);
798  fromNeighbour >> masterPts >> masterFaces;
799  }
800 
801  const pointField& localPts = pp.localPoints();
802  const faceList& localFaces = pp.localFaces();
803 
804  label nMatches = 0;
805 
806  forAll(pp, lFacei)
807  {
808  const face& localFace = localFaces[lFacei];
809  label faceRotation = -1;
810 
811  const scalar absTolSqr = sqr(SMALL);
812 
813  scalar closestMatchDistSqr = sqr(GREAT);
814  scalar matchDistSqr = sqr(GREAT);
815  label closestFaceMatch = -1;
816  label closestFaceRotation = -1;
817 
818  forAll(masterFaces, mFacei)
819  {
820  const face& masterFace = masterFaces[mFacei];
821 
822  faceRotation = matchFace
823  (
824  localFace,
825  localPts,
826  masterFace,
827  masterPts,
828  false,
829  absTolSqr,
830  matchDistSqr
831  );
832 
833  if
834  (
835  faceRotation != -1
836  && matchDistSqr < closestMatchDistSqr
837  )
838  {
839  closestMatchDistSqr = matchDistSqr;
840  closestFaceMatch = mFacei;
841  closestFaceRotation = faceRotation;
842  }
843 
844  if (closestMatchDistSqr == 0)
845  {
846  break;
847  }
848  }
849 
850  if (closestFaceRotation != -1 && closestMatchDistSqr == 0)
851  {
852  faceMap[lFacei] = closestFaceMatch;
853 
854  rotation[lFacei] = closestFaceRotation;
855 
856  if (lFacei != closestFaceMatch || closestFaceRotation > 0)
857  {
858  change = true;
859  }
860 
861  nMatches++;
862  }
863  else
864  {
865  Pout<< "Number of matches = " << nMatches << " / "
866  << pp.size() << endl;
867 
868  pointField pts(localFace.size());
869  forAll(localFace, pI)
870  {
871  const label localPtI = localFace[pI];
872  pts[pI] = localPts[localPtI];
873  }
874 
876  << "No match for face " << localFace << nl << pts
877  << abort(FatalError);
878  }
879  }
880 
881  return change;
882  }
883  else
884  {
885  vectorField masterCtrs;
886  vectorField masterNormals;
887  vectorField masterAnchors;
888  vectorField masterFacePointAverages;
889 
890  // Receive data from neighbour
891  {
892  UIPstream fromNeighbour(neighbProcNo(), pBufs);
893  fromNeighbour >> masterCtrs >> masterNormals
894  >> masterAnchors >> masterFacePointAverages;
895  }
896 
897  // Calculate typical distance from face centre
898  scalarField tols
899  (
900  matchTolerance()*calcFaceTol(pp, pp.points(), pp.faceCentres())
901  );
902 
903  if (debug || masterCtrs.size() != pp.size())
904  {
905  {
906  OFstream nbrStr
907  (
908  boundaryMesh().mesh().time().path()
909  /name() + "_nbrFaceCentres.obj"
910  );
911  Pout<< "processorPolyPatch::order : "
912  << "Dumping neighbour faceCentres to " << nbrStr.name()
913  << endl;
914  forAll(masterCtrs, facei)
915  {
916  writeOBJ(nbrStr, masterCtrs[facei]);
917  }
918  }
919 
920  if (masterCtrs.size() != pp.size())
921  {
923  << "in patch:" << name() << " : "
924  << "Local size of patch is " << pp.size() << " (faces)."
925  << endl
926  << "Received from neighbour " << masterCtrs.size()
927  << " faceCentres!"
928  << abort(FatalError);
929  }
930  }
931 
932  // Geometric match of face centre vectors
933  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
934 
935  // 1. Try existing ordering and transformation
936  bool matchedAll = matchPoints
937  (
938  pp.faceCentres(),
939  masterCtrs,
940  pp.faceNormals(),
941  masterNormals,
942  tols,
943  false,
944  faceMap
945  );
946 
947  // Try using face point average for matching
948  if (!matchedAll)
949  {
950  const pointField& ppPoints = pp.points();
951 
952  pointField facePointAverages(pp.size(), Zero);
953  forAll(pp, fI)
954  {
955  const labelList& facePoints = pp[fI];
956 
957  forAll(facePoints, pI)
958  {
959  facePointAverages[fI] += ppPoints[facePoints[pI]];
960  }
961 
962  facePointAverages[fI] /= facePoints.size();
963  }
964 
965  scalarField tols2
966  (
967  calcFaceTol(pp, pp.points(), facePointAverages)
968  );
969 
970  labelList faceMap2(faceMap.size(), -1);
971  matchedAll = matchPoints
972  (
973  facePointAverages,
974  masterFacePointAverages,
975  pp.faceNormals(),
976  masterNormals,
977  tols2,
978  true,
979  faceMap2
980  );
981 
982  forAll(faceMap, oldFacei)
983  {
984  if (faceMap[oldFacei] == -1)
985  {
986  faceMap[oldFacei] = faceMap2[oldFacei];
987  }
988  }
989 
990  matchedAll = true;
991 
992  forAll(faceMap, oldFacei)
993  {
994  if (faceMap[oldFacei] == -1)
995  {
996  matchedAll = false;
997  }
998  }
999  }
1000 
1001  if (!matchedAll || debug)
1002  {
1003  // Dump faces
1004  fileName str
1005  (
1006  boundaryMesh().mesh().time().path()
1007  /name() + "_faces.obj"
1008  );
1009  Pout<< "processorPolyPatch::order :"
1010  << " Writing faces to OBJ file " << str.name() << endl;
1011  writeOBJ(str, pp, pp.points());
1012 
1013  OFstream ccStr
1014  (
1015  boundaryMesh().mesh().time().path()
1016  /name() + "_faceCentresConnections.obj"
1017  );
1018 
1019  Pout<< "processorPolyPatch::order :"
1020  << " Dumping newly found match as lines between"
1021  << " corresponding face centres to OBJ file "
1022  << ccStr.name()
1023  << endl;
1024 
1025  label vertI = 0;
1026 
1027  forAll(pp.faceCentres(), facei)
1028  {
1029  label masterFacei = faceMap[facei];
1030 
1031  if (masterFacei != -1)
1032  {
1033  const point& c0 = masterCtrs[masterFacei];
1034  const point& c1 = pp.faceCentres()[facei];
1035  writeOBJ(ccStr, c0, c1, vertI);
1036  }
1037  }
1038  }
1039 
1040  if (!matchedAll)
1041  {
1043  << "in patch:" << name() << " : "
1044  << "Cannot match vectors to faces on both sides of patch"
1045  << endl
1046  << " masterCtrs[0]:" << masterCtrs[0] << endl
1047  << " ctrs[0]:" << pp.faceCentres()[0] << endl
1048  << " Check your topology changes or maybe you have"
1049  << " multiple separated (from cyclics) processor patches"
1050  << endl
1051  << " Continuing with incorrect face ordering from now on"
1052  << endl;
1053 
1054  return false;
1055  }
1056 
1057  // Set rotation.
1058  forAll(faceMap, oldFacei)
1059  {
1060  // The face f will be at newFacei (after morphing) and we want
1061  // its anchorPoint (= f[0]) to align with the anchorpoint for
1062  // the corresponding face on the other side.
1063 
1064  label newFacei = faceMap[oldFacei];
1065 
1066  const point& wantedAnchor = masterAnchors[newFacei];
1067 
1068  rotation[newFacei] = getRotation
1069  (
1070  pp.points(),
1071  pp[oldFacei],
1072  wantedAnchor,
1073  tols[oldFacei]
1074  );
1075 
1076  if (rotation[newFacei] == -1)
1077  {
1079  << "in patch " << name()
1080  << " : "
1081  << "Cannot find point on face " << pp[oldFacei]
1082  << " with vertices "
1083  << UIndirectList<point>(pp.points(), pp[oldFacei])()
1084  << " that matches point " << wantedAnchor
1085  << " when matching the halves of processor patch "
1086  << name()
1087  << "Continuing with incorrect face ordering from now on"
1088  << endl;
1089 
1090  return false;
1091  }
1092  }
1093 
1094  forAll(faceMap, facei)
1095  {
1096  if (faceMap[facei] != facei || rotation[facei] != 0)
1097  {
1098  return true;
1099  }
1100  }
1101 
1102  return false;
1103  }
1104  }
1105 }
1106 
1107 
1109 {
1111  os.writeKeyword("myProcNo") << myProcNo_
1112  << token::END_STATEMENT << nl;
1113  os.writeKeyword("neighbProcNo") << neighbProcNo_
1114  << token::END_STATEMENT << nl;
1115 }
1116 
1117 
1118 // ************************************************************************* //
#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
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
const labelList & neighbEdges() const
Return neighbour edge labels. WIP.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
Output to file stream.
Definition: OFstream.H:81
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
bool circulate(const CirculatorBase::direction dir=NONE)
Circulate around the list in the given direction.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:59
static word newName(const label myProcNo, const label neighbProcNo)
Return the name of a processorPolyPatch.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
const Field< PointType > & points() const
Return reference to global points.
difference_type nRotations() const
Return the distance between the iterator and the fulcrum. This is.
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
direction
Direction type enumeration.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
void setIteratorToFulcrum()
Set the iterator to the current position of the fulcrum.
void movePoints(PstreamBuffers &, const pointField &)
Correct patches after moving points.
Determine correspondence between points. See below.
Macros for easy insertion into run-time selection tables.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
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.
Neighbour processor patch.
void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
Input inter-processor communications stream operating on external buffer.
Definition: UIPstream.H:53
A list of faces which address into the list of points.
word name() const
Return file name (part beyond last /)
Definition: fileName.C:179
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
dynamicFvMesh & mesh
const pointField & points
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
virtual ~processorPolyPatch()
Destructor.
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
label readLabel(Istream &is)
Definition: label.H:64
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Foam::polyBoundaryMesh.
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
Output inter-processor communications stream operating on external buffer.
Definition: UOPstream.H:54
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.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)
const Field< PointType > & faceCentres() const
Return face centres for patch.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
static label matchFace(const face &localFace, const pointField &localPts, const face &masterFace, const pointField &masterPts, const bool sameOrientation, const scalar absTolSqr, scalar &matchDistSqr)
Returns rotation.
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
labelList f(nPoints)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Spatial transformation functions for primitive fields.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
void setSize(const label)
Reset size of List.
Definition: List.C:295
Template functions to aid in the implementation of demand driven data.
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:393
const Field< PointType > & faceNormals() const
Return face normals for patch.
vector point
Point is a vector.
Definition: point.H:41
A List with indirect addressing.
Definition: fvMatrix.H:106
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
const labelList & neighbPoints() const
Return neighbour point labels. WIP.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Walks over a container as if it were circular. The container must have the following members defined:...
void setFulcrumToIterator()
Set the fulcrum to the current position of the iterator.
processorPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const int myProcNo, const int neighbProcNo, const transformType transform=UNKNOWN, const word &patchType=typeName)
Construct from components with specified name.
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:118
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451