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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "processorPolyPatch.H"
31 #include "dictionary.H"
32 #include "SubField.H"
33 #include "matchPoints.H"
34 #include "OFstream.H"
35 #include "polyMesh.H"
36 #include "Time.H"
37 #include "transformList.H"
38 #include "PstreamBuffers.H"
39 #include "Circulator.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45  defineTypeNameAndDebug(processorPolyPatch, 0);
47 }
48 
49 
50 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
51 
53 (
54  const label myProcNo,
55  const label neighbProcNo
56 )
57 {
58  return word
59  (
60  "procBoundary"
61  + std::to_string(myProcNo) + "to" + std::to_string(neighbProcNo),
62  false // No stripping needed
63  );
64 }
65 
66 
67 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68 
70 (
71  const word& name,
72  const label size,
73  const label start,
74  const label index,
75  const polyBoundaryMesh& bm,
76  const int myProcNo,
77  const int neighbProcNo,
78  const transformType transform,
79  const word& patchType
80 )
81 :
82  coupledPolyPatch(name, size, start, index, bm, patchType, transform),
83  myProcNo_(myProcNo),
84  neighbProcNo_(neighbProcNo),
85  neighbFaceCentres_(),
86  neighbFaceAreas_(),
87  neighbFaceCellCentres_()
88 {}
89 
90 
92 (
93  const label size,
94  const label start,
95  const label index,
96  const polyBoundaryMesh& bm,
97  const int myProcNo,
98  const int neighbProcNo,
99  const transformType transform,
100  const word& patchType
101 )
102 :
104  (
105  newName(myProcNo, neighbProcNo),
106  size,
107  start,
108  index,
109  bm,
110  patchType,
111  transform
112  ),
113  myProcNo_(myProcNo),
114  neighbProcNo_(neighbProcNo),
115  neighbFaceCentres_(),
116  neighbFaceAreas_(),
117  neighbFaceCellCentres_()
118 {}
119 
120 
122 (
123  const word& name,
124  const dictionary& dict,
125  const label index,
126  const polyBoundaryMesh& bm,
127  const word& patchType
128 )
129 :
130  coupledPolyPatch(name, dict, index, bm, patchType),
131  myProcNo_(dict.get<label>("myProcNo")),
132  neighbProcNo_(dict.get<label>("neighbProcNo")),
133  neighbFaceCentres_(),
134  neighbFaceAreas_(),
135  neighbFaceCellCentres_()
136 {}
137 
138 
140 (
141  const processorPolyPatch& pp,
142  const polyBoundaryMesh& bm
143 )
144 :
145  coupledPolyPatch(pp, bm),
146  myProcNo_(pp.myProcNo_),
147  neighbProcNo_(pp.neighbProcNo_),
148  neighbFaceCentres_(),
149  neighbFaceAreas_(),
150  neighbFaceCellCentres_()
151 {}
152 
153 
155 (
156  const processorPolyPatch& pp,
157  const polyBoundaryMesh& bm,
158  const label index,
159  const label newSize,
160  const label newStart
161 )
162 :
163  coupledPolyPatch(pp, bm, index, newSize, newStart),
164  myProcNo_(pp.myProcNo_),
165  neighbProcNo_(pp.neighbProcNo_),
166  neighbFaceCentres_(),
167  neighbFaceAreas_(),
168  neighbFaceCellCentres_()
169 {}
170 
171 
173 (
174  const processorPolyPatch& pp,
175  const polyBoundaryMesh& bm,
176  const label index,
177  const labelUList& mapAddressing,
178  const label newStart
179 )
180 :
181  coupledPolyPatch(pp, bm, index, mapAddressing, newStart),
182  myProcNo_(pp.myProcNo_),
183  neighbProcNo_(pp.neighbProcNo_),
184  neighbFaceCentres_(),
185  neighbFaceAreas_(),
186  neighbFaceCellCentres_()
187 {}
188 
189 
190 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
191 
193 {
194  neighbPointsPtr_.clear();
195  neighbEdgesPtr_.clear();
196 }
197 
198 
199 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
200 
202 {
203  if (Pstream::parRun())
204  {
205  if (neighbProcNo() >= pBufs.nProcs())
206  {
208  << "On patch " << name()
209  << " trying to access out of range neighbour processor "
210  << neighbProcNo() << ". This can happen if" << nl
211  << " trying to run on an incorrect number of processors"
212  << exit(FatalError);
213  }
214 
215  UOPstream toNeighbProc(neighbProcNo(), pBufs);
216 
217  toNeighbProc
218  << faceCentres()
219  << faceAreas()
220  << faceCellCentres();
221  }
222 }
223 
224 
225 void Foam::processorPolyPatch::calcGeometry(PstreamBuffers& pBufs)
226 {
227  if (Pstream::parRun())
228  {
229  {
230  UIPstream fromNeighbProc(neighbProcNo(), pBufs);
231 
232  fromNeighbProc
233  >> neighbFaceCentres_
234  >> neighbFaceAreas_
235  >> neighbFaceCellCentres_;
236  }
237 
238  // My normals
239  vectorField faceNormals(size());
240 
241  // Neighbour normals
242  vectorField nbrFaceNormals(neighbFaceAreas_.size());
243 
244  // Face match tolerances
245  scalarField tols = calcFaceTol(*this, points(), faceCentres());
246 
247  // Calculate normals from areas and check
248  forAll(faceNormals, facei)
249  {
250  const scalar magSf = mag(faceAreas()[facei]);
251  const scalar nbrMagSf = mag(neighbFaceAreas_[facei]);
252 
253  // For small face area calculation the results of the area
254  // calculation have been found to only be accurate to ~1e-20
255  if (magSf < SMALL || nbrMagSf < SMALL)
256  {
257  // Undetermined normal. Use dummy normal to force separation
258  // check.
259  faceNormals[facei] = point(1, 0, 0);
260  nbrFaceNormals[facei] = -faceNormals[facei];
261  tols[facei] = GREAT;
262  }
263  else if (mag(magSf - nbrMagSf) > matchTolerance()*sqr(tols[facei]))
264  {
265  const scalar avgMagSf = 0.5*(magSf + nbrMagSf);
266 
267  fileName nm
268  (
269  boundaryMesh().mesh().time().path()
270  /name()+"_faces.obj"
271  );
272 
273  Pout<< "processorPolyPatch::calcGeometry : Writing my "
274  << size()
275  << " faces to OBJ file " << nm << endl;
276 
277  writeOBJ(nm, *this, points());
278 
279  OFstream ccStr
280  (
281  boundaryMesh().mesh().time().path()
282  /name() + "_faceCentresConnections.obj"
283  );
284 
285  Pout<< "processorPolyPatch::calcGeometry :"
286  << " Dumping cell centre lines between"
287  << " corresponding face centres to OBJ file" << ccStr.name()
288  << endl;
289 
290  label vertI = 0;
291 
292  const vectorField::subField ownFaceCentres = faceCentres();
293 
294  forAll(ownFaceCentres, facej)
295  {
296  const point& c0 = neighbFaceCentres_[facej];
297  const point& c1 = ownFaceCentres[facej];
298 
299  writeOBJ(ccStr, c0, c1, vertI);
300  }
301 
303  << "face " << facei << " area does not match neighbour by "
304  << 100*mag(magSf - nbrMagSf)/avgMagSf
305  << "% -- possible face ordering problem." << endl
306  << "patch:" << name()
307  << " my area:" << magSf
308  << " neighbour area:" << nbrMagSf
309  << " matching tolerance:"
310  << matchTolerance()*sqr(tols[facei])
311  << endl
312  << "Mesh face:" << start()+facei
313  << " vertices:"
314  << UIndirectList<point>(points(), operator[](facei))
315  << endl
316  << "If you are certain your matching is correct"
317  << " you can increase the 'matchTolerance' setting"
318  << " in the patch dictionary in the boundary file."
319  << endl
320  << "Rerun with processor debug flag set for"
321  << " more information." << exit(FatalError);
322  }
323  else
324  {
325  faceNormals[facei] = faceAreas()[facei]/magSf;
326  nbrFaceNormals[facei] = neighbFaceAreas_[facei]/nbrMagSf;
327  }
328  }
329 
330  calcTransformTensors
331  (
332  faceCentres(),
333  neighbFaceCentres_,
334  faceNormals,
335  nbrFaceNormals,
336  matchTolerance()*tols,
337  matchTolerance(),
339  );
340  }
341 }
342 
343 
345 (
346  PstreamBuffers& pBufs,
347  const pointField& p
348 )
349 {
350  polyPatch::movePoints(pBufs, p);
352 }
353 
354 
356 (
357  PstreamBuffers& pBufs,
358  const pointField&
359 )
360 {
362 }
363 
364 
366 {
368 
369  if (Pstream::parRun())
370  {
371  if (neighbProcNo() >= pBufs.nProcs())
372  {
374  << "On patch " << name()
375  << " trying to access out of range neighbour processor "
376  << neighbProcNo() << ". This can happen if" << nl
377  << " trying to run on an incorrect number of processors"
378  << exit(FatalError);
379  }
380 
381  // Express all points as patch face and index in face.
382  labelList pointFace(nPoints());
383  labelList pointIndex(nPoints());
384 
385  for (label patchPointi = 0; patchPointi < nPoints(); patchPointi++)
386  {
387  label facei = pointFaces()[patchPointi][0];
388 
389  pointFace[patchPointi] = facei;
390 
391  const face& f = localFaces()[facei];
392 
393  pointIndex[patchPointi] = f.find(patchPointi);
394  }
395 
396  // Express all edges as patch face and index in face.
397  labelList edgeFace(nEdges());
398  labelList edgeIndex(nEdges());
399 
400  for (label patchEdgeI = 0; patchEdgeI < nEdges(); patchEdgeI++)
401  {
402  label facei = edgeFaces()[patchEdgeI][0];
403 
404  edgeFace[patchEdgeI] = facei;
405 
406  const labelList& fEdges = faceEdges()[facei];
407 
408  edgeIndex[patchEdgeI] = fEdges.find(patchEdgeI);
409  }
410 
411  UOPstream toNeighbProc(neighbProcNo(), pBufs);
412 
413  toNeighbProc
414  << pointFace
415  << pointIndex
416  << edgeFace
417  << edgeIndex;
418  }
419 }
420 
421 
422 void Foam::processorPolyPatch::updateMesh(PstreamBuffers& pBufs)
423 {
424  // For completeness
425  polyPatch::updateMesh(pBufs);
426 
427  neighbPointsPtr_.clear();
428  neighbEdgesPtr_.clear();
429 
430  if (Pstream::parRun())
431  {
432  labelList nbrPointFace;
433  labelList nbrPointIndex;
434  labelList nbrEdgeFace;
435  labelList nbrEdgeIndex;
436 
437  {
438  // !Note: there is one situation where the opposite points and
439  // do not exactly match and that is while we are distributing
440  // meshes and multiple parts come together from different
441  // processors. This can temporarily create the situation that
442  // we have points which will be merged out later but we still
443  // need the face connectivity to be correct.
444  // So: cannot check here on same points and edges.
445 
446  UIPstream fromNeighbProc(neighbProcNo(), pBufs);
447 
448  fromNeighbProc
449  >> nbrPointFace
450  >> nbrPointIndex
451  >> nbrEdgeFace
452  >> nbrEdgeIndex;
453  }
454 
455  // Convert neighbour faces and indices into face back into
456  // my edges and points.
457 
458  // Convert points.
459  // ~~~~~~~~~~~~~~~
460 
461  neighbPointsPtr_.reset(new labelList(nPoints(), -1));
462  labelList& neighbPoints = neighbPointsPtr_();
463 
464  forAll(nbrPointFace, nbrPointi)
465  {
466  // Find face and index in face on this side.
467  const face& f = localFaces()[nbrPointFace[nbrPointi]];
468 
469  label index = (f.size() - nbrPointIndex[nbrPointi]) % f.size();
470  label patchPointi = f[index];
471 
472  if (neighbPoints[patchPointi] == -1)
473  {
474  // First reference of point
475  neighbPoints[patchPointi] = nbrPointi;
476  }
477  else if (neighbPoints[patchPointi] >= 0)
478  {
479  // Point already visited. Mark as duplicate.
480  neighbPoints[patchPointi] = -2;
481  }
482  }
483 
484  // Reset all duplicate entries to -1.
485  forAll(neighbPoints, patchPointi)
486  {
487  if (neighbPoints[patchPointi] == -2)
488  {
489  neighbPoints[patchPointi] = -1;
490  }
491  }
492 
493  // Convert edges.
494  // ~~~~~~~~~~~~~~
495 
496  neighbEdgesPtr_.reset(new labelList(nEdges(), -1));
497  labelList& neighbEdges = neighbEdgesPtr_();
498 
499  forAll(nbrEdgeFace, nbrEdgeI)
500  {
501  // Find face and index in face on this side.
502  const labelList& f = faceEdges()[nbrEdgeFace[nbrEdgeI]];
503  label index = (f.size() - nbrEdgeIndex[nbrEdgeI] - 1) % f.size();
504  label patchEdgeI = f[index];
505 
506  if (neighbEdges[patchEdgeI] == -1)
507  {
508  // First reference of edge
509  neighbEdges[patchEdgeI] = nbrEdgeI;
510  }
511  else if (neighbEdges[patchEdgeI] >= 0)
512  {
513  // Edge already visited. Mark as duplicate.
514  neighbEdges[patchEdgeI] = -2;
515  }
516  }
517 
518  // Reset all duplicate entries to -1.
519  forAll(neighbEdges, patchEdgeI)
520  {
521  if (neighbEdges[patchEdgeI] == -2)
522  {
523  neighbEdges[patchEdgeI] = -1;
524  }
525  }
526 
527  // Remove any addressing used for shared points/edges calculation
528  // since mostly not needed.
530  }
531 }
532 
533 
535 {
536  if (!neighbPointsPtr_)
537  {
539  << "No extended addressing calculated for patch " << name()
540  << abort(FatalError);
541  }
542  return *neighbPointsPtr_;
543 }
544 
545 
547 {
548  if (!neighbEdgesPtr_)
549  {
551  << "No extended addressing calculated for patch " << name()
553  }
554  return *neighbEdgesPtr_;
555 }
556 
557 
559 (
560  PstreamBuffers& pBufs,
561  const primitivePatch& pp
562 ) const
563 {
564  if
565  (
566  !Pstream::parRun()
567  || transform() == NOORDERING
568  )
569  {
570  return;
571  }
572 
573  if (debug)
574  {
575  fileName nm
576  (
577  boundaryMesh().mesh().time().path()
578  /name()+"_faces.obj"
579  );
580  Pout<< "processorPolyPatch::order : Writing my " << pp.size()
581  << " faces to OBJ file " << nm << endl;
582  writeOBJ(nm, pp, pp.points());
583 
584  // Calculate my face centres
585  const pointField& fc = pp.faceCentres();
586 
587  OFstream localStr
588  (
589  boundaryMesh().mesh().time().path()
590  /name() + "_localFaceCentres.obj"
591  );
592  Pout<< "processorPolyPatch::order : "
593  << "Dumping " << fc.size()
594  << " local faceCentres to " << localStr.name() << endl;
595 
596  forAll(fc, facei)
597  {
598  writeOBJ(localStr, fc[facei]);
599  }
600  }
601 
602  if (owner())
603  {
604  if (transform() == COINCIDENTFULLMATCH)
605  {
606  // Pass the patch points and faces across
607  UOPstream toNeighbour(neighbProcNo(), pBufs);
608  toNeighbour << pp.localPoints()
609  << pp.localFaces();
610  }
611  else
612  {
613  const pointField& ppPoints = pp.points();
614 
615  pointField anchors(getAnchorPoints(pp, ppPoints, transform()));
616 
617  // Get the average of the points of each face. This is needed in
618  // case the face centroid calculation is incorrect due to the face
619  // having a very high aspect ratio.
620  pointField facePointAverages(pp.size(), Zero);
621  forAll(pp, fI)
622  {
623  const labelList& facePoints = pp[fI];
624 
625  forAll(facePoints, pI)
626  {
627  facePointAverages[fI] += ppPoints[facePoints[pI]];
628  }
629 
630  facePointAverages[fI] /= facePoints.size();
631  }
632 
633  // Now send all info over to the neighbour
634  UOPstream toNeighbour(neighbProcNo(), pBufs);
635  toNeighbour << pp.faceCentres() << pp.faceNormals()
636  << anchors << facePointAverages;
637  }
638  }
639 }
640 
641 
643 (
644  const face& a,
645  const pointField& aPts,
646  const face& b,
647  const pointField& bPts,
648  const bool sameOrientation,
649  const scalar absTolSqr,
650  scalar& matchDistSqr
651 )
652 {
653  if (a.size() != b.size())
654  {
655  return -1;
656  }
657 
658  enum CirculatorBase::direction circulateDirection
660 
661  if (!sameOrientation)
662  {
663  circulateDirection = CirculatorBase::ANTICLOCKWISE;
664  }
665 
666  label matchFp = -1;
667 
668  scalar closestMatchDistSqr = sqr(GREAT);
669 
670  ConstCirculator<face> aCirc(a);
671  ConstCirculator<face> bCirc(b);
672 
673  do
674  {
675  const scalar diffSqr = magSqr(aPts[*aCirc] - bPts[*bCirc]);
676 
677  if (diffSqr < absTolSqr)
678  {
679  // Found a matching point. Set the fulcrum of b to the iterator
680  ConstCirculator<face> bCirc2(bCirc);
681  ++aCirc;
682 
683  bCirc2.setFulcrumToIterator();
684 
685  if (!sameOrientation)
686  {
687  --bCirc2;
688  }
689  else
690  {
691  ++bCirc2;
692  }
693 
694  matchDistSqr = diffSqr;
695 
696  do
697  {
698  const scalar diffSqr2 = magSqr(aPts[*aCirc] - bPts[*bCirc2]);
699 
700  if (diffSqr2 > absTolSqr)
701  {
702  // No match.
703  break;
704  }
705 
706  matchDistSqr += diffSqr2;
707  }
708  while
709  (
710  aCirc.circulate(CirculatorBase::CLOCKWISE),
711  bCirc2.circulate(circulateDirection)
712  );
713 
714  if (!aCirc.circulate())
715  {
716  if (matchDistSqr < closestMatchDistSqr)
717  {
718  closestMatchDistSqr = matchDistSqr;
719 
720  if (!sameOrientation)
721  {
722  matchFp = a.size() - bCirc.nRotations();
723  }
724  else
725  {
726  matchFp = bCirc.nRotations();
727  }
728 
729  if (closestMatchDistSqr == 0)
730  {
731  break;
732  }
733  }
734  }
735 
736  // Reset aCirc
737  aCirc.setIteratorToFulcrum();
738  }
739 
740  } while (bCirc.circulate(circulateDirection));
741 
742  matchDistSqr = closestMatchDistSqr;
743 
744  return matchFp;
745 }
746 
747 
749 (
750  PstreamBuffers& pBufs,
751  const primitivePatch& pp,
753  labelList& rotation
754 ) const
755 {
756  // Note: we only get the faces that originate from internal faces.
757 
758  if
759  (
760  !Pstream::parRun()
761  || transform() == NOORDERING
762  )
763  {
764  return false;
765  }
766 
767  faceMap.setSize(pp.size());
768  faceMap = -1;
769 
770  rotation.setSize(pp.size());
771  rotation = 0;
772 
773  bool change = false;
774 
775  if (owner())
776  {
777  // Do nothing (i.e. identical mapping, zero rotation).
778  // See comment at top.
779  forAll(faceMap, patchFacei)
780  {
781  faceMap[patchFacei] = patchFacei;
782  }
783 
784  if (transform() != COINCIDENTFULLMATCH)
785  {
786  const pointField& ppPoints = pp.points();
787 
788  pointField anchors(getAnchorPoints(pp, ppPoints, transform()));
789 
790  // Calculate typical distance from face centre
791  scalarField tols
792  (
793  matchTolerance()*calcFaceTol(pp, pp.points(), pp.faceCentres())
794  );
795 
796  forAll(faceMap, patchFacei)
797  {
798  const point& wantedAnchor = anchors[patchFacei];
799 
800  rotation[patchFacei] = getRotation
801  (
802  ppPoints,
803  pp[patchFacei],
804  wantedAnchor,
805  tols[patchFacei]
806  );
807 
808  if (rotation[patchFacei] > 0)
809  {
810  change = true;
811  }
812  }
813  }
814 
815  return change;
816  }
817  else
818  {
819  // Calculate the absolute matching tolerance
820  scalarField tols
821  (
822  matchTolerance()*calcFaceTol(pp, pp.points(), pp.faceCentres())
823  );
824 
825  if (transform() == COINCIDENTFULLMATCH)
826  {
827  vectorField masterPts;
828  faceList masterFaces;
829 
830  {
831  UIPstream fromNeighbour(neighbProcNo(), pBufs);
832  fromNeighbour >> masterPts >> masterFaces;
833  }
834 
835  const pointField& localPts = pp.localPoints();
836  const faceList& localFaces = pp.localFaces();
837 
838  label nMatches = 0;
839 
840  forAll(pp, lFacei)
841  {
842  const face& localFace = localFaces[lFacei];
843  label faceRotation = -1;
844 
845  const scalar absTolSqr = sqr(tols[lFacei]);
846 
847  scalar closestMatchDistSqr = sqr(GREAT);
848  scalar matchDistSqr = sqr(GREAT);
849  label closestFaceMatch = -1;
850  label closestFaceRotation = -1;
851 
852  forAll(masterFaces, mFacei)
853  {
854  const face& masterFace = masterFaces[mFacei];
855 
856  faceRotation = matchFace
857  (
858  localFace,
859  localPts,
860  masterFace,
861  masterPts,
862  false,
863  absTolSqr,
864  matchDistSqr
865  );
866 
867  if
868  (
869  faceRotation != -1
870  && matchDistSqr < closestMatchDistSqr
871  )
872  {
873  closestMatchDistSqr = matchDistSqr;
874  closestFaceMatch = mFacei;
875  closestFaceRotation = faceRotation;
876  }
877 
878  if (closestMatchDistSqr == 0)
879  {
880  break;
881  }
882  }
883 
884  if
885  (
886  closestFaceRotation != -1
887  && closestMatchDistSqr < absTolSqr
888  )
889  {
890  faceMap[lFacei] = closestFaceMatch;
891 
892  rotation[lFacei] = closestFaceRotation;
893 
894  if (lFacei != closestFaceMatch || closestFaceRotation > 0)
895  {
896  change = true;
897  }
898 
899  nMatches++;
900  }
901  else
902  {
903  Pout<< "Number of matches = " << nMatches << " / "
904  << pp.size() << endl;
905 
906  pointField pts(localFace.size());
907  forAll(localFace, pI)
908  {
909  const label localPtI = localFace[pI];
910  pts[pI] = localPts[localPtI];
911  }
912 
914  << "No match for face " << localFace << nl << pts
915  << abort(FatalError);
916  }
917  }
918 
919  return change;
920  }
921  else
922  {
923  vectorField masterCtrs;
924  vectorField masterNormals;
925  vectorField masterAnchors;
926  vectorField masterFacePointAverages;
927 
928  // Receive data from neighbour
929  {
930  UIPstream fromNeighbour(neighbProcNo(), pBufs);
931  fromNeighbour >> masterCtrs >> masterNormals
932  >> masterAnchors >> masterFacePointAverages;
933  }
934 
935  if (debug || masterCtrs.size() != pp.size())
936  {
937  {
938  OFstream nbrStr
939  (
940  boundaryMesh().mesh().time().path()
941  /name() + "_nbrFaceCentres.obj"
942  );
943  Pout<< "processorPolyPatch::order : "
944  << "Dumping neighbour faceCentres to " << nbrStr.name()
945  << endl;
946  forAll(masterCtrs, facei)
947  {
948  writeOBJ(nbrStr, masterCtrs[facei]);
949  }
950  }
951 
952  if (masterCtrs.size() != pp.size())
953  {
955  << "in patch:" << name() << " : "
956  << "Local size of patch is " << pp.size() << " (faces)."
957  << endl
958  << "Received from neighbour " << masterCtrs.size()
959  << " faceCentres!"
960  << abort(FatalError);
961  }
962  }
963 
964  // Geometric match of face centre vectors
965  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
966 
967  // 1. Try existing ordering and transformation
968  bool matchedAll = matchPoints
969  (
970  pp.faceCentres(),
971  masterCtrs,
972  pp.faceNormals(),
973  masterNormals,
974  tols,
975  false,
976  faceMap
977  );
978 
979  // Fallback: try using face point average for matching
980  if (!matchedAll)
981  {
982  const pointField& ppPoints = pp.points();
983 
984  pointField facePointAverages(pp.size(), Zero);
985  forAll(pp, fI)
986  {
987  const labelList& facePoints = pp[fI];
988 
989  forAll(facePoints, pI)
990  {
991  facePointAverages[fI] += ppPoints[facePoints[pI]];
992  }
993 
994  facePointAverages[fI] /= facePoints.size();
995  }
996 
997  scalarField tols2
998  (
999  matchTolerance()
1000  *calcFaceTol(pp, pp.points(), facePointAverages)
1001  );
1002 
1003  // Note that we do not use the faceNormals anymore for
1004  // comparison. Since we're
1005  // having problems with the face centres (e.g. due to extreme
1006  // aspect ratios) we will probably also have problems with
1007  // reliable normals calculation
1008  labelList faceMap2(faceMap.size(), -1);
1009  matchPoints
1010  (
1011  facePointAverages,
1012  masterFacePointAverages,
1013  tols2,
1014  true,
1015  faceMap2
1016  );
1017 
1018  matchedAll = true;
1019 
1020  forAll(faceMap, oldFacei)
1021  {
1022  if (faceMap[oldFacei] == -1)
1023  {
1024  faceMap[oldFacei] = faceMap2[oldFacei];
1025 
1026  if (faceMap[oldFacei] == -1)
1027  {
1028  matchedAll = false;
1029  }
1030  }
1031  }
1032  }
1033 
1034  if (!matchedAll || debug)
1035  {
1036  // Dump faces
1037  fileName str
1038  (
1039  boundaryMesh().mesh().time().path()
1040  /name() + "_faces.obj"
1041  );
1042  Pout<< "processorPolyPatch::order :"
1043  << " Writing faces to OBJ file " << str.name() << endl;
1044  writeOBJ(str, pp, pp.points());
1045 
1046  OFstream ccStr
1047  (
1048  boundaryMesh().mesh().time().path()
1049  /name() + "_faceCentresConnections.obj"
1050  );
1051 
1052  Pout<< "processorPolyPatch::order :"
1053  << " Dumping newly found match as lines between"
1054  << " corresponding face centres to OBJ file "
1055  << ccStr.name()
1056  << endl;
1057 
1058  label vertI = 0;
1059 
1060  forAll(pp.faceCentres(), facei)
1061  {
1062  label masterFacei = faceMap[facei];
1063 
1064  if (masterFacei != -1)
1065  {
1066  const point& c0 = masterCtrs[masterFacei];
1067  const point& c1 = pp.faceCentres()[facei];
1068  writeOBJ(ccStr, c0, c1, vertI);
1069  }
1070  }
1071  }
1072 
1073  if (!matchedAll)
1074  {
1076  << "in patch:" << name() << " : "
1077  << "Cannot match vectors to faces on both sides of patch"
1078  << endl
1079  << " masterCtrs[0]:" << masterCtrs[0] << endl
1080  << " ctrs[0]:" << pp.faceCentres()[0] << endl
1081  << " Check your topology changes or maybe you have"
1082  << " multiple separated (from cyclics) processor patches"
1083  << endl
1084  << " Continuing with incorrect face ordering from now on"
1085  << endl;
1086 
1087  return false;
1088  }
1089 
1090  // Set rotation.
1091  forAll(faceMap, oldFacei)
1092  {
1093  // The face f will be at newFacei (after morphing) and we want
1094  // its anchorPoint (= f[0]) to align with the anchorpoint for
1095  // the corresponding face on the other side.
1096 
1097  label newFacei = faceMap[oldFacei];
1098 
1099  const point& wantedAnchor = masterAnchors[newFacei];
1100 
1101  rotation[newFacei] = getRotation
1102  (
1103  pp.points(),
1104  pp[oldFacei],
1105  wantedAnchor,
1106  tols[oldFacei]
1107  );
1108 
1109  if (rotation[newFacei] == -1)
1110  {
1112  << "in patch " << name()
1113  << " : "
1114  << "Cannot find point on face " << pp[oldFacei]
1115  << " with vertices "
1116  << UIndirectList<point>(pp.points(), pp[oldFacei])
1117  << " that matches point " << wantedAnchor
1118  << " when matching the halves of processor patch "
1119  << name()
1120  << "Continuing with incorrect face ordering from now on"
1121  << endl;
1122 
1123  return false;
1124  }
1125  }
1126 
1127  forAll(faceMap, facei)
1128  {
1129  if (faceMap[facei] != facei || rotation[facei] != 0)
1130  {
1131  return true;
1132  }
1133  }
1135  return false;
1136  }
1137  }
1138 }
1139 
1140 
1141 void Foam::processorPolyPatch::write(Ostream& os) const
1142 {
1144  os.writeEntry("myProcNo", myProcNo_);
1145  os.writeEntry("neighbProcNo", neighbProcNo_);
1146 }
1147 
1148 
1149 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const labelList & neighbPoints() const
Return neighbour point labels. WIP.
virtual void write(Ostream &os) const
Write the polyPatch data as a dictionary.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const Field< point_type > & localPoints() const
Return pointField of points in patch.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: polyPatch.H:140
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
static word newName(const label myProcNo, const label neighbProcNo)
Return the name of a processorPolyPatch ("procBoundary..") constructed from the pair of processor IDs...
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1061
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
direction
Direction type enumeration.
Definition: Circulator.H:86
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:196
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:134
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:53
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
Neighbour processor patch.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
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:29
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
void setSize(const label n)
Alias for resize()
Definition: List.H:320
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
A class for handling words, derived from Foam::string.
Definition: word.H:63
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const Field< point_type > & faceCentres() const
Return face centres for patch.
virtual ~processorPolyPatch()
Destructor.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: polyPatch.C:59
const Field< point_type > & points() const noexcept
Return reference to global points.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
label find(const T &val) const
Find index of the first occurrence of the value.
Definition: UList.C:173
void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
Output inter-processor communications stream using MPI send/recv etc. - operating on external buffer...
Definition: UOPstream.H:394
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
defineTypeNameAndDebug(combustionModel, 0)
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.
labelList f(nPoints)
Buffers for inter-processor communications streams (UOPstream, UIPstream).
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field...
Spatial transformation functions for list of values and primitive fields.
int neighbProcNo() const noexcept
Return neighbour processor number.
vector point
Point is a vector.
Definition: point.H:37
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
SubField< vector > subField
Declare type of subField.
Definition: Field.H:128
Field< vector > vectorField
Specialisation of Field<T> for vector.
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
const labelList & neighbEdges() const
Return neighbour edge labels. WIP.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
int myProcNo() const noexcept
Return processor number.
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.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127