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  forAll(faceCentres(), facej)
293  {
294  const point& c0 = neighbFaceCentres_[facej];
295  const point& c1 = faceCentres()[facej];
296 
297  writeOBJ(ccStr, c0, c1, vertI);
298  }
299 
301  << "face " << facei << " area does not match neighbour by "
302  << 100*mag(magSf - nbrMagSf)/avgMagSf
303  << "% -- possible face ordering problem." << endl
304  << "patch:" << name()
305  << " my area:" << magSf
306  << " neighbour area:" << nbrMagSf
307  << " matching tolerance:"
308  << matchTolerance()*sqr(tols[facei])
309  << endl
310  << "Mesh face:" << start()+facei
311  << " vertices:"
312  << UIndirectList<point>(points(), operator[](facei))
313  << endl
314  << "If you are certain your matching is correct"
315  << " you can increase the 'matchTolerance' setting"
316  << " in the patch dictionary in the boundary file."
317  << endl
318  << "Rerun with processor debug flag set for"
319  << " more information." << exit(FatalError);
320  }
321  else
322  {
323  faceNormals[facei] = faceAreas()[facei]/magSf;
324  nbrFaceNormals[facei] = neighbFaceAreas_[facei]/nbrMagSf;
325  }
326  }
327 
328  calcTransformTensors
329  (
330  faceCentres(),
331  neighbFaceCentres_,
332  faceNormals,
333  nbrFaceNormals,
334  matchTolerance()*tols,
335  matchTolerance(),
337  );
338  }
339 }
340 
341 
343 (
344  PstreamBuffers& pBufs,
345  const pointField& p
346 )
347 {
348  polyPatch::movePoints(pBufs, p);
350 }
351 
352 
354 (
355  PstreamBuffers& pBufs,
356  const pointField&
357 )
358 {
360 }
361 
362 
364 {
366 
367  if (Pstream::parRun())
368  {
369  if (neighbProcNo() >= pBufs.nProcs())
370  {
372  << "On patch " << name()
373  << " trying to access out of range neighbour processor "
374  << neighbProcNo() << ". This can happen if" << nl
375  << " trying to run on an incorrect number of processors"
376  << exit(FatalError);
377  }
378 
379  // Express all points as patch face and index in face.
380  labelList pointFace(nPoints());
381  labelList pointIndex(nPoints());
382 
383  for (label patchPointi = 0; patchPointi < nPoints(); patchPointi++)
384  {
385  label facei = pointFaces()[patchPointi][0];
386 
387  pointFace[patchPointi] = facei;
388 
389  const face& f = localFaces()[facei];
390 
391  pointIndex[patchPointi] = f.find(patchPointi);
392  }
393 
394  // Express all edges as patch face and index in face.
395  labelList edgeFace(nEdges());
396  labelList edgeIndex(nEdges());
397 
398  for (label patchEdgeI = 0; patchEdgeI < nEdges(); patchEdgeI++)
399  {
400  label facei = edgeFaces()[patchEdgeI][0];
401 
402  edgeFace[patchEdgeI] = facei;
403 
404  const labelList& fEdges = faceEdges()[facei];
405 
406  edgeIndex[patchEdgeI] = fEdges.find(patchEdgeI);
407  }
408 
409  UOPstream toNeighbProc(neighbProcNo(), pBufs);
410 
411  toNeighbProc
412  << pointFace
413  << pointIndex
414  << edgeFace
415  << edgeIndex;
416  }
417 }
418 
419 
420 void Foam::processorPolyPatch::updateMesh(PstreamBuffers& pBufs)
421 {
422  // For completeness
423  polyPatch::updateMesh(pBufs);
424 
425  neighbPointsPtr_.clear();
426  neighbEdgesPtr_.clear();
427 
428  if (Pstream::parRun())
429  {
430  labelList nbrPointFace;
431  labelList nbrPointIndex;
432  labelList nbrEdgeFace;
433  labelList nbrEdgeIndex;
434 
435  {
436  // !Note: there is one situation where the opposite points and
437  // do not exactly match and that is while we are distributing
438  // meshes and multiple parts come together from different
439  // processors. This can temporarily create the situation that
440  // we have points which will be merged out later but we still
441  // need the face connectivity to be correct.
442  // So: cannot check here on same points and edges.
443 
444  UIPstream fromNeighbProc(neighbProcNo(), pBufs);
445 
446  fromNeighbProc
447  >> nbrPointFace
448  >> nbrPointIndex
449  >> nbrEdgeFace
450  >> nbrEdgeIndex;
451  }
452 
453  // Convert neighbour faces and indices into face back into
454  // my edges and points.
455 
456  // Convert points.
457  // ~~~~~~~~~~~~~~~
458 
459  neighbPointsPtr_.reset(new labelList(nPoints(), -1));
460  labelList& neighbPoints = neighbPointsPtr_();
461 
462  forAll(nbrPointFace, nbrPointi)
463  {
464  // Find face and index in face on this side.
465  const face& f = localFaces()[nbrPointFace[nbrPointi]];
466 
467  label index = (f.size() - nbrPointIndex[nbrPointi]) % f.size();
468  label patchPointi = f[index];
469 
470  if (neighbPoints[patchPointi] == -1)
471  {
472  // First reference of point
473  neighbPoints[patchPointi] = nbrPointi;
474  }
475  else if (neighbPoints[patchPointi] >= 0)
476  {
477  // Point already visited. Mark as duplicate.
478  neighbPoints[patchPointi] = -2;
479  }
480  }
481 
482  // Reset all duplicate entries to -1.
483  forAll(neighbPoints, patchPointi)
484  {
485  if (neighbPoints[patchPointi] == -2)
486  {
487  neighbPoints[patchPointi] = -1;
488  }
489  }
490 
491  // Convert edges.
492  // ~~~~~~~~~~~~~~
493 
494  neighbEdgesPtr_.reset(new labelList(nEdges(), -1));
495  labelList& neighbEdges = neighbEdgesPtr_();
496 
497  forAll(nbrEdgeFace, nbrEdgeI)
498  {
499  // Find face and index in face on this side.
500  const labelList& f = faceEdges()[nbrEdgeFace[nbrEdgeI]];
501  label index = (f.size() - nbrEdgeIndex[nbrEdgeI] - 1) % f.size();
502  label patchEdgeI = f[index];
503 
504  if (neighbEdges[patchEdgeI] == -1)
505  {
506  // First reference of edge
507  neighbEdges[patchEdgeI] = nbrEdgeI;
508  }
509  else if (neighbEdges[patchEdgeI] >= 0)
510  {
511  // Edge already visited. Mark as duplicate.
512  neighbEdges[patchEdgeI] = -2;
513  }
514  }
515 
516  // Reset all duplicate entries to -1.
517  forAll(neighbEdges, patchEdgeI)
518  {
519  if (neighbEdges[patchEdgeI] == -2)
520  {
521  neighbEdges[patchEdgeI] = -1;
522  }
523  }
524 
525  // Remove any addressing used for shared points/edges calculation
526  // since mostly not needed.
528  }
529 }
530 
531 
533 {
534  if (!neighbPointsPtr_)
535  {
537  << "No extended addressing calculated for patch " << name()
538  << abort(FatalError);
539  }
540  return *neighbPointsPtr_;
541 }
542 
543 
545 {
546  if (!neighbEdgesPtr_)
547  {
549  << "No extended addressing calculated for patch " << name()
551  }
552  return *neighbEdgesPtr_;
553 }
554 
555 
557 (
558  PstreamBuffers& pBufs,
559  const primitivePatch& pp
560 ) const
561 {
562  if
563  (
564  !Pstream::parRun()
565  || transform() == NOORDERING
566  )
567  {
568  return;
569  }
570 
571  if (debug)
572  {
573  fileName nm
574  (
575  boundaryMesh().mesh().time().path()
576  /name()+"_faces.obj"
577  );
578  Pout<< "processorPolyPatch::order : Writing my " << pp.size()
579  << " faces to OBJ file " << nm << endl;
580  writeOBJ(nm, pp, pp.points());
581 
582  // Calculate my face centres
583  const pointField& fc = pp.faceCentres();
584 
585  OFstream localStr
586  (
587  boundaryMesh().mesh().time().path()
588  /name() + "_localFaceCentres.obj"
589  );
590  Pout<< "processorPolyPatch::order : "
591  << "Dumping " << fc.size()
592  << " local faceCentres to " << localStr.name() << endl;
593 
594  forAll(fc, facei)
595  {
596  writeOBJ(localStr, fc[facei]);
597  }
598  }
599 
600  if (owner())
601  {
602  if (transform() == COINCIDENTFULLMATCH)
603  {
604  // Pass the patch points and faces across
605  UOPstream toNeighbour(neighbProcNo(), pBufs);
606  toNeighbour << pp.localPoints()
607  << pp.localFaces();
608  }
609  else
610  {
611  const pointField& ppPoints = pp.points();
612 
613  pointField anchors(getAnchorPoints(pp, ppPoints, transform()));
614 
615  // Get the average of the points of each face. This is needed in
616  // case the face centroid calculation is incorrect due to the face
617  // having a very high aspect ratio.
618  pointField facePointAverages(pp.size(), Zero);
619  forAll(pp, fI)
620  {
621  const labelList& facePoints = pp[fI];
622 
623  forAll(facePoints, pI)
624  {
625  facePointAverages[fI] += ppPoints[facePoints[pI]];
626  }
627 
628  facePointAverages[fI] /= facePoints.size();
629  }
630 
631  // Now send all info over to the neighbour
632  UOPstream toNeighbour(neighbProcNo(), pBufs);
633  toNeighbour << pp.faceCentres() << pp.faceNormals()
634  << anchors << facePointAverages;
635  }
636  }
637 }
638 
639 
641 (
642  const face& a,
643  const pointField& aPts,
644  const face& b,
645  const pointField& bPts,
646  const bool sameOrientation,
647  const scalar absTolSqr,
648  scalar& matchDistSqr
649 )
650 {
651  if (a.size() != b.size())
652  {
653  return -1;
654  }
655 
656  enum CirculatorBase::direction circulateDirection
658 
659  if (!sameOrientation)
660  {
661  circulateDirection = CirculatorBase::ANTICLOCKWISE;
662  }
663 
664  label matchFp = -1;
665 
666  scalar closestMatchDistSqr = sqr(GREAT);
667 
668  ConstCirculator<face> aCirc(a);
669  ConstCirculator<face> bCirc(b);
670 
671  do
672  {
673  const scalar diffSqr = magSqr(aPts[*aCirc] - bPts[*bCirc]);
674 
675  if (diffSqr < absTolSqr)
676  {
677  // Found a matching point. Set the fulcrum of b to the iterator
678  ConstCirculator<face> bCirc2(bCirc);
679  ++aCirc;
680 
681  bCirc2.setFulcrumToIterator();
682 
683  if (!sameOrientation)
684  {
685  --bCirc2;
686  }
687  else
688  {
689  ++bCirc2;
690  }
691 
692  matchDistSqr = diffSqr;
693 
694  do
695  {
696  const scalar diffSqr2 = magSqr(aPts[*aCirc] - bPts[*bCirc2]);
697 
698  if (diffSqr2 > absTolSqr)
699  {
700  // No match.
701  break;
702  }
703 
704  matchDistSqr += diffSqr2;
705  }
706  while
707  (
708  aCirc.circulate(CirculatorBase::CLOCKWISE),
709  bCirc2.circulate(circulateDirection)
710  );
711 
712  if (!aCirc.circulate())
713  {
714  if (matchDistSqr < closestMatchDistSqr)
715  {
716  closestMatchDistSqr = matchDistSqr;
717 
718  if (!sameOrientation)
719  {
720  matchFp = a.size() - bCirc.nRotations();
721  }
722  else
723  {
724  matchFp = bCirc.nRotations();
725  }
726 
727  if (closestMatchDistSqr == 0)
728  {
729  break;
730  }
731  }
732  }
733 
734  // Reset aCirc
735  aCirc.setIteratorToFulcrum();
736  }
737 
738  } while (bCirc.circulate(circulateDirection));
739 
740  matchDistSqr = closestMatchDistSqr;
741 
742  return matchFp;
743 }
744 
745 
747 (
748  PstreamBuffers& pBufs,
749  const primitivePatch& pp,
751  labelList& rotation
752 ) const
753 {
754  // Note: we only get the faces that originate from internal faces.
755 
756  if
757  (
758  !Pstream::parRun()
759  || transform() == NOORDERING
760  )
761  {
762  return false;
763  }
764 
765  faceMap.setSize(pp.size());
766  faceMap = -1;
767 
768  rotation.setSize(pp.size());
769  rotation = 0;
770 
771  bool change = false;
772 
773  if (owner())
774  {
775  // Do nothing (i.e. identical mapping, zero rotation).
776  // See comment at top.
777  forAll(faceMap, patchFacei)
778  {
779  faceMap[patchFacei] = patchFacei;
780  }
781 
782  if (transform() != COINCIDENTFULLMATCH)
783  {
784  const pointField& ppPoints = pp.points();
785 
786  pointField anchors(getAnchorPoints(pp, ppPoints, transform()));
787 
788  // Calculate typical distance from face centre
789  scalarField tols
790  (
791  matchTolerance()*calcFaceTol(pp, pp.points(), pp.faceCentres())
792  );
793 
794  forAll(faceMap, patchFacei)
795  {
796  const point& wantedAnchor = anchors[patchFacei];
797 
798  rotation[patchFacei] = getRotation
799  (
800  ppPoints,
801  pp[patchFacei],
802  wantedAnchor,
803  tols[patchFacei]
804  );
805 
806  if (rotation[patchFacei] > 0)
807  {
808  change = true;
809  }
810  }
811  }
812 
813  return change;
814  }
815  else
816  {
817  // Calculate the absolute matching tolerance
818  scalarField tols
819  (
820  matchTolerance()*calcFaceTol(pp, pp.points(), pp.faceCentres())
821  );
822 
823  if (transform() == COINCIDENTFULLMATCH)
824  {
825  vectorField masterPts;
826  faceList masterFaces;
827 
828  {
829  UIPstream fromNeighbour(neighbProcNo(), pBufs);
830  fromNeighbour >> masterPts >> masterFaces;
831  }
832 
833  const pointField& localPts = pp.localPoints();
834  const faceList& localFaces = pp.localFaces();
835 
836  label nMatches = 0;
837 
838  forAll(pp, lFacei)
839  {
840  const face& localFace = localFaces[lFacei];
841  label faceRotation = -1;
842 
843  const scalar absTolSqr = sqr(tols[lFacei]);
844 
845  scalar closestMatchDistSqr = sqr(GREAT);
846  scalar matchDistSqr = sqr(GREAT);
847  label closestFaceMatch = -1;
848  label closestFaceRotation = -1;
849 
850  forAll(masterFaces, mFacei)
851  {
852  const face& masterFace = masterFaces[mFacei];
853 
854  faceRotation = matchFace
855  (
856  localFace,
857  localPts,
858  masterFace,
859  masterPts,
860  false,
861  absTolSqr,
862  matchDistSqr
863  );
864 
865  if
866  (
867  faceRotation != -1
868  && matchDistSqr < closestMatchDistSqr
869  )
870  {
871  closestMatchDistSqr = matchDistSqr;
872  closestFaceMatch = mFacei;
873  closestFaceRotation = faceRotation;
874  }
875 
876  if (closestMatchDistSqr == 0)
877  {
878  break;
879  }
880  }
881 
882  if
883  (
884  closestFaceRotation != -1
885  && closestMatchDistSqr < absTolSqr
886  )
887  {
888  faceMap[lFacei] = closestFaceMatch;
889 
890  rotation[lFacei] = closestFaceRotation;
891 
892  if (lFacei != closestFaceMatch || closestFaceRotation > 0)
893  {
894  change = true;
895  }
896 
897  nMatches++;
898  }
899  else
900  {
901  Pout<< "Number of matches = " << nMatches << " / "
902  << pp.size() << endl;
903 
904  pointField pts(localFace.size());
905  forAll(localFace, pI)
906  {
907  const label localPtI = localFace[pI];
908  pts[pI] = localPts[localPtI];
909  }
910 
912  << "No match for face " << localFace << nl << pts
913  << abort(FatalError);
914  }
915  }
916 
917  return change;
918  }
919  else
920  {
921  vectorField masterCtrs;
922  vectorField masterNormals;
923  vectorField masterAnchors;
924  vectorField masterFacePointAverages;
925 
926  // Receive data from neighbour
927  {
928  UIPstream fromNeighbour(neighbProcNo(), pBufs);
929  fromNeighbour >> masterCtrs >> masterNormals
930  >> masterAnchors >> masterFacePointAverages;
931  }
932 
933  if (debug || masterCtrs.size() != pp.size())
934  {
935  {
936  OFstream nbrStr
937  (
938  boundaryMesh().mesh().time().path()
939  /name() + "_nbrFaceCentres.obj"
940  );
941  Pout<< "processorPolyPatch::order : "
942  << "Dumping neighbour faceCentres to " << nbrStr.name()
943  << endl;
944  forAll(masterCtrs, facei)
945  {
946  writeOBJ(nbrStr, masterCtrs[facei]);
947  }
948  }
949 
950  if (masterCtrs.size() != pp.size())
951  {
953  << "in patch:" << name() << " : "
954  << "Local size of patch is " << pp.size() << " (faces)."
955  << endl
956  << "Received from neighbour " << masterCtrs.size()
957  << " faceCentres!"
958  << abort(FatalError);
959  }
960  }
961 
962  // Geometric match of face centre vectors
963  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
964 
965  // 1. Try existing ordering and transformation
966  bool matchedAll = matchPoints
967  (
968  pp.faceCentres(),
969  masterCtrs,
970  pp.faceNormals(),
971  masterNormals,
972  tols,
973  false,
974  faceMap
975  );
976 
977  // Fallback: try using face point average for matching
978  if (!matchedAll)
979  {
980  const pointField& ppPoints = pp.points();
981 
982  pointField facePointAverages(pp.size(), Zero);
983  forAll(pp, fI)
984  {
985  const labelList& facePoints = pp[fI];
986 
987  forAll(facePoints, pI)
988  {
989  facePointAverages[fI] += ppPoints[facePoints[pI]];
990  }
991 
992  facePointAverages[fI] /= facePoints.size();
993  }
994 
995  scalarField tols2
996  (
997  matchTolerance()
998  *calcFaceTol(pp, pp.points(), facePointAverages)
999  );
1000 
1001  // Note that we do not use the faceNormals anymore for
1002  // comparison. Since we're
1003  // having problems with the face centres (e.g. due to extreme
1004  // aspect ratios) we will probably also have problems with
1005  // reliable normals calculation
1006  labelList faceMap2(faceMap.size(), -1);
1007  matchPoints
1008  (
1009  facePointAverages,
1010  masterFacePointAverages,
1011  tols2,
1012  true,
1013  faceMap2
1014  );
1015 
1016  matchedAll = true;
1017 
1018  forAll(faceMap, oldFacei)
1019  {
1020  if (faceMap[oldFacei] == -1)
1021  {
1022  faceMap[oldFacei] = faceMap2[oldFacei];
1023 
1024  if (faceMap[oldFacei] == -1)
1025  {
1026  matchedAll = false;
1027  }
1028  }
1029  }
1030  }
1031 
1032  if (!matchedAll || debug)
1033  {
1034  // Dump faces
1035  fileName str
1036  (
1037  boundaryMesh().mesh().time().path()
1038  /name() + "_faces.obj"
1039  );
1040  Pout<< "processorPolyPatch::order :"
1041  << " Writing faces to OBJ file " << str.name() << endl;
1042  writeOBJ(str, pp, pp.points());
1043 
1044  OFstream ccStr
1045  (
1046  boundaryMesh().mesh().time().path()
1047  /name() + "_faceCentresConnections.obj"
1048  );
1049 
1050  Pout<< "processorPolyPatch::order :"
1051  << " Dumping newly found match as lines between"
1052  << " corresponding face centres to OBJ file "
1053  << ccStr.name()
1054  << endl;
1055 
1056  label vertI = 0;
1057 
1058  forAll(pp.faceCentres(), facei)
1059  {
1060  label masterFacei = faceMap[facei];
1061 
1062  if (masterFacei != -1)
1063  {
1064  const point& c0 = masterCtrs[masterFacei];
1065  const point& c1 = pp.faceCentres()[facei];
1066  writeOBJ(ccStr, c0, c1, vertI);
1067  }
1068  }
1069  }
1070 
1071  if (!matchedAll)
1072  {
1074  << "in patch:" << name() << " : "
1075  << "Cannot match vectors to faces on both sides of patch"
1076  << endl
1077  << " masterCtrs[0]:" << masterCtrs[0] << endl
1078  << " ctrs[0]:" << pp.faceCentres()[0] << endl
1079  << " Check your topology changes or maybe you have"
1080  << " multiple separated (from cyclics) processor patches"
1081  << endl
1082  << " Continuing with incorrect face ordering from now on"
1083  << endl;
1084 
1085  return false;
1086  }
1087 
1088  // Set rotation.
1089  forAll(faceMap, oldFacei)
1090  {
1091  // The face f will be at newFacei (after morphing) and we want
1092  // its anchorPoint (= f[0]) to align with the anchorpoint for
1093  // the corresponding face on the other side.
1094 
1095  label newFacei = faceMap[oldFacei];
1096 
1097  const point& wantedAnchor = masterAnchors[newFacei];
1098 
1099  rotation[newFacei] = getRotation
1100  (
1101  pp.points(),
1102  pp[oldFacei],
1103  wantedAnchor,
1104  tols[oldFacei]
1105  );
1106 
1107  if (rotation[newFacei] == -1)
1108  {
1110  << "in patch " << name()
1111  << " : "
1112  << "Cannot find point on face " << pp[oldFacei]
1113  << " with vertices "
1114  << UIndirectList<point>(pp.points(), pp[oldFacei])
1115  << " that matches point " << wantedAnchor
1116  << " when matching the halves of processor patch "
1117  << name()
1118  << "Continuing with incorrect face ordering from now on"
1119  << endl;
1120 
1121  return false;
1122  }
1123  }
1124 
1125  forAll(faceMap, facei)
1126  {
1127  if (faceMap[facei] != facei || rotation[facei] != 0)
1128  {
1129  return true;
1130  }
1131  }
1133  return false;
1134  }
1135  }
1136 }
1137 
1138 
1139 void Foam::processorPolyPatch::write(Ostream& os) const
1140 {
1142  os.writeEntry("myProcNo", myProcNo_);
1143  os.writeEntry("neighbProcNo", neighbProcNo_);
1144 }
1145 
1146 
1147 // ************************************************************************* //
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:598
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: polyPatch.H:135
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:1049
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:128
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:316
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:395
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.
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