FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

slidingInterface.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "slidingInterface.H"
00027 #include <dynamicMesh/polyTopoChanger.H>
00028 #include <OpenFOAM/polyMesh.H>
00029 #include <dynamicMesh/polyTopoChange.H>
00030 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00031 #include <OpenFOAM/plane.H>
00032 
00033 // Index of debug signs:
00034 // p - adjusting a projection point
00035 // * - adjusting edge intersection
00036 
00037 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00038 
00039 namespace Foam
00040 {
00041     defineTypeNameAndDebug(slidingInterface, 0);
00042     addToRunTimeSelectionTable
00043     (
00044         polyMeshModifier,
00045         slidingInterface,
00046         dictionary
00047     );
00048 }
00049 
00050 
00051 template<>
00052 const char* Foam::NamedEnum<Foam::slidingInterface::typeOfMatch, 2>::names[] =
00053 {
00054     "integral",
00055     "partial"
00056 };
00057 
00058 
00059 const Foam::NamedEnum<Foam::slidingInterface::typeOfMatch, 2>
00060 Foam::slidingInterface::typeOfMatchNames_;
00061 
00062 
00063 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00064 
00065 void Foam::slidingInterface::checkDefinition()
00066 {
00067     const polyMesh& mesh = topoChanger().mesh();
00068 
00069     if
00070     (
00071         !masterFaceZoneID_.active()
00072      || !slaveFaceZoneID_.active()
00073      || !cutPointZoneID_.active()
00074      || !cutFaceZoneID_.active()
00075      || !masterPatchID_.active()
00076      || !slavePatchID_.active()
00077     )
00078     {
00079         FatalErrorIn
00080         (
00081             "void slidingInterface::checkDefinition()"
00082         )   << "Not all zones and patches needed in the definition "
00083             << "have been found.  Please check your mesh definition."
00084             << abort(FatalError);
00085     }
00086 
00087     // Check the sizes and set up state
00088     if
00089     (
00090         mesh.faceZones()[masterFaceZoneID_.index()].empty()
00091      || mesh.faceZones()[slaveFaceZoneID_.index()].empty()
00092     )
00093     {
00094         FatalErrorIn("void slidingInterface::checkDefinition()")
00095             << "Master or slave face zone contain no faces.  "
00096             << "Please check your mesh definition."
00097             << abort(FatalError);
00098     }
00099 
00100     if (debug)
00101     {
00102         Pout<< "Sliding interface object " << name() << " :" << nl
00103             << "    master face zone: " << masterFaceZoneID_.index() << nl
00104             << "    slave face zone: " << slaveFaceZoneID_.index() << endl;
00105     }
00106 }
00107 
00108 
00109 void Foam::slidingInterface::clearOut() const
00110 {
00111     clearPointProjection();
00112     clearAttachedAddressing();
00113     clearAddressing();
00114 }
00115 
00116 
00117 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00118 
00119 
00120 // Construct from components
00121 Foam::slidingInterface::slidingInterface
00122 (
00123     const word& name,
00124     const label index,
00125     const polyTopoChanger& mme,
00126     const word& masterFaceZoneName,
00127     const word& slaveFaceZoneName,
00128     const word& cutPointZoneName,
00129     const word& cutFaceZoneName,
00130     const word& masterPatchName,
00131     const word& slavePatchName,
00132     const typeOfMatch tom,
00133     const bool coupleDecouple,
00134     const intersection::algorithm algo
00135 )
00136 :
00137     polyMeshModifier(name, index, mme, true),
00138     masterFaceZoneID_
00139     (
00140         masterFaceZoneName,
00141         mme.mesh().faceZones()
00142     ),
00143     slaveFaceZoneID_
00144     (
00145         slaveFaceZoneName,
00146         mme.mesh().faceZones()
00147     ),
00148     cutPointZoneID_
00149     (
00150         cutPointZoneName,
00151         mme.mesh().pointZones()
00152     ),
00153     cutFaceZoneID_
00154     (
00155         cutFaceZoneName,
00156         mme.mesh().faceZones()
00157     ),
00158     masterPatchID_
00159     (
00160         masterPatchName,
00161         mme.mesh().boundaryMesh()
00162     ),
00163     slavePatchID_
00164     (
00165         slavePatchName,
00166         mme.mesh().boundaryMesh()
00167     ),
00168     matchType_(tom),
00169     coupleDecouple_(coupleDecouple),
00170     attached_(false),
00171     projectionAlgo_(algo),
00172     trigger_(false),
00173     pointMergeTol_(pointMergeTolDefault_),
00174     edgeMergeTol_(edgeMergeTolDefault_),
00175     nFacesPerSlaveEdge_(nFacesPerSlaveEdgeDefault_),
00176     edgeFaceEscapeLimit_(edgeFaceEscapeLimitDefault_),
00177     integralAdjTol_(integralAdjTolDefault_),
00178     edgeMasterCatchFraction_(edgeMasterCatchFractionDefault_),
00179     edgeCoPlanarTol_(edgeCoPlanarTolDefault_),
00180     edgeEndCutoffTol_(edgeEndCutoffTolDefault_),
00181     cutFaceMasterPtr_(NULL),
00182     cutFaceSlavePtr_(NULL),
00183     masterFaceCellsPtr_(NULL),
00184     slaveFaceCellsPtr_(NULL),
00185     masterStickOutFacesPtr_(NULL),
00186     slaveStickOutFacesPtr_(NULL),
00187     retiredPointMapPtr_(NULL),
00188     cutPointEdgePairMapPtr_(NULL),
00189     slavePointPointHitsPtr_(NULL),
00190     slavePointEdgeHitsPtr_(NULL),
00191     slavePointFaceHitsPtr_(NULL),
00192     masterPointEdgeHitsPtr_(NULL),
00193     projectedSlavePointsPtr_(NULL)
00194 {
00195     checkDefinition();
00196 
00197     if (attached_)
00198     {
00199         FatalErrorIn
00200         (
00201             "Foam::slidingInterface::slidingInterface\n"
00202             "(\n"
00203             "    const word& name,\n"
00204             "    const label index,\n"
00205             "    const polyTopoChanger& mme,\n"
00206             "    const word& masterFaceZoneName,\n"
00207             "    const word& slaveFaceZoneName,\n"
00208             "    const word& cutFaceZoneName,\n"
00209             "    const word& cutPointZoneName,\n"
00210             "    const word& masterPatchName,\n"
00211             "    const word& slavePatchName,\n"
00212             "    const typeOfMatch tom,\n"
00213             "    const bool coupleDecouple\n"
00214             ")"
00215         )   << "Creation of a sliding interface from components "
00216             << "in attached state not supported."
00217             << abort(FatalError);
00218     }
00219     else
00220     {
00221         calcAttachedAddressing();
00222     }
00223 }
00224 
00225 
00226 // Construct from components
00227 Foam::slidingInterface::slidingInterface
00228 (
00229     const word& name,
00230     const dictionary& dict,
00231     const label index,
00232     const polyTopoChanger& mme
00233 )
00234 :
00235     polyMeshModifier(name, index, mme, Switch(dict.lookup("active"))),
00236     masterFaceZoneID_
00237     (
00238         dict.lookup("masterFaceZoneName"),
00239         mme.mesh().faceZones()
00240     ),
00241     slaveFaceZoneID_
00242     (
00243         dict.lookup("slaveFaceZoneName"),
00244         mme.mesh().faceZones()
00245     ),
00246     cutPointZoneID_
00247     (
00248         dict.lookup("cutPointZoneName"),
00249         mme.mesh().pointZones()
00250     ),
00251     cutFaceZoneID_
00252     (
00253         dict.lookup("cutFaceZoneName"),
00254         mme.mesh().faceZones()
00255     ),
00256     masterPatchID_
00257     (
00258         dict.lookup("masterPatchName"),
00259         mme.mesh().boundaryMesh()
00260     ),
00261     slavePatchID_
00262     (
00263         dict.lookup("slavePatchName"),
00264         mme.mesh().boundaryMesh()
00265     ),
00266     matchType_(typeOfMatchNames_.read((dict.lookup("typeOfMatch")))),
00267     coupleDecouple_(dict.lookup("coupleDecouple")),
00268     attached_(dict.lookup("attached")),
00269     projectionAlgo_
00270     (
00271         intersection::algorithmNames_.read(dict.lookup("projection"))
00272     ),
00273     trigger_(false),
00274     cutFaceMasterPtr_(NULL),
00275     cutFaceSlavePtr_(NULL),
00276     masterFaceCellsPtr_(NULL),
00277     slaveFaceCellsPtr_(NULL),
00278     masterStickOutFacesPtr_(NULL),
00279     slaveStickOutFacesPtr_(NULL),
00280     retiredPointMapPtr_(NULL),
00281     cutPointEdgePairMapPtr_(NULL),
00282     slavePointPointHitsPtr_(NULL),
00283     slavePointEdgeHitsPtr_(NULL),
00284     slavePointFaceHitsPtr_(NULL),
00285     masterPointEdgeHitsPtr_(NULL),
00286     projectedSlavePointsPtr_(NULL)
00287 {
00288     // Optionally default tolerances from dictionary
00289     setTolerances(dict);
00290 
00291     checkDefinition();
00292 
00293     // If the interface is attached, the master and slave face zone addressing
00294     // needs to be read in; otherwise it will be created
00295     if (attached_)
00296     {
00297         if (debug)
00298         {
00299             Pout<< "slidingInterface::slidingInterface(...) "
00300                 << " for object " << name << " : "
00301                 << "Interface attached.  Reading master and slave face zones "
00302                 << "and retired point lookup." << endl;
00303         }
00304 
00305         // The face zone addressing is written out in the definition dictionary
00306         masterFaceCellsPtr_ = new labelList(dict.lookup("masterFaceCells"));
00307         slaveFaceCellsPtr_ = new labelList(dict.lookup("slaveFaceCells"));
00308 
00309         masterStickOutFacesPtr_ =
00310             new labelList(dict.lookup("masterStickOutFaces"));
00311         slaveStickOutFacesPtr_ =
00312             new labelList(dict.lookup("slaveStickOutFaces"));
00313 
00314         retiredPointMapPtr_ = new Map<label>(dict.lookup("retiredPointMap"));
00315         cutPointEdgePairMapPtr_ =
00316             new Map<Pair<edge> >(dict.lookup("cutPointEdgePairMap"));
00317     }
00318     else
00319     {
00320         calcAttachedAddressing();
00321     }
00322 }
00323 
00324 
00325 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00326 
00327 Foam::slidingInterface::~slidingInterface()
00328 {
00329     clearOut();
00330 }
00331 
00332 
00333 void Foam::slidingInterface::clearAddressing() const
00334 {
00335     deleteDemandDrivenData(cutFaceMasterPtr_);
00336     deleteDemandDrivenData(cutFaceSlavePtr_);
00337 }
00338 
00339 
00340 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00341 
00342 const Foam::faceZoneID& Foam::slidingInterface::masterFaceZoneID() const
00343 {
00344     return masterFaceZoneID_;
00345 }
00346 
00347 
00348 const Foam::faceZoneID& Foam::slidingInterface::slaveFaceZoneID() const
00349 {
00350     return slaveFaceZoneID_;
00351 }
00352 
00353 
00354 bool Foam::slidingInterface::changeTopology() const
00355 {
00356     if (coupleDecouple_)
00357     {
00358         // Always changes.  If not attached, project points
00359         if (debug)
00360         {
00361             Pout<< "bool slidingInterface::changeTopology() const "
00362                 << "for object " << name() << " : "
00363                 << "Couple-decouple mode." << endl;
00364         }
00365 
00366         if (!attached_)
00367         {
00368             projectPoints();
00369         }
00370         else
00371         {
00372         }
00373 
00374         return true;
00375     }
00376 
00377     if
00378     (
00379         attached_
00380      && !topoChanger().mesh().changing()
00381     )
00382     {
00383         // If the mesh is not moving or morphing and the interface is
00384         // already attached, the topology will not change
00385         return false;
00386     }
00387     else
00388     {
00389         // Check if the motion changes point projection
00390         return projectPoints();
00391     }
00392 }
00393 
00394 
00395 void Foam::slidingInterface::setRefinement(polyTopoChange& ref) const
00396 {
00397     if (coupleDecouple_)
00398     {
00399         if (attached_)
00400         {
00401             // Attached, coupling
00402             decoupleInterface(ref);
00403         }
00404         else
00405         {
00406             // Detached, coupling
00407             coupleInterface(ref);
00408         }
00409 
00410         return;
00411     }
00412 
00413     if (trigger_)
00414     {
00415         if (attached_)
00416         {
00417             // Clear old coupling data
00418             clearCouple(ref);
00419         }
00420 
00421         coupleInterface(ref);
00422         
00423         trigger_ = false;
00424     }
00425 }
00426 
00427 
00428 void Foam::slidingInterface::modifyMotionPoints(pointField& motionPoints) const
00429 {
00430     if (debug)
00431     {
00432         Pout<< "void slidingInterface::modifyMotionPoints(" 
00433             << "pointField& motionPoints) const for object " << name() << " : "
00434             << "Adjusting motion points." << endl;
00435     }
00436 
00437     const polyMesh& mesh = topoChanger().mesh();
00438 
00439     // Get point from the cut zone
00440     const labelList& cutPoints = mesh.pointZones()[cutPointZoneID_.index()];
00441 
00442     if (cutPoints.size() && !projectedSlavePointsPtr_)
00443     {
00444         return;
00445     }
00446     else
00447     {
00448         const pointField& projectedSlavePoints = *projectedSlavePointsPtr_;
00449 
00450         const Map<label>& rpm = retiredPointMap();
00451 
00452         const Map<Pair<edge> >& cpepm = cutPointEdgePairMap();
00453 
00454         const Map<label>& slaveZonePointMap =
00455             mesh.faceZones()[slaveFaceZoneID_.index()]().meshPointMap();
00456 
00457         const primitiveFacePatch& masterPatch =
00458             mesh.faceZones()[masterFaceZoneID_.index()]();
00459         const edgeList& masterEdges = masterPatch.edges();
00460         const pointField& masterLocalPoints = masterPatch.localPoints();
00461 
00462         const primitiveFacePatch& slavePatch =
00463             mesh.faceZones()[slaveFaceZoneID_.index()]();
00464         const edgeList& slaveEdges = slavePatch.edges();
00465         const pointField& slaveLocalPoints = slavePatch.localPoints();
00466         const vectorField& slavePointNormals = slavePatch.pointNormals();
00467 
00468         forAll (cutPoints, pointI)
00469         {
00470             // Try to find the cut point in retired points
00471             Map<label>::const_iterator rpmIter = rpm.find(cutPoints[pointI]);
00472 
00473             if (rpmIter != rpm.end())
00474             {
00475                 if (debug)
00476                 {
00477                     Pout << "p";
00478                 }
00479 
00480                 // Cut point is a retired point
00481                 motionPoints[cutPoints[pointI]] =
00482                     projectedSlavePoints[slaveZonePointMap.find(rpmIter())()];
00483             }
00484             else
00485             {
00486                 // A cut point is not a projected slave point.  Therefore, it
00487                 // must be an edge-to-edge intersection.  
00488 
00489                 Map<Pair<edge> >::const_iterator cpepmIter =
00490                     cpepm.find(cutPoints[pointI]);
00491 
00492                 if (cpepmIter != cpepm.end())
00493                 {
00494 //                     Pout << "Need to re-create hit for point " << cutPoints[pointI] << " lookup: " << cpepmIter() << endl;
00495 
00496                     // Note.
00497                     // The edge cutting code is repeated in
00498                     // slidingInterface::coupleInterface.  This is done for
00499                     // efficiency reasons and avoids multiple creation of
00500                     // cutting planes.  Please update both simultaneously.
00501                     // 
00502                     const edge& globalMasterEdge = cpepmIter().first();
00503 
00504                     const label curMasterEdgeIndex =
00505                         masterPatch.whichEdge
00506                         (
00507                             edge
00508                             (
00509                                 masterPatch.whichPoint
00510                                 (
00511                                     globalMasterEdge.start()
00512                                 ),
00513                                 masterPatch.whichPoint
00514                                 (
00515                                     globalMasterEdge.end()
00516                                 )
00517                             )
00518                         );
00519 
00520                     const edge& cme = masterEdges[curMasterEdgeIndex];
00521 //                     Pout << "curMasterEdgeIndex: " << curMasterEdgeIndex << " cme: " << cme << endl;
00522                     const edge& globalSlaveEdge = cpepmIter().second();
00523 
00524                     const label curSlaveEdgeIndex =
00525                         slavePatch.whichEdge
00526                         (
00527                             edge
00528                             (
00529                                 slavePatch.whichPoint
00530                                 (
00531                                     globalSlaveEdge.start()
00532                                 ),
00533                                 slavePatch.whichPoint
00534                                 (
00535                                     globalSlaveEdge.end()
00536                                 )
00537                             )
00538                         );
00539 
00540                     const edge& curSlaveEdge = slaveEdges[curSlaveEdgeIndex];
00541 //                     Pout << "curSlaveEdgeIndex: " << curSlaveEdgeIndex << " curSlaveEdge: " << curSlaveEdge << endl;
00542                     const point& a = projectedSlavePoints[curSlaveEdge.start()];
00543                     const point& b = projectedSlavePoints[curSlaveEdge.end()];
00544 
00545                     point c =
00546                         0.5*
00547                         (
00548                             slaveLocalPoints[curSlaveEdge.start()]
00549                           + slavePointNormals[curSlaveEdge.start()]
00550                           + slaveLocalPoints[curSlaveEdge.end()]
00551                           + slavePointNormals[curSlaveEdge.end()]
00552                         );
00553 
00554                     // Create the plane
00555                     plane cutPlane(a, b, c);
00556 
00557                     linePointRef curSlaveLine =
00558                         curSlaveEdge.line(slaveLocalPoints);
00559                     const scalar curSlaveLineMag = curSlaveLine.mag();
00560 
00561                     scalar cutOnMaster =
00562                         cutPlane.lineIntersect
00563                         (
00564                             cme.line(masterLocalPoints)
00565                         );
00566 
00567                     if
00568                     (
00569                         cutOnMaster > edgeEndCutoffTol_
00570                      && cutOnMaster < 1.0 - edgeEndCutoffTol_
00571                     )
00572                     {
00573                         // Master is cut, check the slave
00574                         point masterCutPoint =
00575                             masterLocalPoints[cme.start()]
00576                           + cutOnMaster*cme.vec(masterLocalPoints);
00577 
00578                         pointHit slaveCut =
00579                             curSlaveLine.nearestDist(masterCutPoint);
00580 
00581                         if (slaveCut.hit())
00582                         {
00583                             // Strict checking of slave cut to avoid capturing
00584                             // end points.  
00585                             scalar cutOnSlave =
00586                                 (
00587                                     (
00588                                         slaveCut.hitPoint()
00589                                       - curSlaveLine.start()
00590                                     ) & curSlaveLine.vec()
00591                                 )/sqr(curSlaveLineMag);
00592 
00593                             // Calculate merge tolerance from the
00594                             // target edge length
00595                             scalar mergeTol =
00596                                 edgeCoPlanarTol_*mag(b - a);
00597 
00598                             if
00599                             (
00600                                 cutOnSlave > edgeEndCutoffTol_
00601                              && cutOnSlave < 1.0 - edgeEndCutoffTol_
00602                              && slaveCut.distance() < mergeTol
00603                             )
00604                             {
00605                                 // Cut both master and slave.
00606                                 motionPoints[cutPoints[pointI]] =
00607                                     masterCutPoint;
00608                             }
00609                         }
00610                         else
00611                         {
00612                             Pout<< "Missed slave edge!!!  This is an error.  "
00613                                 << "Master edge: "
00614                                 << cme.line(masterLocalPoints)
00615                                 << " slave edge: " << curSlaveLine
00616                                 << " point: " << masterCutPoint
00617                                 << " weight: " << 
00618                                 (
00619                                     (
00620                                         slaveCut.missPoint()
00621                                       - curSlaveLine.start()
00622                                     ) & curSlaveLine.vec()
00623                                 )/sqr(curSlaveLineMag)
00624                                 << endl;
00625                         }
00626                     }
00627                     else
00628                     {
00629                         Pout<< "Missed master edge!!!  This is an error"
00630                             << endl;
00631                     }
00632                 }
00633                 else
00634                 {
00635                     FatalErrorIn
00636                     (
00637                         "void slidingInterface::modifyMotionPoints"
00638                         "(pointField&) const"
00639                     )   << "Cut point " << cutPoints[pointI]
00640                         << " not recognised as either the projected "
00641                         << "or as intersection point.  Error in point "
00642                         << "projection or data mapping"
00643                         << abort(FatalError);
00644                 }
00645             }
00646         }
00647         if (debug)
00648         {
00649             Pout << endl;
00650         }
00651     }
00652 }
00653 
00654 
00655 void Foam::slidingInterface::updateMesh(const mapPolyMesh& m)
00656 {
00657     if (debug)
00658     {
00659         Pout<< "void slidingInterface::updateMesh(const mapPolyMesh& m)" 
00660             << " const for object " << name() << " : "
00661             << "Updating topology." << endl;
00662     }
00663 
00664     // Mesh has changed topologically.  Update local topological data
00665     const polyMesh& mesh = topoChanger().mesh();
00666 
00667     masterFaceZoneID_.update(mesh.faceZones());
00668     slaveFaceZoneID_.update(mesh.faceZones());
00669     cutPointZoneID_.update(mesh.pointZones());
00670     cutFaceZoneID_.update(mesh.faceZones());
00671 
00672     masterPatchID_.update(mesh.boundaryMesh());
00673     slavePatchID_.update(mesh.boundaryMesh());
00674 
00675 //MJ:Disabled updating
00676 //    if (!attached())
00677 //    {
00678 //        calcAttachedAddressing();
00679 //    }
00680 //    else
00681 //    {
00682 //        renumberAttachedAddressing(m);
00683 //    }
00684 }
00685 
00686 
00687 const Foam::pointField& Foam::slidingInterface::pointProjection() const
00688 {
00689     if (!projectedSlavePointsPtr_)
00690     {
00691         projectPoints();
00692     }
00693 
00694     return *projectedSlavePointsPtr_;
00695 }
00696 
00697 void Foam::slidingInterface::setTolerances(const dictionary&dict, bool report)
00698 {
00699     pointMergeTol_ = dict.lookupOrDefault<scalar>
00700     (
00701         "pointMergeTol",
00702         pointMergeTol_
00703     );
00704     edgeMergeTol_ = dict.lookupOrDefault<scalar>
00705     (
00706         "edgeMergeTol",
00707         edgeMergeTol_
00708     );
00709     nFacesPerSlaveEdge_ = dict.lookupOrDefault<label>
00710     (
00711         "nFacesPerSlaveEdge",
00712         nFacesPerSlaveEdge_
00713     );
00714     edgeFaceEscapeLimit_ = dict.lookupOrDefault<label>
00715     (
00716         "edgeFaceEscapeLimit",
00717         edgeFaceEscapeLimit_
00718     );
00719     integralAdjTol_ = dict.lookupOrDefault<scalar>
00720     (
00721         "integralAdjTol",
00722         integralAdjTol_
00723     );
00724     edgeMasterCatchFraction_ = dict.lookupOrDefault<scalar>
00725     (
00726         "edgeMasterCatchFraction",
00727         edgeMasterCatchFraction_
00728     );
00729     edgeCoPlanarTol_ = dict.lookupOrDefault<scalar>
00730     (
00731         "edgeCoPlanarTol",
00732         edgeCoPlanarTol_
00733     );
00734     edgeEndCutoffTol_ = dict.lookupOrDefault<scalar>
00735     (
00736         "edgeEndCutoffTol",
00737         edgeEndCutoffTol_
00738     );
00739 
00740     if (report)
00741     {
00742         Info<< "Sliding interface parameters:" << nl
00743             << "pointMergeTol            : " << pointMergeTol_ << nl
00744             << "edgeMergeTol             : " << edgeMergeTol_ << nl
00745             << "nFacesPerSlaveEdge       : " << nFacesPerSlaveEdge_ << nl
00746             << "edgeFaceEscapeLimit      : " << edgeFaceEscapeLimit_ << nl
00747             << "integralAdjTol           : " << integralAdjTol_ << nl
00748             << "edgeMasterCatchFraction  : " << edgeMasterCatchFraction_ << nl
00749             << "edgeCoPlanarTol          : " << edgeCoPlanarTol_ << nl
00750             << "edgeEndCutoffTol         : " << edgeEndCutoffTol_ << endl;
00751     }
00752 }
00753 
00754 
00755 void Foam::slidingInterface::write(Ostream& os) const
00756 {
00757     os  << nl << type() << nl
00758         << name()<< nl
00759         << masterFaceZoneID_.name() << nl
00760         << slaveFaceZoneID_.name() << nl
00761         << cutPointZoneID_.name() << nl
00762         << cutFaceZoneID_.name() << nl
00763         << masterPatchID_.name() << nl
00764         << slavePatchID_.name() << nl
00765         << typeOfMatchNames_[matchType_] << nl
00766         << coupleDecouple_ << nl
00767         << attached_ << endl;
00768 }
00769 
00770 
00771 // To write out all those tolerances
00772 #define WRITE_NON_DEFAULT(name) \
00773     if( name ## _ != name ## Default_ )\
00774     { \
00775         os << "    " #name " " <<  name ## _ << token::END_STATEMENT << nl; \
00776     }
00777 
00778 
00779 void Foam::slidingInterface::writeDict(Ostream& os) const
00780 {
00781     os  << nl << name() << nl << token::BEGIN_BLOCK << nl
00782         << "    type " << type() << token::END_STATEMENT << nl
00783         << "    masterFaceZoneName " << masterFaceZoneID_.name()
00784         << token::END_STATEMENT << nl
00785         << "    slaveFaceZoneName " << slaveFaceZoneID_.name()
00786         << token::END_STATEMENT << nl
00787         << "    cutPointZoneName " << cutPointZoneID_.name()
00788         << token::END_STATEMENT << nl
00789         << "    cutFaceZoneName " << cutFaceZoneID_.name()
00790         << token::END_STATEMENT << nl
00791         << "    masterPatchName " << masterPatchID_.name()
00792         << token::END_STATEMENT << nl
00793         << "    slavePatchName " << slavePatchID_.name()
00794         << token::END_STATEMENT << nl
00795         << "    typeOfMatch " << typeOfMatchNames_[matchType_]
00796         << token::END_STATEMENT << nl
00797         << "    coupleDecouple " << coupleDecouple_
00798         << token::END_STATEMENT << nl
00799         << "    projection " << intersection::algorithmNames_[projectionAlgo_]
00800         << token::END_STATEMENT << nl
00801         << "    attached " << attached_
00802         << token::END_STATEMENT << nl
00803         << "    active " << active()
00804         << token::END_STATEMENT << nl;
00805 
00806     if (attached_)
00807     {
00808         masterFaceCellsPtr_->writeEntry("masterFaceCells", os);
00809         slaveFaceCellsPtr_->writeEntry("slaveFaceCells", os);
00810         masterStickOutFacesPtr_->writeEntry("masterStickOutFaces", os);
00811         slaveStickOutFacesPtr_->writeEntry("slaveStickOutFaces", os);
00812 
00813          os << "    retiredPointMap " << retiredPointMap()
00814             << token::END_STATEMENT << nl
00815             << "    cutPointEdgePairMap " << cutPointEdgePairMap()
00816             << token::END_STATEMENT << nl;
00817     }
00818 
00819     WRITE_NON_DEFAULT(pointMergeTol)
00820     WRITE_NON_DEFAULT(edgeMergeTol)
00821     WRITE_NON_DEFAULT(nFacesPerSlaveEdge)
00822     WRITE_NON_DEFAULT(edgeFaceEscapeLimit)
00823     WRITE_NON_DEFAULT(integralAdjTol)
00824     WRITE_NON_DEFAULT(edgeMasterCatchFraction)
00825     WRITE_NON_DEFAULT(edgeCoPlanarTol)
00826     WRITE_NON_DEFAULT(edgeEndCutoffTol)
00827 
00828     os  << token::END_BLOCK << endl;
00829 }
00830 
00831 
00832 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines