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