cyclicPolyPatch.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 "cyclicPolyPatch.H"
28 #include "polyBoundaryMesh.H"
29 #include "polyMesh.H"
30 #include "demandDrivenData.H"
31 #include "OFstream.H"
32 #include "matchPoints.H"
33 #include "EdgeMap.H"
34 #include "Time.H"
35 #include "transformField.H"
36 #include "SubField.H"
37 #include "unitConversion.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(cyclicPolyPatch, 0);
44 
45  addToRunTimeSelectionTable(polyPatch, cyclicPolyPatch, word);
46  addToRunTimeSelectionTable(polyPatch, cyclicPolyPatch, dictionary);
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 Foam::label Foam::cyclicPolyPatch::findMaxArea
53 (
54  const pointField& points,
55  const faceList& faces
56 )
57 {
58  label maxI = -1;
59  scalar maxAreaSqr = -GREAT;
60 
61  forAll(faces, facei)
62  {
63  scalar areaSqr = magSqr(faces[facei].normal(points));
64 
65  if (areaSqr > maxAreaSqr)
66  {
67  maxAreaSqr = areaSqr;
68  maxI = facei;
69  }
70  }
71  return maxI;
72 }
73 
74 
76 {
77  if (size())
78  {
79  // Half0
80  const cyclicPolyPatch& half0 = *this;
81  vectorField half0Areas(half0.size());
82  forAll(half0, facei)
83  {
84  half0Areas[facei] = half0[facei].normal(half0.points());
85  }
86 
87  // Half1
88  const cyclicPolyPatch& half1 = neighbPatch();
89  vectorField half1Areas(half1.size());
90  forAll(half1, facei)
91  {
92  half1Areas[facei] = half1[facei].normal(half1.points());
93  }
94 
95  calcTransforms
96  (
97  half0,
98  half0.faceCentres(),
99  half0Areas,
100  half1.faceCentres(),
101  half1Areas
102  );
103  }
104 }
105 
106 
108 (
109  const primitivePatch& half0,
110  const pointField& half0Ctrs,
111  const vectorField& half0Areas,
112  const pointField& half1Ctrs,
113  const vectorField& half1Areas
114 )
115 {
116  if (debug && owner())
117  {
118  fileName casePath(boundaryMesh().mesh().time().path());
119  {
120  fileName nm0(casePath/name()+"_faces.obj");
121  Pout<< "cyclicPolyPatch::calcTransforms : Writing " << name()
122  << " faces to OBJ file " << nm0 << endl;
123  writeOBJ(nm0, half0, half0.points());
124  }
125  const cyclicPolyPatch& half1 = neighbPatch();
126  {
127  fileName nm1(casePath/half1.name()+"_faces.obj");
128  Pout<< "cyclicPolyPatch::calcTransforms : Writing " << half1.name()
129  << " faces to OBJ file " << nm1 << endl;
130  writeOBJ(nm1, half1, half1.points());
131  }
132  {
133  OFstream str(casePath/name()+"_to_" + half1.name() + ".obj");
134  label vertI = 0;
135  Pout<< "cyclicPolyPatch::calcTransforms :"
136  << " Writing coupled face centres as lines to " << str.name()
137  << endl;
138  forAll(half0Ctrs, i)
139  {
140  const point& p0 = half0Ctrs[i];
141  str << "v " << p0.x() << ' ' << p0.y() << ' ' << p0.z() << nl;
142  vertI++;
143  const point& p1 = half1Ctrs[i];
144  str << "v " << p1.x() << ' ' << p1.y() << ' ' << p1.z() << nl;
145  vertI++;
146  str << "l " << vertI-1 << ' ' << vertI << nl;
147  }
148  }
149  }
150 
151 
152  // Some sanity checks
153 
154  if (half0Ctrs.size() != half1Ctrs.size())
155  {
157  << "For patch " << name()
158  << " there are " << half0Ctrs.size()
159  << " face centres, for the neighbour patch " << neighbPatch().name()
160  << " there are " << half1Ctrs.size()
161  << exit(FatalError);
162  }
163 
164  if (transform() != neighbPatch().transform())
165  {
167  << "Patch " << name()
168  << " has transform type " << transformTypeNames[transform()]
169  << ", neighbour patch " << neighbPatchName()
170  << " has transform type "
171  << neighbPatch().transformTypeNames[neighbPatch().transform()]
172  << exit(FatalError);
173  }
174 
175 
176  // Calculate transformation tensors
177 
178  if (half0Ctrs.size() > 0)
179  {
180  vectorField half0Normals(half0Areas.size());
181  vectorField half1Normals(half1Areas.size());
182 
183  scalar maxAreaDiff = -GREAT;
184  label maxAreaFacei = -1;
185 
186  forAll(half0, facei)
187  {
188  scalar magSf = mag(half0Areas[facei]);
189  scalar nbrMagSf = mag(half1Areas[facei]);
190  scalar avSf = (magSf + nbrMagSf)/2.0;
191 
192  if (magSf < ROOTVSMALL && nbrMagSf < ROOTVSMALL)
193  {
194  // Undetermined normal. Use dummy normal to force separation
195  // check. (note use of sqrt(VSMALL) since that is how mag
196  // scales)
197  half0Normals[facei] = point(1, 0, 0);
198  half1Normals[facei] = half0Normals[facei];
199  }
200  else
201  {
202  scalar areaDiff = mag(magSf - nbrMagSf)/avSf;
203 
204  if (areaDiff > maxAreaDiff)
205  {
206  maxAreaDiff = areaDiff;
207  maxAreaFacei = facei;
208  }
209 
210  if (areaDiff > matchTolerance())
211  {
213  << "face " << facei
214  << " area does not match neighbour by "
215  << 100*areaDiff
216  << "% -- possible face ordering problem." << endl
217  << "patch:" << name()
218  << " my area:" << magSf
219  << " neighbour area:" << nbrMagSf
220  << " matching tolerance:" << matchTolerance()
221  << endl
222  << "Mesh face:" << start()+facei
223  << " fc:" << half0Ctrs[facei]
224  << endl
225  << "Neighbour fc:" << half1Ctrs[facei]
226  << endl
227  << "If you are certain your matching is correct"
228  << " you can increase the 'matchTolerance' setting"
229  << " in the patch dictionary in the boundary file."
230  << endl
231  << "Rerun with cyclic debug flag set"
232  << " for more information." << exit(FatalError);
233  }
234  else
235  {
236  half0Normals[facei] = half0Areas[facei] / magSf;
237  half1Normals[facei] = half1Areas[facei] / nbrMagSf;
238  }
239  }
240  }
241 
242 
243  // Print area match
244  if (debug)
245  {
246  Pout<< "cyclicPolyPatch::calcTransforms :"
247  << " patch:" << name()
248  << " Max area error:" << 100*maxAreaDiff << "% at face:"
249  << maxAreaFacei << " at:" << half0Ctrs[maxAreaFacei]
250  << " coupled face at:" << half1Ctrs[maxAreaFacei]
251  << endl;
252  }
253 
254 
255  // Calculate transformation tensors
256 
257  if (transform() == ROTATIONAL)
258  {
259  // Calculate using the given rotation axis and centre. Do not
260  // use calculated normals.
261  vector n0 = findFaceMaxRadius(half0Ctrs);
262  vector n1 = -findFaceMaxRadius(half1Ctrs);
263  n0 /= mag(n0) + VSMALL;
264  n1 /= mag(n1) + VSMALL;
265 
266  if (debug)
267  {
268  scalar theta = radToDeg(acos(n0 & n1));
269 
270  Pout<< "cyclicPolyPatch::calcTransforms :"
271  << " patch:" << name()
272  << " Specified rotation :"
273  << " n0:" << n0 << " n1:" << n1
274  << " swept angle: " << theta << " [deg]"
275  << endl;
276  }
277 
278  // Extended tensor from two local coordinate systems calculated
279  // using normal and rotation axis
280  const tensor E0
281  (
282  rotationAxis_,
283  (n0 ^ rotationAxis_),
284  n0
285  );
286  const tensor E1
287  (
288  rotationAxis_,
289  (-n1 ^ rotationAxis_),
290  -n1
291  );
292  const tensor revT(E1.T() & E0);
293 
294  const_cast<tensorField&>(forwardT()) = tensorField(1, revT.T());
295  const_cast<tensorField&>(reverseT()) = tensorField(1, revT);
296  const_cast<vectorField&>(separation()).setSize(0);
297  const_cast<boolList&>(collocated()) = boolList(1, false);
298  }
299  else
300  {
301  scalarField half0Tols
302  (
303  matchTolerance()
304  *calcFaceTol
305  (
306  half0,
307  half0.points(),
308  static_cast<const pointField&>(half0Ctrs)
309  )
310  );
311 
312  calcTransformTensors
313  (
314  static_cast<const pointField&>(half0Ctrs),
315  static_cast<const pointField&>(half1Ctrs),
316  half0Normals,
317  half1Normals,
318  half0Tols,
319  matchTolerance(),
320  transform()
321  );
322 
323 
324  if (transform() == TRANSLATIONAL)
325  {
326  if (debug)
327  {
328  Pout<< "cyclicPolyPatch::calcTransforms :"
329  << " patch:" << name()
330  << " Specified separation vector : "
331  << separationVector_ << endl;
332  }
333 
334  // Check that separation vectors are same.
335  const scalar avgTol = average(half0Tols);
336  if
337  (
338  mag(separationVector_ + neighbPatch().separationVector_)
339  > avgTol
340  )
341  {
343  << "Specified separation vector " << separationVector_
344  << " differs by that of neighbouring patch "
345  << neighbPatch().separationVector_
346  << " by more than tolerance " << avgTol << endl
347  << "patch:" << name()
348  << " neighbour:" << neighbPatchName()
349  << endl;
350  }
351 
352 
353  // Override computed transform with specified.
354  if
355  (
356  separation().size() != 1
357  || mag(separation()[0] - separationVector_) > avgTol
358  )
359  {
361  << "Specified separationVector " << separationVector_
362  << " differs from computed separation vector "
363  << separation() << endl
364  << "This probably means your geometry is not consistent"
365  << " with the specified separation and might lead"
366  << " to problems." << endl
367  << "Continuing with specified separation vector "
368  << separationVector_ << endl
369  << "patch:" << name()
370  << " neighbour:" << neighbPatchName()
371  << endl;
372  }
373 
374  // Set tensors
375  const_cast<tensorField&>(forwardT()).clear();
376  const_cast<tensorField&>(reverseT()).clear();
377  const_cast<vectorField&>(separation()) = vectorField
378  (
379  1,
380  separationVector_
381  );
382  const_cast<boolList&>(collocated()) = boolList(1, false);
383  }
384  }
385  }
386 }
387 
388 
389 void Foam::cyclicPolyPatch::getCentresAndAnchors
390 (
391  const primitivePatch& pp0,
392  const primitivePatch& pp1,
393 
394  pointField& half0Ctrs,
395  pointField& half1Ctrs,
396  pointField& anchors0,
397  scalarField& tols
398 ) const
399 {
400  // Get geometric data on both halves.
401  half0Ctrs = pp0.faceCentres();
402  anchors0 = getAnchorPoints(pp0, pp0.points(), transform());
403  half1Ctrs = pp1.faceCentres();
404 
405  if (debug)
406  {
407  Pout<< "cyclicPolyPatch::getCentresAndAnchors :"
408  << " patch:" << name() << nl
409  << "half0 untransformed faceCentres (avg) : "
410  << gAverage(half0Ctrs) << nl
411  << "half1 untransformed faceCentres (avg) : "
412  << gAverage(half1Ctrs) << endl;
413  }
414 
415  if (half0Ctrs.size())
416  {
417  switch (transform())
418  {
419  case ROTATIONAL:
420  {
421  vector n0 = findFaceMaxRadius(half0Ctrs);
422  vector n1 = -findFaceMaxRadius(half1Ctrs);
423  n0 /= mag(n0) + VSMALL;
424  n1 /= mag(n1) + VSMALL;
425 
426  if (debug)
427  {
428  scalar theta = radToDeg(acos(n0 & n1));
429 
430  Pout<< "cyclicPolyPatch::getCentresAndAnchors :"
431  << " patch:" << name()
432  << " Specified rotation :"
433  << " n0:" << n0 << " n1:" << n1
434  << " swept angle: " << theta << " [deg]"
435  << endl;
436  }
437 
438  // Extended tensor from two local coordinate systems calculated
439  // using normal and rotation axis
440  const tensor E0
441  (
442  rotationAxis_,
443  (n0 ^ rotationAxis_),
444  n0
445  );
446  const tensor E1
447  (
448  rotationAxis_,
449  (-n1 ^ rotationAxis_),
450  -n1
451  );
452  const tensor revT(E1.T() & E0);
453 
454  // Rotation
455  forAll(half0Ctrs, facei)
456  {
457  half0Ctrs[facei] =
459  (
460  revT,
461  half0Ctrs[facei] - rotationCentre_
462  )
463  + rotationCentre_;
464  anchors0[facei] =
466  (
467  revT,
468  anchors0[facei] - rotationCentre_
469  )
470  + rotationCentre_;
471  }
472 
473  break;
474  }
475  case TRANSLATIONAL:
476  {
477  // Transform 0 points.
478 
479  if (debug)
480  {
481  Pout<< "cyclicPolyPatch::getCentresAndAnchors :"
482  << " patch:" << name()
483  << "Specified translation : " << separationVector_
484  << endl;
485  }
486 
487  // Note: getCentresAndAnchors gets called on the slave side
488  // so separationVector is owner-slave points.
489 
490  half0Ctrs -= separationVector_;
491  anchors0 -= separationVector_;
492  break;
493  }
494  default:
495  {
496  // Assumes that cyclic is rotational. This is also the initial
497  // condition for patches without faces.
498 
499  // Determine the face with max area on both halves. These
500  // two faces are used to determine the transformation tensors
501  label max0I = findMaxArea(pp0.points(), pp0);
502  vector n0 = pp0[max0I].normal(pp0.points());
503  n0 /= mag(n0) + VSMALL;
504 
505  label max1I = findMaxArea(pp1.points(), pp1);
506  vector n1 = pp1[max1I].normal(pp1.points());
507  n1 /= mag(n1) + VSMALL;
508 
509  if (mag(n0 & n1) < 1-matchTolerance())
510  {
511  if (debug)
512  {
513  Pout<< "cyclicPolyPatch::getCentresAndAnchors :"
514  << " patch:" << name()
515  << " Detected rotation :"
516  << " n0:" << n0 << " n1:" << n1 << endl;
517  }
518 
519  // Rotation (around origin)
520  const tensor revT(rotationTensor(n0, -n1));
521 
522  // Rotation
523  forAll(half0Ctrs, facei)
524  {
525  half0Ctrs[facei] = Foam::transform
526  (
527  revT,
528  half0Ctrs[facei]
529  );
530  anchors0[facei] = Foam::transform
531  (
532  revT,
533  anchors0[facei]
534  );
535  }
536  }
537  else
538  {
539  // Parallel translation. Get average of all used points.
540 
541  const point ctr0(sum(pp0.localPoints())/pp0.nPoints());
542  const point ctr1(sum(pp1.localPoints())/pp1.nPoints());
543 
544  if (debug)
545  {
546  Pout<< "cyclicPolyPatch::getCentresAndAnchors :"
547  << " patch:" << name()
548  << " Detected translation :"
549  << " n0:" << n0 << " n1:" << n1
550  << " ctr0:" << ctr0 << " ctr1:" << ctr1 << endl;
551  }
552 
553  half0Ctrs += ctr1 - ctr0;
554  anchors0 += ctr1 - ctr0;
555  }
556  break;
557  }
558  }
559  }
560 
561  // Calculate typical distance per face
562  tols = matchTolerance()*calcFaceTol(pp1, pp1.points(), half1Ctrs);
563 }
564 
565 
566 Foam::vector Foam::cyclicPolyPatch::findFaceMaxRadius
567 (
568  const pointField& faceCentres
569 ) const
570 {
571  // Determine a face furthest away from the axis
572 
573  const vectorField n((faceCentres - rotationCentre_) ^ rotationAxis_);
574 
575  const scalarField magRadSqr(magSqr(n));
576 
577  label facei = findMax(magRadSqr);
578 
579  if (debug)
580  {
581  Info<< "findFaceMaxRadius(const pointField&) : patch: " << name() << nl
582  << " rotFace = " << facei << nl
583  << " point = " << faceCentres[facei] << nl
584  << " distance = " << Foam::sqrt(magRadSqr[facei])
585  << endl;
586  }
587 
588  return n[facei];
589 }
590 
591 
592 // * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * * * //
593 
595 (
596  const word& name,
597  const label size,
598  const label start,
599  const label index,
600  const polyBoundaryMesh& bm,
601  const word& patchType,
602  const transformType transform
603 )
604 :
605  coupledPolyPatch(name, size, start, index, bm, patchType, transform),
606  neighbPatchName_(word::null),
607  neighbPatchID_(-1),
608  rotationAxis_(Zero),
609  rotationCentre_(Zero),
610  separationVector_(Zero),
611  coupledPointsPtr_(nullptr),
612  coupledEdgesPtr_(nullptr)
613 {
614  // Neighbour patch might not be valid yet so no transformation
615  // calculation possible.
616 }
617 
618 
620 (
621  const word& name,
622  const label size,
623  const label start,
624  const label index,
625  const polyBoundaryMesh& bm,
626  const word& neighbPatchName,
627  const transformType transform,
628  const vector& rotationAxis,
629  const point& rotationCentre,
630  const vector& separationVector
631 )
632 :
633  coupledPolyPatch(name, size, start, index, bm, typeName, transform),
634  neighbPatchName_(neighbPatchName),
635  neighbPatchID_(-1),
636  rotationAxis_(rotationAxis),
637  rotationCentre_(rotationCentre),
638  separationVector_(separationVector),
639  coupledPointsPtr_(nullptr),
640  coupledEdgesPtr_(nullptr)
641 {
642  // Neighbour patch might not be valid yet so no transformation
643  // calculation possible.
644 }
645 
646 
648 (
649  const word& name,
650  const dictionary& dict,
651  const label index,
652  const polyBoundaryMesh& bm,
653  const word& patchType
654 )
655 :
656  coupledPolyPatch(name, dict, index, bm, patchType),
657  neighbPatchName_(dict.lookupOrDefault("neighbourPatch", word::null)),
658  coupleGroup_(dict),
659  neighbPatchID_(-1),
660  rotationAxis_(Zero),
661  rotationCentre_(Zero),
662  separationVector_(Zero),
663  coupledPointsPtr_(nullptr),
664  coupledEdgesPtr_(nullptr)
665 {
666  if (neighbPatchName_ == word::null && !coupleGroup_.valid())
667  {
669  (
670  dict
671  ) << "No \"neighbourPatch\" provided." << endl
672  << "Is your mesh uptodate with split cyclics?" << endl
673  << "Run foamUpgradeCyclics to convert mesh and fields"
674  << " to split cyclics." << exit(FatalIOError);
675  }
676 
677  if (neighbPatchName_ == name)
678  {
680  << "Neighbour patch name " << neighbPatchName_
681  << " cannot be the same as this patch " << name
682  << exit(FatalIOError);
683  }
684 
685  switch (transform())
686  {
687  case ROTATIONAL:
688  {
689  dict.lookup("rotationAxis") >> rotationAxis_;
690  dict.lookup("rotationCentre") >> rotationCentre_;
691 
692  scalar magRot = mag(rotationAxis_);
693  if (magRot < SMALL)
694  {
696  << "Illegal rotationAxis " << rotationAxis_ << endl
697  << "Please supply a non-zero vector."
698  << exit(FatalIOError);
699  }
700  rotationAxis_ /= magRot;
701 
702  break;
703  }
704  case TRANSLATIONAL:
705  {
706  dict.lookup("separationVector") >> separationVector_;
707  break;
708  }
709  default:
710  {
711  // no additional info required
712  }
713  }
714 
715  // Neighbour patch might not be valid yet so no transformation
716  // calculation possible.
717 }
718 
719 
721 (
722  const cyclicPolyPatch& pp,
723  const polyBoundaryMesh& bm
724 )
725 :
726  coupledPolyPatch(pp, bm),
727  neighbPatchName_(pp.neighbPatchName_),
728  coupleGroup_(pp.coupleGroup_),
729  neighbPatchID_(-1),
730  rotationAxis_(pp.rotationAxis_),
731  rotationCentre_(pp.rotationCentre_),
732  separationVector_(pp.separationVector_),
733  coupledPointsPtr_(nullptr),
734  coupledEdgesPtr_(nullptr)
735 {
736  // Neighbour patch might not be valid yet so no transformation
737  // calculation possible.
738 }
739 
740 
742 (
743  const cyclicPolyPatch& pp,
744  const polyBoundaryMesh& bm,
745  const label index,
746  const label newSize,
747  const label newStart,
748  const word& neighbName
749 )
750 :
751  coupledPolyPatch(pp, bm, index, newSize, newStart),
752  neighbPatchName_(neighbName),
753  coupleGroup_(pp.coupleGroup_),
754  neighbPatchID_(-1),
755  rotationAxis_(pp.rotationAxis_),
756  rotationCentre_(pp.rotationCentre_),
757  separationVector_(pp.separationVector_),
758  coupledPointsPtr_(nullptr),
759  coupledEdgesPtr_(nullptr)
760 {
761  if (neighbName == name())
762  {
764  << "Neighbour patch name " << neighbName
765  << " cannot be the same as this patch " << name()
766  << exit(FatalError);
767  }
768 
769  // Neighbour patch might not be valid yet so no transformation
770  // calculation possible.
771 }
772 
773 
775 (
776  const cyclicPolyPatch& pp,
777  const polyBoundaryMesh& bm,
778  const label index,
779  const labelUList& mapAddressing,
780  const label newStart
781 )
782 :
783  coupledPolyPatch(pp, bm, index, mapAddressing, newStart),
784  neighbPatchName_(pp.neighbPatchName_),
785  coupleGroup_(pp.coupleGroup_),
786  neighbPatchID_(-1),
787  rotationAxis_(pp.rotationAxis_),
788  rotationCentre_(pp.rotationCentre_),
789  separationVector_(pp.separationVector_),
790  coupledPointsPtr_(nullptr),
791  coupledEdgesPtr_(nullptr)
792 {}
793 
794 
795 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
796 
798 {
799  deleteDemandDrivenData(coupledPointsPtr_);
800  deleteDemandDrivenData(coupledEdgesPtr_);
801 }
802 
803 
804 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
805 
807 {
808  if (neighbPatchName_.empty())
809  {
810  // Try and use patchGroup to find samplePatch and sampleRegion
811  label patchID = coupleGroup_.findOtherPatchID(*this);
812 
813  neighbPatchName_ = boundaryMesh()[patchID].name();
814  }
815  return neighbPatchName_;
816 }
817 
818 
820 {
821  if (neighbPatchID_ == -1)
822  {
823  neighbPatchID_ = this->boundaryMesh().findPatchID(neighbPatchName());
824 
825  if (neighbPatchID_ == -1)
826  {
828  << "Illegal neighbourPatch name " << neighbPatchName()
829  << endl << "Valid patch names are "
830  << this->boundaryMesh().names()
831  << exit(FatalError);
832  }
833 
834  // Check that it is a cyclic
835  const cyclicPolyPatch& nbrPatch = refCast<const cyclicPolyPatch>
836  (
837  this->boundaryMesh()[neighbPatchID_]
838  );
839 
840  if (nbrPatch.neighbPatchName() != name())
841  {
843  << "Patch " << name()
844  << " specifies neighbour patch " << neighbPatchName()
845  << endl << " but that in return specifies "
846  << nbrPatch.neighbPatchName()
847  << endl;
848  }
849  }
850  return neighbPatchID_;
851 }
852 
853 
855 {
856  if (!parallel())
857  {
858  if (transform() == ROTATIONAL)
859  {
860  l =
861  Foam::transform(forwardT(), l-rotationCentre_)
862  + rotationCentre_;
863  }
864  else
865  {
866  l = Foam::transform(forwardT(), l);
867  }
868  }
869  else if (separated())
870  {
871  // transformPosition gets called on the receiving side,
872  // separation gets calculated on the sending side so subtract.
873 
874  const vectorField& s = separation();
875  if (s.size() == 1)
876  {
877  forAll(l, i)
878  {
879  l[i] -= s[0];
880  }
881  }
882  else
883  {
884  l -= s;
885  }
886  }
887 }
888 
889 
891 {
892  if (!parallel())
893  {
894  const tensor& T =
895  (
896  forwardT().size() == 1
897  ? forwardT()[0]
898  : forwardT()[facei]
899  );
900 
901  if (transform() == ROTATIONAL)
902  {
903  l = Foam::transform(T, l-rotationCentre_) + rotationCentre_;
904  }
905  else
906  {
907  l = Foam::transform(T, l);
908  }
909  }
910  else if (separated())
911  {
912  const vector& s =
913  (
914  separation().size() == 1
915  ? separation()[0]
916  : separation()[facei]
917  );
918 
919  l -= s;
920  }
921 }
922 
923 
925 {
927 }
928 
929 
931 (
932  const primitivePatch& referPatch,
933  pointField& nbrCtrs,
934  vectorField& nbrAreas,
935  pointField& nbrCc
936 )
937 {}
938 
939 
941 (
942  const primitivePatch& referPatch,
943  const pointField& thisCtrs,
944  const vectorField& thisAreas,
945  const pointField& thisCc,
946  const pointField& nbrCtrs,
947  const vectorField& nbrAreas,
948  const pointField& nbrCc
949 )
950 {
951  calcTransforms
952  (
953  referPatch,
954  thisCtrs,
955  thisAreas,
956  nbrCtrs,
957  nbrAreas
958  );
959 }
960 
961 
963 {
964  calcGeometry
965  (
966  *this,
967  faceCentres(),
968  faceAreas(),
969  faceCellCentres(),
970  neighbPatch().faceCentres(),
971  neighbPatch().faceAreas(),
972  neighbPatch().faceCellCentres()
973  );
974 }
975 
976 
978 (
979  PstreamBuffers& pBufs,
980  const pointField& p
981 )
982 {
983  polyPatch::initMovePoints(pBufs, p);
984 }
985 
986 
988 (
989  PstreamBuffers& pBufs,
990  const pointField& p
991 )
992 {
993  polyPatch::movePoints(pBufs, p);
994  calcTransforms();
995 }
996 
997 
999 {
1001 }
1002 
1003 
1005 {
1006  polyPatch::updateMesh(pBufs);
1007  deleteDemandDrivenData(coupledPointsPtr_);
1008  deleteDemandDrivenData(coupledEdgesPtr_);
1009 }
1010 
1011 
1013 {
1014  if (!coupledPointsPtr_)
1015  {
1016  const faceList& nbrLocalFaces = neighbPatch().localFaces();
1017  const labelList& nbrMeshPoints = neighbPatch().meshPoints();
1018 
1019  // Now all we know is that relative face index in *this is same
1020  // as coupled face in nbrPatch and also that the 0th vertex
1021  // corresponds.
1022 
1023  // From local point to nbrPatch or -1.
1024  labelList coupledPoint(nPoints(), -1);
1025 
1026  forAll(*this, patchFacei)
1027  {
1028  const face& fA = localFaces()[patchFacei];
1029  const face& fB = nbrLocalFaces[patchFacei];
1030 
1031  forAll(fA, indexA)
1032  {
1033  label patchPointA = fA[indexA];
1034 
1035  if (coupledPoint[patchPointA] == -1)
1036  {
1037  label indexB = (fB.size() - indexA) % fB.size();
1038 
1039  // Filter out points on wedge axis
1040  if (meshPoints()[patchPointA] != nbrMeshPoints[fB[indexB]])
1041  {
1042  coupledPoint[patchPointA] = fB[indexB];
1043  }
1044  }
1045  }
1046  }
1047 
1048  coupledPointsPtr_ = new edgeList(nPoints());
1049  edgeList& connected = *coupledPointsPtr_;
1050 
1051  // Extract coupled points.
1052  label connectedI = 0;
1053 
1054  forAll(coupledPoint, i)
1055  {
1056  if (coupledPoint[i] != -1)
1057  {
1058  connected[connectedI++] = edge(i, coupledPoint[i]);
1059  }
1060  }
1061 
1062  connected.setSize(connectedI);
1063 
1064  if (debug)
1065  {
1066  OFstream str
1067  (
1068  boundaryMesh().mesh().time().path()
1069  /name() + "_coupledPoints.obj"
1070  );
1071  label vertI = 0;
1072 
1073  Pout<< "Writing file " << str.name() << " with coordinates of "
1074  << "coupled points" << endl;
1075 
1076  forAll(connected, i)
1077  {
1078  const point& a = points()[meshPoints()[connected[i][0]]];
1079  const point& b = points()[nbrMeshPoints[connected[i][1]]];
1080 
1081  str<< "v " << a.x() << ' ' << a.y() << ' ' << a.z() << nl;
1082  str<< "v " << b.x() << ' ' << b.y() << ' ' << b.z() << nl;
1083  vertI += 2;
1084 
1085  str<< "l " << vertI-1 << ' ' << vertI << nl;
1086  }
1087  }
1088  }
1089  return *coupledPointsPtr_;
1090 }
1091 
1092 
1094 {
1095  if (!coupledEdgesPtr_)
1096  {
1097  const edgeList& pointCouples = coupledPoints();
1098 
1099  // Build map from points on *this (A) to points on neighbourpatch (B)
1100  Map<label> aToB(2*pointCouples.size());
1101 
1102  forAll(pointCouples, i)
1103  {
1104  const edge& e = pointCouples[i];
1105 
1106  aToB.insert(e[0], e[1]);
1107  }
1108 
1109  // Map from edge on A to points (in B indices)
1110  EdgeMap<label> edgeMap(nEdges());
1111 
1112  forAll(*this, patchFacei)
1113  {
1114  const labelList& fEdges = faceEdges()[patchFacei];
1115 
1116  forAll(fEdges, i)
1117  {
1118  label edgeI = fEdges[i];
1119 
1120  const edge& e = edges()[edgeI];
1121 
1122  // Convert edge end points to corresponding points on B side.
1123  Map<label>::const_iterator fnd0 = aToB.find(e[0]);
1124  if (fnd0 != aToB.end())
1125  {
1126  Map<label>::const_iterator fnd1 = aToB.find(e[1]);
1127  if (fnd1 != aToB.end())
1128  {
1129  edgeMap.insert(edge(fnd0(), fnd1()), edgeI);
1130  }
1131  }
1132  }
1133  }
1134 
1135  // Use the edgeMap to get the edges on the B side.
1136 
1137  const cyclicPolyPatch& neighbPatch = this->neighbPatch();
1138  const labelList& nbrMp = neighbPatch.meshPoints();
1139  const labelList& mp = meshPoints();
1140 
1141 
1142 
1143  coupledEdgesPtr_ = new edgeList(edgeMap.size());
1144  edgeList& coupledEdges = *coupledEdgesPtr_;
1145  label coupleI = 0;
1146 
1147  forAll(neighbPatch, patchFacei)
1148  {
1149  const labelList& fEdges = neighbPatch.faceEdges()[patchFacei];
1150 
1151  forAll(fEdges, i)
1152  {
1153  label edgeI = fEdges[i];
1154 
1155  const edge& e = neighbPatch.edges()[edgeI];
1156 
1157  // Look up A edge from HashTable.
1158  EdgeMap<label>::iterator iter = edgeMap.find(e);
1159 
1160  if (iter != edgeMap.end())
1161  {
1162  label edgeA = iter();
1163  const edge& eA = edges()[edgeA];
1164 
1165  // Store correspondence. Filter out edges on wedge axis.
1166  if
1167  (
1168  edge(mp[eA[0]], mp[eA[1]])
1169  != edge(nbrMp[e[0]], nbrMp[e[1]])
1170  )
1171  {
1172  coupledEdges[coupleI++] = edge(edgeA, edgeI);
1173  }
1174 
1175  // Remove so we build unique list
1176  edgeMap.erase(iter);
1177  }
1178  }
1179  }
1180  coupledEdges.setSize(coupleI);
1181 
1182 
1183  // Some checks
1184 
1185  forAll(coupledEdges, i)
1186  {
1187  const edge& e = coupledEdges[i];
1188 
1189  if (e[0] < 0 || e[1] < 0)
1190  {
1192  << "Problem : at position " << i
1193  << " illegal couple:" << e
1194  << abort(FatalError);
1195  }
1196  }
1197 
1198  if (debug)
1199  {
1200  OFstream str
1201  (
1202  boundaryMesh().mesh().time().path()
1203  /name() + "_coupledEdges.obj"
1204  );
1205  label vertI = 0;
1206 
1207  Pout<< "Writing file " << str.name() << " with centres of "
1208  << "coupled edges" << endl;
1209 
1210  forAll(coupledEdges, i)
1211  {
1212  const edge& e = coupledEdges[i];
1213 
1214  const point& a = edges()[e[0]].centre(localPoints());
1215  const point& b = neighbPatch.edges()[e[1]].centre
1216  (
1217  neighbPatch.localPoints()
1218  );
1219 
1220  str<< "v " << a.x() << ' ' << a.y() << ' ' << a.z() << nl;
1221  str<< "v " << b.x() << ' ' << b.y() << ' ' << b.z() << nl;
1222  vertI += 2;
1223 
1224  str<< "l " << vertI-1 << ' ' << vertI << nl;
1225  }
1226  }
1227  }
1228  return *coupledEdgesPtr_;
1229 }
1230 
1231 
1234  PstreamBuffers&,
1235  const primitivePatch& pp
1236 ) const
1237 {
1238  if (owner())
1239  {
1240  // Save patch for use in non-owner side ordering. Equivalent to
1241  // processorPolyPatch using OPstream.
1242  ownerPatchPtr_.reset
1243  (
1244  new primitivePatch
1245  (
1246  pp,
1247  pp.points()
1248  )
1249  );
1250  }
1251 }
1252 
1253 
1256  PstreamBuffers& pBufs,
1257  const primitivePatch& pp,
1258  labelList& faceMap,
1259  labelList& rotation
1260 ) const
1261 {
1262  if (debug)
1263  {
1264  Pout<< "order : of " << pp.size()
1265  << " faces of patch:" << name()
1266  << " neighbour:" << neighbPatchName()
1267  << endl;
1268  }
1269  faceMap.setSize(pp.size());
1270  faceMap = -1;
1271 
1272  rotation.setSize(pp.size());
1273  rotation = 0;
1274 
1275  if (transform() == NOORDERING)
1276  {
1277  // No faces, nothing to change.
1278  return false;
1279  }
1280 
1281  if (owner())
1282  {
1283  // Do nothing (i.e. identical mapping, zero rotation).
1284  // See comment at top.
1285  forAll(faceMap, patchFacei)
1286  {
1287  faceMap[patchFacei] = patchFacei;
1288  }
1289 
1290  return false;
1291  }
1292  else
1293  {
1294  // Get stored geometry from initOrder invocation of owner.
1295  const primitivePatch& pp0 = neighbPatch().ownerPatchPtr_();
1296 
1297  // Get geometric quantities
1298  pointField half0Ctrs, half1Ctrs, anchors0;
1299  scalarField tols;
1300  getCentresAndAnchors
1301  (
1302  pp0,
1303  pp,
1304 
1305  half0Ctrs,
1306  half1Ctrs,
1307  anchors0,
1308  tols
1309  );
1310 
1311  if (debug)
1312  {
1313  Pout<< "half0 transformed faceCentres (avg) : "
1314  << gAverage(half0Ctrs) << nl
1315  << "half1 untransformed faceCentres (avg) : "
1316  << gAverage(half1Ctrs) << endl;
1317  }
1318 
1319  // Geometric match of face centre vectors
1320  bool matchedAll = matchPoints
1321  (
1322  half1Ctrs,
1323  half0Ctrs,
1324  tols,
1325  true,
1326  faceMap
1327  );
1328 
1329  if (!matchedAll || debug)
1330  {
1331  // Dump halves
1332  fileName nm0
1333  (
1334  boundaryMesh().mesh().time().path()
1335  /neighbPatch().name()+"_faces.obj"
1336  );
1337  Pout<< "cyclicPolyPatch::order : Writing neighbour"
1338  << " faces to OBJ file " << nm0 << endl;
1339  writeOBJ(nm0, pp0, pp0.points());
1340 
1341  fileName nm1
1342  (
1343  boundaryMesh().mesh().time().path()
1344  /name()+"_faces.obj"
1345  );
1346  Pout<< "cyclicPolyPatch::order : Writing my"
1347  << " faces to OBJ file " << nm1 << endl;
1348  writeOBJ(nm1, pp, pp.points());
1349 
1350  OFstream ccStr
1351  (
1352  boundaryMesh().mesh().time().path()
1353  /name() + "_faceCentres.obj"
1354  );
1355  Pout<< "cyclicPolyPatch::order : "
1356  << "Dumping currently found cyclic match as lines between"
1357  << " corresponding face centres to file " << ccStr.name()
1358  << endl;
1359 
1360  // Recalculate untransformed face centres
1361  //pointField rawHalf0Ctrs =
1362  // calcFaceCentres(half0Faces, pp.points());
1363  label vertI = 0;
1364 
1365  forAll(half1Ctrs, i)
1366  {
1367  if (faceMap[i] != -1)
1368  {
1369  // Write edge between c1 and c0
1370  const point& c0 = half0Ctrs[faceMap[i]];
1371  const point& c1 = half1Ctrs[i];
1372  writeOBJ(ccStr, c0, c1, vertI);
1373  }
1374  }
1375  }
1376 
1377  if (!matchedAll)
1378  {
1380  << "Patch:" << name() << " : "
1381  << "Cannot match vectors to faces on both sides of patch"
1382  << endl
1383  << " Perhaps your faces do not match?"
1384  << " The obj files written contain the current match." << endl
1385  << " Continuing with incorrect face ordering from now on!"
1386  << endl;
1387 
1388  return false;
1389  }
1390 
1391 
1392  // Set rotation.
1393  forAll(faceMap, oldFacei)
1394  {
1395  // The face f will be at newFacei (after morphing) and we want its
1396  // anchorPoint (= f[0]) to align with the anchorpoint for the
1397  // corresponding face on the other side.
1398 
1399  label newFacei = faceMap[oldFacei];
1400 
1401  const point& wantedAnchor = anchors0[newFacei];
1402 
1403  rotation[newFacei] = getRotation
1404  (
1405  pp.points(),
1406  pp[oldFacei],
1407  wantedAnchor,
1408  tols[oldFacei]
1409  );
1410 
1411  if (rotation[newFacei] == -1)
1412  {
1414  << "in patch " << name()
1415  << " : "
1416  << "Cannot find point on face " << pp[oldFacei]
1417  << " with vertices "
1418  << IndirectList<point>(pp.points(), pp[oldFacei])()
1419  << " that matches point " << wantedAnchor
1420  << " when matching the halves of processor patch " << name()
1421  << "Continuing with incorrect face ordering from now on!"
1422  << endl;
1423 
1424  return false;
1425  }
1426  }
1427 
1428  ownerPatchPtr_.clear();
1429 
1430  // Return false if no change neccesary, true otherwise.
1431 
1432  forAll(faceMap, facei)
1433  {
1434  if (faceMap[facei] != facei || rotation[facei] != 0)
1435  {
1436  return true;
1437  }
1438  }
1439 
1440  return false;
1441  }
1442 }
1443 
1444 
1446 {
1448  if (!neighbPatchName_.empty())
1449  {
1450  os.writeKeyword("neighbourPatch") << neighbPatchName_
1451  << token::END_STATEMENT << nl;
1452  }
1453  coupleGroup_.write(os);
1454  switch (transform())
1455  {
1456  case ROTATIONAL:
1457  {
1458  os.writeKeyword("rotationAxis") << rotationAxis_
1459  << token::END_STATEMENT << nl;
1460  os.writeKeyword("rotationCentre") << rotationCentre_
1461  << token::END_STATEMENT << nl;
1462  break;
1463  }
1464  case TRANSLATIONAL:
1465  {
1466  os.writeKeyword("separationVector") << separationVector_
1467  << token::END_STATEMENT << nl;
1468  break;
1469  }
1470  case NOORDERING:
1471  {
1472  break;
1473  }
1474  default:
1475  {
1476  // no additional info to write
1477  }
1478  }
1479 }
1480 
1481 
1482 // ************************************************************************* //
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
label nPoints() const
Return number of points supporting patch faces.
static bool valid(char)
Is this character valid for a word.
Definition: wordI.H:115
dimensionedScalar acos(const dimensionedScalar &ds)
#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
scalar radToDeg(const scalar rad)
Conversion from radians to degrees.
A class for handling file names.
Definition: fileName.H:69
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const edgeList & coupledPoints() const
Return connected points (from patch local to neighbour patch local)
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:110
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
const double e
Elementary charge.
Definition: doubleFloat.H:78
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
tUEqn clear()
Unit conversion functions.
Output to file stream.
Definition: OFstream.H:82
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Tensor< Cmpt > T() const
Return transpose.
Definition: TensorI.H:321
List< face > faceList
Definition: faceListFwd.H:43
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:59
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
cyclicPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform=UNKNOWN)
Construct from components.
const Field< PointType > & faceCentres() const
Return face centres for patch.
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
const Cmpt & z() const
Definition: VectorI.H:87
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:47
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
Definition: polyPatch.H:108
const edgeList & coupledEdges() const
Return connected edges (from patch local to neighbour patch local).
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:120
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
Determine correspondence between points. See below.
Macros for easy insertion into run-time selection tables.
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
const Cmpt & y() const
Definition: VectorI.H:81
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
bool erase(const iterator &)
Erase a hashedEntry specified by given iterator.
Definition: HashTable.C:371
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.
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
const Vector< Cmpt > & centre(const Foam::List< Vector< Cmpt >> &) const
Return *this (used for point which is a typedef to Vector<scalar>.
Definition: VectorI.H:116
points setSize(newPointi)
Spatial transformation functions for primitive fields.
A list of faces which address into the list of points.
bool matchPoints(const UList< point > &pts0, const UList< point > &pts1, const UList< scalar > &matchDistance, const bool verbose, labelList &from0To1, const point &origin=point::zero)
Determine correspondence between pointFields. Gets passed.
Definition: matchPoints.C:33
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dynamicFvMesh & mesh
const pointField & points
List< edge > edgeList
Definition: edgeList.H:38
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
Cyclic plane patch.
word name() const
Return file name (part beyond last /)
Definition: fileName.C:180
const Field< PointType > & points() const
Return reference to global points.
static const word null
An empty word.
Definition: word.H:77
virtual void movePoints(PstreamBuffers &, const pointField &)
Correct patches after moving points.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
static const zero Zero
Definition: zero.H:91
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: polyPatch.C:62
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
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:61
const Cmpt & x() const
Definition: VectorI.H:75
Foam::polyBoundaryMesh.
label findMax(const ListType &, const label start=0)
Find index of max element (and larger than given element).
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
const Field< PointType > & localPoints() const
Return pointField of points in patch.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:262
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
defineTypeNameAndDebug(combustionModel, 0)
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: polyPatch.H:100
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
void setSize(const label)
Reset size of List.
Definition: List.C:281
Template functions to aid in the implementation of demand driven data.
virtual void transformPosition(pointField &l) const
Transform a patch-based position from other side to this side.
virtual ~cyclicPolyPatch()
Destructor.
vector point
Point is a vector.
Definition: point.H:41
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
Type gAverage(const FieldField< Field, Type > &f)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
virtual label neighbPatchID() const
Neighbour patchID.
const word & neighbPatchName() const
Neighbour patch name.
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
messageStream Info
const labelListList & faceEdges() const
Return face-edge addressing.
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
void deleteDemandDrivenData(DataPtr &dataPtr)
A List with indirect addressing.
Definition: IndirectList.H:102
const dimensionedScalar mp
Proton mass.
virtual void calcTransforms()
Recalculate the transformation tensors.
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
IOerror FatalIOError
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:284