slidingInterface.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 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 "slidingInterface.H"
27 #include "polyTopoChanger.H"
28 #include "polyMesh.H"
29 #include "polyTopoChange.H"
31 #include "plane.H"
32 
33 // Index of debug signs:
34 // p - adjusting a projection point
35 // * - adjusting edge intersection
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(slidingInterface, 0);
43  (
44  polyMeshModifier,
45  slidingInterface,
46  dictionary
47  );
48 
49  template<>
50  const char* Foam::NamedEnum
51  <
53  2
54  >::names[] =
55  {
56  "integral",
57  "partial"
58  };
59 }
60 
61 
64 
65 
66 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
67 
68 void Foam::slidingInterface::checkDefinition()
69 {
70  const polyMesh& mesh = topoChanger().mesh();
71 
72  if
73  (
74  !masterFaceZoneID_.active()
75  || !slaveFaceZoneID_.active()
76  || !cutPointZoneID_.active()
77  || !cutFaceZoneID_.active()
78  || !masterPatchID_.active()
79  || !slavePatchID_.active()
80  )
81  {
83  (
84  "void slidingInterface::checkDefinition()"
85  ) << "Not all zones and patches needed in the definition "
86  << "have been found. Please check your mesh definition."
87  << abort(FatalError);
88  }
89 
90  // Check the sizes and set up state
91  if
92  (
93  mesh.faceZones()[masterFaceZoneID_.index()].empty()
94  || mesh.faceZones()[slaveFaceZoneID_.index()].empty()
95  )
96  {
97  FatalErrorIn("void slidingInterface::checkDefinition()")
98  << "Master or slave face zone contain no faces. "
99  << "Please check your mesh definition."
100  << abort(FatalError);
101  }
102 
103  if (debug)
104  {
105  Pout<< "Sliding interface object " << name() << " :" << nl
106  << " master face zone: " << masterFaceZoneID_.index() << nl
107  << " slave face zone: " << slaveFaceZoneID_.index() << endl;
108  }
109 }
110 
111 
112 void Foam::slidingInterface::clearOut() const
113 {
114  clearPointProjection();
115  clearAttachedAddressing();
116  clearAddressing();
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
122 
123 // Construct from components
124 Foam::slidingInterface::slidingInterface
125 (
126  const word& name,
127  const label index,
128  const polyTopoChanger& mme,
129  const word& masterFaceZoneName,
130  const word& slaveFaceZoneName,
131  const word& cutPointZoneName,
132  const word& cutFaceZoneName,
133  const word& masterPatchName,
134  const word& slavePatchName,
135  const typeOfMatch tom,
136  const bool coupleDecouple,
137  const intersection::algorithm algo
138 )
139 :
140  polyMeshModifier(name, index, mme, true),
141  masterFaceZoneID_
142  (
143  masterFaceZoneName,
144  mme.mesh().faceZones()
145  ),
146  slaveFaceZoneID_
147  (
148  slaveFaceZoneName,
149  mme.mesh().faceZones()
150  ),
151  cutPointZoneID_
152  (
153  cutPointZoneName,
154  mme.mesh().pointZones()
155  ),
156  cutFaceZoneID_
157  (
158  cutFaceZoneName,
159  mme.mesh().faceZones()
160  ),
161  masterPatchID_
162  (
163  masterPatchName,
164  mme.mesh().boundaryMesh()
165  ),
166  slavePatchID_
167  (
168  slavePatchName,
169  mme.mesh().boundaryMesh()
170  ),
171  matchType_(tom),
172  coupleDecouple_(coupleDecouple),
173  attached_(false),
174  projectionAlgo_(algo),
175  trigger_(false),
176  pointMergeTol_(pointMergeTolDefault_),
177  edgeMergeTol_(edgeMergeTolDefault_),
178  nFacesPerSlaveEdge_(nFacesPerSlaveEdgeDefault_),
179  edgeFaceEscapeLimit_(edgeFaceEscapeLimitDefault_),
180  integralAdjTol_(integralAdjTolDefault_),
181  edgeMasterCatchFraction_(edgeMasterCatchFractionDefault_),
182  edgeCoPlanarTol_(edgeCoPlanarTolDefault_),
183  edgeEndCutoffTol_(edgeEndCutoffTolDefault_),
184  cutFaceMasterPtr_(NULL),
185  cutFaceSlavePtr_(NULL),
186  masterFaceCellsPtr_(NULL),
187  slaveFaceCellsPtr_(NULL),
188  masterStickOutFacesPtr_(NULL),
189  slaveStickOutFacesPtr_(NULL),
190  retiredPointMapPtr_(NULL),
191  cutPointEdgePairMapPtr_(NULL),
192  slavePointPointHitsPtr_(NULL),
193  slavePointEdgeHitsPtr_(NULL),
194  slavePointFaceHitsPtr_(NULL),
195  masterPointEdgeHitsPtr_(NULL),
196  projectedSlavePointsPtr_(NULL)
197 {
198  checkDefinition();
199 
200  if (attached_)
201  {
203  (
204  "Foam::slidingInterface::slidingInterface\n"
205  "(\n"
206  " const word& name,\n"
207  " const label index,\n"
208  " const polyTopoChanger& mme,\n"
209  " const word& masterFaceZoneName,\n"
210  " const word& slaveFaceZoneName,\n"
211  " const word& cutFaceZoneName,\n"
212  " const word& cutPointZoneName,\n"
213  " const word& masterPatchName,\n"
214  " const word& slavePatchName,\n"
215  " const typeOfMatch tom,\n"
216  " const bool coupleDecouple\n"
217  ")"
218  ) << "Creation of a sliding interface from components "
219  << "in attached state not supported."
220  << abort(FatalError);
221  }
222  else
223  {
224  calcAttachedAddressing();
225  }
226 }
227 
228 
229 // Construct from components
230 Foam::slidingInterface::slidingInterface
231 (
232  const word& name,
233  const dictionary& dict,
234  const label index,
235  const polyTopoChanger& mme
236 )
237 :
238  polyMeshModifier(name, index, mme, Switch(dict.lookup("active"))),
239  masterFaceZoneID_
240  (
241  dict.lookup("masterFaceZoneName"),
242  mme.mesh().faceZones()
243  ),
244  slaveFaceZoneID_
245  (
246  dict.lookup("slaveFaceZoneName"),
247  mme.mesh().faceZones()
248  ),
249  cutPointZoneID_
250  (
251  dict.lookup("cutPointZoneName"),
252  mme.mesh().pointZones()
253  ),
254  cutFaceZoneID_
255  (
256  dict.lookup("cutFaceZoneName"),
257  mme.mesh().faceZones()
258  ),
259  masterPatchID_
260  (
261  dict.lookup("masterPatchName"),
262  mme.mesh().boundaryMesh()
263  ),
264  slavePatchID_
265  (
266  dict.lookup("slavePatchName"),
267  mme.mesh().boundaryMesh()
268  ),
269  matchType_(typeOfMatchNames_.read((dict.lookup("typeOfMatch")))),
270  coupleDecouple_(dict.lookup("coupleDecouple")),
271  attached_(dict.lookup("attached")),
272  projectionAlgo_
273  (
274  intersection::algorithmNames_.read(dict.lookup("projection"))
275  ),
276  trigger_(false),
277  cutFaceMasterPtr_(NULL),
278  cutFaceSlavePtr_(NULL),
279  masterFaceCellsPtr_(NULL),
280  slaveFaceCellsPtr_(NULL),
281  masterStickOutFacesPtr_(NULL),
282  slaveStickOutFacesPtr_(NULL),
283  retiredPointMapPtr_(NULL),
284  cutPointEdgePairMapPtr_(NULL),
285  slavePointPointHitsPtr_(NULL),
286  slavePointEdgeHitsPtr_(NULL),
287  slavePointFaceHitsPtr_(NULL),
288  masterPointEdgeHitsPtr_(NULL),
289  projectedSlavePointsPtr_(NULL)
290 {
291  // Optionally default tolerances from dictionary
292  setTolerances(dict);
293 
294  checkDefinition();
295 
296  // If the interface is attached, the master and slave face zone addressing
297  // needs to be read in; otherwise it will be created
298  if (attached_)
299  {
300  if (debug)
301  {
302  Pout<< "slidingInterface::slidingInterface(...) "
303  << " for object " << name << " : "
304  << "Interface attached. Reading master and slave face zones "
305  << "and retired point lookup." << endl;
306  }
307 
308  // The face zone addressing is written out in the definition dictionary
309  masterFaceCellsPtr_ = new labelList(dict.lookup("masterFaceCells"));
310  slaveFaceCellsPtr_ = new labelList(dict.lookup("slaveFaceCells"));
311 
312  masterStickOutFacesPtr_ =
313  new labelList(dict.lookup("masterStickOutFaces"));
314  slaveStickOutFacesPtr_ =
315  new labelList(dict.lookup("slaveStickOutFaces"));
316 
317  retiredPointMapPtr_ = new Map<label>(dict.lookup("retiredPointMap"));
318  cutPointEdgePairMapPtr_ =
319  new Map<Pair<edge> >(dict.lookup("cutPointEdgePairMap"));
320  }
321  else
322  {
323  calcAttachedAddressing();
324  }
325 }
326 
327 
328 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
329 
331 {
332  clearOut();
333 }
334 
335 
336 void Foam::slidingInterface::clearAddressing() const
337 {
338  deleteDemandDrivenData(cutFaceMasterPtr_);
339  deleteDemandDrivenData(cutFaceSlavePtr_);
340 }
341 
342 
343 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
344 
346 {
347  return masterFaceZoneID_;
348 }
349 
350 
352 {
353  return slaveFaceZoneID_;
354 }
355 
356 
358 {
359  if (coupleDecouple_)
360  {
361  // Always changes. If not attached, project points
362  if (debug)
363  {
364  Pout<< "bool slidingInterface::changeTopology() const "
365  << "for object " << name() << " : "
366  << "Couple-decouple mode." << endl;
367  }
368 
369  if (!attached_)
370  {
371  projectPoints();
372  }
373  else
374  {
375  }
376 
377  return true;
378  }
379 
380  if
381  (
382  attached_
383  && !topoChanger().mesh().changing()
384  )
385  {
386  // If the mesh is not moving or morphing and the interface is
387  // already attached, the topology will not change
388  return false;
389  }
390  else
391  {
392  // Check if the motion changes point projection
393  return projectPoints();
394  }
395 }
396 
397 
399 {
400  if (coupleDecouple_)
401  {
402  if (attached_)
403  {
404  // Attached, coupling
405  decoupleInterface(ref);
406  }
407  else
408  {
409  // Detached, coupling
410  coupleInterface(ref);
411  }
412 
413  return;
414  }
415 
416  if (trigger_)
417  {
418  if (attached_)
419  {
420  // Clear old coupling data
421  clearCouple(ref);
422  }
423 
424  coupleInterface(ref);
425 
426  trigger_ = false;
427  }
428 }
429 
430 
432 {
433  if (debug)
434  {
435  Pout<< "void slidingInterface::modifyMotionPoints("
436  << "pointField& motionPoints) const for object " << name() << " : "
437  << "Adjusting motion points." << endl;
438  }
439 
440  const polyMesh& mesh = topoChanger().mesh();
441 
442  // Get point from the cut zone
443  const labelList& cutPoints = mesh.pointZones()[cutPointZoneID_.index()];
444 
445  if (cutPoints.size() && !projectedSlavePointsPtr_)
446  {
447  return;
448  }
449  else
450  {
451  const pointField& projectedSlavePoints = *projectedSlavePointsPtr_;
452 
453  const Map<label>& rpm = retiredPointMap();
454 
455  const Map<Pair<edge> >& cpepm = cutPointEdgePairMap();
456 
457  const Map<label>& slaveZonePointMap =
458  mesh.faceZones()[slaveFaceZoneID_.index()]().meshPointMap();
459 
460  const primitiveFacePatch& masterPatch =
461  mesh.faceZones()[masterFaceZoneID_.index()]();
462  const edgeList& masterEdges = masterPatch.edges();
463  const pointField& masterLocalPoints = masterPatch.localPoints();
464 
465  const primitiveFacePatch& slavePatch =
466  mesh.faceZones()[slaveFaceZoneID_.index()]();
467  const edgeList& slaveEdges = slavePatch.edges();
468  const pointField& slaveLocalPoints = slavePatch.localPoints();
469  const vectorField& slavePointNormals = slavePatch.pointNormals();
470 
471  forAll(cutPoints, pointI)
472  {
473  // Try to find the cut point in retired points
474  Map<label>::const_iterator rpmIter = rpm.find(cutPoints[pointI]);
475 
476  if (rpmIter != rpm.end())
477  {
478  if (debug)
479  {
480  Pout<< "p";
481  }
482 
483  // Cut point is a retired point
484  motionPoints[cutPoints[pointI]] =
485  projectedSlavePoints[slaveZonePointMap.find(rpmIter())()];
486  }
487  else
488  {
489  // A cut point is not a projected slave point. Therefore, it
490  // must be an edge-to-edge intersection.
491 
492  Map<Pair<edge> >::const_iterator cpepmIter =
493  cpepm.find(cutPoints[pointI]);
494 
495  if (cpepmIter != cpepm.end())
496  {
497  // Pout<< "Need to re-create hit for point "
498  // << cutPoints[pointI]
499  // << " lookup: " << cpepmIter()
500  // << endl;
501 
502  // Note.
503  // The edge cutting code is repeated in
504  // slidingInterface::coupleInterface. This is done for
505  // efficiency reasons and avoids multiple creation of
506  // cutting planes. Please update both simultaneously.
507  //
508  const edge& globalMasterEdge = cpepmIter().first();
509 
510  const label curMasterEdgeIndex =
511  masterPatch.whichEdge
512  (
513  edge
514  (
515  masterPatch.whichPoint
516  (
517  globalMasterEdge.start()
518  ),
519  masterPatch.whichPoint
520  (
521  globalMasterEdge.end()
522  )
523  )
524  );
525 
526  const edge& cme = masterEdges[curMasterEdgeIndex];
527 
528  // Pout<< "curMasterEdgeIndex: " << curMasterEdgeIndex
529  // << " cme: " << cme
530  // << endl;
531 
532  const edge& globalSlaveEdge = cpepmIter().second();
533 
534  const label curSlaveEdgeIndex =
535  slavePatch.whichEdge
536  (
537  edge
538  (
539  slavePatch.whichPoint
540  (
541  globalSlaveEdge.start()
542  ),
543  slavePatch.whichPoint
544  (
545  globalSlaveEdge.end()
546  )
547  )
548  );
549 
550  const edge& curSlaveEdge = slaveEdges[curSlaveEdgeIndex];
551  // Pout<< "curSlaveEdgeIndex: " << curSlaveEdgeIndex
552  // << " curSlaveEdge: " << curSlaveEdge
553  // << endl;
554  const point& a = projectedSlavePoints[curSlaveEdge.start()];
555  const point& b = projectedSlavePoints[curSlaveEdge.end()];
556 
557  point c =
558  0.5*
559  (
560  slaveLocalPoints[curSlaveEdge.start()]
561  + slavePointNormals[curSlaveEdge.start()]
562  + slaveLocalPoints[curSlaveEdge.end()]
563  + slavePointNormals[curSlaveEdge.end()]
564  );
565 
566  // Create the plane
567  plane cutPlane(a, b, c);
568 
569  linePointRef curSlaveLine =
570  curSlaveEdge.line(slaveLocalPoints);
571  const scalar curSlaveLineMag = curSlaveLine.mag();
572 
573  scalar cutOnMaster =
574  cutPlane.lineIntersect
575  (
576  cme.line(masterLocalPoints)
577  );
578 
579  if
580  (
581  cutOnMaster > edgeEndCutoffTol_
582  && cutOnMaster < 1.0 - edgeEndCutoffTol_
583  )
584  {
585  // Master is cut, check the slave
586  point masterCutPoint =
587  masterLocalPoints[cme.start()]
588  + cutOnMaster*cme.vec(masterLocalPoints);
589 
590  pointHit slaveCut =
591  curSlaveLine.nearestDist(masterCutPoint);
592 
593  if (slaveCut.hit())
594  {
595  // Strict checking of slave cut to avoid capturing
596  // end points.
597  scalar cutOnSlave =
598  (
599  (
600  slaveCut.hitPoint()
601  - curSlaveLine.start()
602  ) & curSlaveLine.vec()
603  )/sqr(curSlaveLineMag);
604 
605  // Calculate merge tolerance from the
606  // target edge length
607  scalar mergeTol =
608  edgeCoPlanarTol_*mag(b - a);
609 
610  if
611  (
612  cutOnSlave > edgeEndCutoffTol_
613  && cutOnSlave < 1.0 - edgeEndCutoffTol_
614  && slaveCut.distance() < mergeTol
615  )
616  {
617  // Cut both master and slave.
618  motionPoints[cutPoints[pointI]] =
619  masterCutPoint;
620  }
621  }
622  else
623  {
624  Pout<< "Missed slave edge!!! This is an error. "
625  << "Master edge: "
626  << cme.line(masterLocalPoints)
627  << " slave edge: " << curSlaveLine
628  << " point: " << masterCutPoint
629  << " weight: " <<
630  (
631  (
632  slaveCut.missPoint()
633  - curSlaveLine.start()
634  ) & curSlaveLine.vec()
635  )/sqr(curSlaveLineMag)
636  << endl;
637  }
638  }
639  else
640  {
641  Pout<< "Missed master edge!!! This is an error"
642  << endl;
643  }
644  }
645  else
646  {
648  (
649  "void slidingInterface::modifyMotionPoints"
650  "(pointField&) const"
651  ) << "Cut point " << cutPoints[pointI]
652  << " not recognised as either the projected "
653  << "or as intersection point. Error in point "
654  << "projection or data mapping"
655  << abort(FatalError);
656  }
657  }
658  }
659  if (debug)
660  {
661  Pout<< endl;
662  }
663  }
664 }
665 
666 
668 {
669  if (debug)
670  {
671  Pout<< "void slidingInterface::updateMesh(const mapPolyMesh& m)"
672  << " const for object " << name() << " : "
673  << "Updating topology." << endl;
674  }
675 
676  // Mesh has changed topologically. Update local topological data
677  const polyMesh& mesh = topoChanger().mesh();
678 
679  masterFaceZoneID_.update(mesh.faceZones());
680  slaveFaceZoneID_.update(mesh.faceZones());
681  cutPointZoneID_.update(mesh.pointZones());
682  cutFaceZoneID_.update(mesh.faceZones());
683 
684  masterPatchID_.update(mesh.boundaryMesh());
685  slavePatchID_.update(mesh.boundaryMesh());
686 
687 //MJ:Disabled updating
688 // if (!attached())
689 // {
690 // calcAttachedAddressing();
691 // }
692 // else
693 // {
694 // renumberAttachedAddressing(m);
695 // }
696 }
697 
698 
700 {
701  if (!projectedSlavePointsPtr_)
702  {
703  projectPoints();
704  }
705 
706  return *projectedSlavePointsPtr_;
707 }
708 
710 {
711  pointMergeTol_ = dict.lookupOrDefault<scalar>
712  (
713  "pointMergeTol",
714  pointMergeTol_
715  );
716  edgeMergeTol_ = dict.lookupOrDefault<scalar>
717  (
718  "edgeMergeTol",
719  edgeMergeTol_
720  );
721  nFacesPerSlaveEdge_ = dict.lookupOrDefault<label>
722  (
723  "nFacesPerSlaveEdge",
724  nFacesPerSlaveEdge_
725  );
726  edgeFaceEscapeLimit_ = dict.lookupOrDefault<label>
727  (
728  "edgeFaceEscapeLimit",
729  edgeFaceEscapeLimit_
730  );
731  integralAdjTol_ = dict.lookupOrDefault<scalar>
732  (
733  "integralAdjTol",
734  integralAdjTol_
735  );
736  edgeMasterCatchFraction_ = dict.lookupOrDefault<scalar>
737  (
738  "edgeMasterCatchFraction",
739  edgeMasterCatchFraction_
740  );
741  edgeCoPlanarTol_ = dict.lookupOrDefault<scalar>
742  (
743  "edgeCoPlanarTol",
744  edgeCoPlanarTol_
745  );
746  edgeEndCutoffTol_ = dict.lookupOrDefault<scalar>
747  (
748  "edgeEndCutoffTol",
749  edgeEndCutoffTol_
750  );
751 
752  if (report)
753  {
754  Info<< "Sliding interface parameters:" << nl
755  << "pointMergeTol : " << pointMergeTol_ << nl
756  << "edgeMergeTol : " << edgeMergeTol_ << nl
757  << "nFacesPerSlaveEdge : " << nFacesPerSlaveEdge_ << nl
758  << "edgeFaceEscapeLimit : " << edgeFaceEscapeLimit_ << nl
759  << "integralAdjTol : " << integralAdjTol_ << nl
760  << "edgeMasterCatchFraction : " << edgeMasterCatchFraction_ << nl
761  << "edgeCoPlanarTol : " << edgeCoPlanarTol_ << nl
762  << "edgeEndCutoffTol : " << edgeEndCutoffTol_ << endl;
763  }
764 }
765 
766 
768 {
769  os << nl << type() << nl
770  << name()<< nl
771  << masterFaceZoneID_.name() << nl
772  << slaveFaceZoneID_.name() << nl
773  << cutPointZoneID_.name() << nl
774  << cutFaceZoneID_.name() << nl
775  << masterPatchID_.name() << nl
776  << slavePatchID_.name() << nl
777  << typeOfMatchNames_[matchType_] << nl
778  << coupleDecouple_ << nl
779  << attached_ << endl;
780 }
781 
782 
783 // To write out all those tolerances
784 #define WRITE_NON_DEFAULT(name) \
785  if ( name ## _ != name ## Default_ )\
786  { \
787  os << " " #name " " << name ## _ << token::END_STATEMENT << nl; \
788  }
789 
790 
792 {
793  os << nl << name() << nl << token::BEGIN_BLOCK << nl
794  << " type " << type() << token::END_STATEMENT << nl
795  << " masterFaceZoneName " << masterFaceZoneID_.name()
797  << " slaveFaceZoneName " << slaveFaceZoneID_.name()
799  << " cutPointZoneName " << cutPointZoneID_.name()
801  << " cutFaceZoneName " << cutFaceZoneID_.name()
803  << " masterPatchName " << masterPatchID_.name()
805  << " slavePatchName " << slavePatchID_.name()
807  << " typeOfMatch " << typeOfMatchNames_[matchType_]
809  << " coupleDecouple " << coupleDecouple_
811  << " projection " << intersection::algorithmNames_[projectionAlgo_]
813  << " attached " << attached_
815  << " active " << active()
816  << token::END_STATEMENT << nl;
817 
818  if (attached_)
819  {
820  masterFaceCellsPtr_->writeEntry("masterFaceCells", os);
821  slaveFaceCellsPtr_->writeEntry("slaveFaceCells", os);
822  masterStickOutFacesPtr_->writeEntry("masterStickOutFaces", os);
823  slaveStickOutFacesPtr_->writeEntry("slaveStickOutFaces", os);
824 
825  os << " retiredPointMap " << retiredPointMap()
826  << token::END_STATEMENT << nl
827  << " cutPointEdgePairMap " << cutPointEdgePairMap()
828  << token::END_STATEMENT << nl;
829  }
830 
831  WRITE_NON_DEFAULT(pointMergeTol)
832  WRITE_NON_DEFAULT(edgeMergeTol)
833  WRITE_NON_DEFAULT(nFacesPerSlaveEdge)
834  WRITE_NON_DEFAULT(edgeFaceEscapeLimit)
835  WRITE_NON_DEFAULT(integralAdjTol)
836  WRITE_NON_DEFAULT(edgeMasterCatchFraction)
837  WRITE_NON_DEFAULT(edgeCoPlanarTol)
838  WRITE_NON_DEFAULT(edgeEndCutoffTol)
839 
840  os << token::END_BLOCK << endl;
841 }
842 
843 
844 // ************************************************************************* //
PointHit< Point > nearestDist(const Point &p) const
Return nearest distance to line from a given point.
Definition: lineI.H:95
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
PointRef start() const
Return first vertex.
Definition: lineI.H:60
virtual bool update()=0
Update the mesh for both mesh motion and topology change.
const Field< PointType > & pointNormals() const
Return point normals for patch.
bool empty() const
Return true if the PtrList is empty (ie, size() is zero).
Definition: PtrListI.H:39
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
dimensioned< scalar > mag(const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void write(Ostream &) const
Write.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
const Point & missPoint() const
Return miss point.
Definition: PointHit.H:145
virtual void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Virtual base class for mesh modifiers.
const Point & hitPoint() const
Return hit point.
Definition: PointHit.H:126
label whichPoint(const label gp) const
Given a global point index, return the local point index.
typeOfMatch
Type of match.
void deleteDemandDrivenData(DataPtr &dataPtr)
A class for handling words, derived from string.
Definition: word.H:59
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
const faceZoneID & slaveFaceZoneID() const
Return slave face zone ID.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
const pointZoneMesh & pointZones() const
Return point zone mesh.
Definition: polyMesh.H:457
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
messageStream Info
const polyMesh & mesh() const
Return the mesh reference.
static const NamedEnum< algorithm, 3 > algorithmNames_
Projection algorithm names.
Definition: intersection.H:82
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void setRefinement(polyTopoChange &) const
Insert the layer addition/removal instructions.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Namespace for OpenFOAM.
A line primitive.
Definition: line.H:56
static const NamedEnum< typeOfMatch, 2 > typeOfMatchNames_
Direction names.
virtual void modifyMotionPoints(pointField &motionPoints) const
Modify motion points to comply with the topological change.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
linePointRef line(const pointField &) const
Return edge line.
Definition: edgeI.H:169
void setTolerances(const dictionary &, bool report=false)
Set the tolerances from the values in a dictionary.
#define WRITE_NON_DEFAULT(name)
#define forAll(list, i)
Definition: UList.H:421
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const faceZoneID & masterFaceZoneID() const
Return master face zone ID.
Macros for easy insertion into run-time selection tables.
List of mesh modifiers defining the mesh dynamics.
virtual ~slidingInterface()
Destructor.
vector vec(const pointField &) const
Return the vector (end - start)
Definition: edgeI.H:157
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const pointField & pointProjection() const
Return projected points for a slave patch.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
A list of faces which address into the list of points.
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
label whichEdge(const edge &) const
Given an edge in local point labels, return its.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
error FatalError
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
label end() const
Return end vertex label.
Definition: edgeI.H:92
List< label > labelList
A List of labels.
Definition: labelList.H:56
bool hit() const
Is there a hit.
Definition: PointHit.H:120
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
scalar mag() const
Return scalar magnitude.
Definition: lineI.H:80
const dimensionedScalar c
Speed of light in a vacuum.
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Direct mesh changes based on v1.3 polyTopoChange syntax.
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:106
scalar lineIntersect(const line< Point, PointRef > &l) const
Return the cutting point between the plane and.
Definition: plane.H:180
label start() const
Return start vertex label.
Definition: edgeI.H:81
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:139
virtual bool changeTopology() const
Check for topology change.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual void writeDict(Ostream &) const
Write dictionary.
Point vec() const
Return start-end vector.
Definition: lineI.H:87
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53