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

createCoupleMatches.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 Description
00025     Create coupled match faces and add them to the cells
00026 
00027 \*---------------------------------------------------------------------------*/
00028 
00029 #include "starMesh.H"
00030 #include <OpenFOAM/boolList.H>
00031 #include <OpenFOAM/pointHit.H>
00032 #include <OpenFOAM/IOmanip.H>
00033 #include <OpenFOAM/boundBox.H>
00034 #include <OpenFOAM/Map.H>
00035 #include <OpenFOAM/mathematicalConstants.H>
00036 
00037 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00038 
00039 void starMesh::createCoupleMatches()
00040 {
00041     // Loop through all couples and create intersection faces. Add all points
00042     // of intersection faces to the couple points lists. The numbering of
00043     // the list is set such that the list can be appended to the
00044     // existing points list
00045 
00046     // Estimate the number of cells affected by couple matches
00047     const label cellMapSize = min
00048     (
00049         cellShapes_.size()/10,
00050         couples_.size()*2
00051     );
00052 
00053     // Store newly created faces for each cell
00054     Map<SLList<face> > cellAddedFaces(cellMapSize);
00055 
00056     Map<SLList<label> > cellRemovedFaces(cellMapSize);
00057 
00058     // In order to remove often allocation, remember the number of live points.
00059     // If you run out of space in point creation, increase it by the number of
00060     // couples (good scale) and resize at the end;
00061     label nLivePoints = points_.size();
00062 
00063     const label infoJump = max(1000, couples_.size()/20);
00064 
00065     forAll (couples_, coupleI)
00066     {
00067         if (coupleI % infoJump == 0)
00068         {
00069             Info << "Doing couple " << coupleI << ". STAR couple ID: "
00070                 << couples_[coupleI].coupleID() << endl;
00071         }
00072 
00073         // Initialise cell edges for master and slave cells
00074         const coupledFacePair& fp = couples_[coupleI];
00075         const face& masterFace = cellFaces_[fp.masterCell()][fp.masterFace()];
00076         const face& slaveFace = cellFaces_[fp.slaveCell()][fp.slaveFace()];
00077 
00078 #       ifdef DEBUG_COUPLE
00079         Info<< "coupleI: " << coupleI << endl
00080             << "masterFace: " << masterFace << endl
00081             << "master points: " << masterFace.points(points_) << endl
00082             << "slaveFace: " << slaveFace << endl
00083             << "slave points: " << slaveFace.points(points_)
00084             << endl << endl;
00085 #       endif
00086 
00087         // check the angle of face area vectors
00088          scalar faceAreaAngle =
00089              mag
00090              (
00091                  -(masterFace.normal(points_) & slaveFace.normal(points_))/
00092                  (masterFace.mag(points_)*slaveFace.mag(points_) + VSMALL)
00093              );
00094 
00095          if (faceAreaAngle < 0.94)
00096          {
00097              Info<< "Couple direction mismatch in the couple match "
00098                  << coupleI << ". STAR couple ID: "
00099                  << couples_[coupleI].coupleID() << endl
00100                  << "The angle between face normals is "
00101                  << Foam::acos(faceAreaAngle)/mathematicalConstant::pi*180
00102                  << " deg." << endl
00103                  << "master cell: " << fp.masterCell()
00104                  << " STAR number: " << starCellID_[fp.masterCell()]
00105                  << " type: " << cellShapes_[fp.masterCell()].model().name()
00106                  << " face: " << fp.masterFace() << endl
00107                  << "slave cell : " << fp.slaveCell()
00108                  << " STAR number: " << starCellID_[fp.slaveCell()]
00109                  << " type: " << cellShapes_[fp.slaveCell()].model().name()
00110                  << " face: " << fp.slaveFace() << endl;
00111          }
00112 
00113         // Deal with integral patches
00114         if (fp.integralMatch())
00115         {
00116             // Master face is replaced by a set of slave faces
00117 
00118             Map<SLList<label> >::iterator crfIter =
00119                 cellRemovedFaces.find(fp.masterCell());
00120 
00121             if (crfIter == cellRemovedFaces.end())
00122             {
00123                 cellRemovedFaces.insert
00124                 (
00125                     fp.masterCell(),
00126                     fp.masterFace()
00127                 );
00128             }
00129             else
00130             {
00131                 crfIter().append(fp.masterFace());
00132             }
00133 
00134             Map<SLList<face> >::iterator cafIter =
00135                 cellAddedFaces.find(fp.masterCell());
00136             if (cafIter == cellAddedFaces.end())
00137             {
00138                 cellAddedFaces.insert
00139                 (
00140                     fp.masterCell(),
00141                     SLList<face>(slaveFace.reverseFace())
00142                 );
00143             }
00144             else
00145             {
00146                 cafIter().append(slaveFace.reverseFace());
00147             }
00148         }
00149         else
00150         {
00151             // Create cut faces, which replace both master and slave faces
00152 
00153             // Store newly created points
00154             SLList<point> coupleFacePoints;
00155 
00156             // Master data
00157             edgeList masterEdges = masterFace.edges();
00158             List<SLList<label> > masterEdgePoints(masterEdges.size());
00159 
00160             // Slave data
00161             edgeList slaveEdges = slaveFace.edges();
00162             List<SLList<label> > slaveEdgePoints(slaveEdges.size());
00163 
00164             // Find common plane
00165             vector n = masterFace.normal(points_);
00166             n /= mag(n) + VSMALL;
00167 
00168             // Loop through all edges of the master face. For every edge,
00169             // intersect it with all edges of the cutting face.
00170             forAll (masterEdges, masterEdgeI)
00171             {
00172                 const edge& curMasterEdge = masterEdges[masterEdgeI];
00173 
00174                 point P = points_[curMasterEdge.start()];
00175 
00176                 // get d and return it into plane
00177                 vector d = curMasterEdge.vec(points_);
00178                 d -= n*(n & d);
00179 
00180 #               ifdef DEBUG_COUPLE_INTERSECTION
00181                 Info << "curMasterEdge: " << curMasterEdge << endl
00182                     << "P: " << P << endl << "d: " << d << endl;
00183 #               endif
00184 
00185                 // go through all slave edges and try to get an intersection.
00186                 // The point is created along the original master edge rather
00187                 // than its corrected direction.
00188                 forAll (slaveEdges, slaveEdgeI)
00189                 {
00190                     const edge& curSlaveEdge = slaveEdges[slaveEdgeI];
00191 
00192                     point S = points_[curSlaveEdge.start()];
00193 
00194                     // get e and return it into plane
00195                     vector e = curSlaveEdge.vec(points_);
00196                     e -= n*(n & e);
00197                     scalar det = -(e & (n ^ d));
00198 
00199 #                   ifdef DEBUG_COUPLE_INTERSECTION
00200                     Info << "curSlaveEdge: " << curSlaveEdge << endl
00201                         << "S: " << S << endl
00202                         << "e: " << e << endl;
00203 #                   endif
00204 
00205                     if (mag(det) > SMALL)
00206                     {
00207                         // non-singular matrix. Look for intersection
00208                         scalar beta = ((S - P) & (n ^ d))/det;
00209 
00210 #                       ifdef DEBUG_COUPLE_INTERSECTION
00211                         Info << " beta: " << beta << endl;
00212 #                       endif
00213 
00214                         if (beta > -smallMergeTol_ && beta < 1 + smallMergeTol_)
00215                         {
00216                             // slave intersection OK. Try master intersection
00217                             scalar alpha =
00218                                 (((S - P) & d) + beta*(d & e))/magSqr(d);
00219 
00220 #                           ifdef DEBUG_COUPLE_INTERSECTION
00221                             Info << " alpha: " << alpha << endl;
00222 #                           endif
00223 
00224                             if
00225                             (
00226                                 alpha > -smallMergeTol_
00227                              && alpha < 1 + smallMergeTol_
00228                             )
00229                             {
00230                                 // intersection of non-parallel edges
00231 #                               ifdef DEBUG_COUPLE_INTERSECTION
00232                                 Info<< "intersection of non-parallel edges"
00233                                     << endl;
00234 #                               endif
00235 
00236 
00237                                 // check for insertion of start-end
00238                                 // points in the middle of the other
00239                                 // edge
00240                                 if (alpha < smallMergeTol_)
00241                                 {
00242                                     // inserting the start of master edge
00243                                     if
00244                                     (
00245                                         beta > smallMergeTol_
00246                                      && beta < 1 - smallMergeTol_
00247                                     )
00248                                     {
00249                                         slaveEdgePoints[slaveEdgeI].append
00250                                         (
00251                                             curMasterEdge.start()
00252                                         );
00253                                     }
00254                                 }
00255                                 else if (alpha > 1 - smallMergeTol_)
00256                                 {
00257                                     // inserting the end of master edge
00258                                     if
00259                                     (
00260                                         beta > smallMergeTol_
00261                                      && beta < 1 - smallMergeTol_
00262                                     )
00263                                     {
00264                                         slaveEdgePoints[slaveEdgeI].append
00265                                         (
00266                                             curMasterEdge.end()
00267                                         );
00268                                     }
00269                                 }
00270                                 else if (beta < smallMergeTol_)
00271                                 {
00272                                     // inserting the start of the slave edge
00273                                     if
00274                                     (
00275                                         alpha > smallMergeTol_
00276                                      && alpha < 1 - smallMergeTol_
00277                                     )
00278                                     {
00279                                         masterEdgePoints[masterEdgeI].append
00280                                         (
00281                                             curSlaveEdge.start()
00282                                         );
00283                                     }
00284                                 }
00285                                 else if (beta > 1 - smallMergeTol_)
00286                                 {
00287                                     // inserting the start of the slave edge
00288                                     if
00289                                     (
00290                                         alpha > smallMergeTol_
00291                                      && alpha < 1 - smallMergeTol_
00292                                     )
00293                                     {
00294                                         masterEdgePoints[masterEdgeI].append
00295                                         (
00296                                             curSlaveEdge.end()
00297                                         );
00298                                     }
00299                                 }
00300                                 else
00301                                 {
00302                                     masterEdgePoints[masterEdgeI].append
00303                                     (
00304                                         nLivePoints + coupleFacePoints.size()
00305                                     );
00306 
00307                                     slaveEdgePoints[slaveEdgeI].append
00308                                     (
00309                                         nLivePoints + coupleFacePoints.size()
00310                                     );
00311 
00312 #                                   ifdef DEBUG_COUPLE_INTERSECTION
00313                                     Info<< "regular intersection. "
00314                                         << "Adding point: "
00315                                         << coupleFacePoints.size()
00316                                         << " which is "
00317                                         << P + alpha*curMasterEdge.vec(points_)
00318                                         << endl;
00319 #                                   endif
00320 
00321                                     // A new point is created. Warning:
00322                                     // using original edge for accuracy.
00323                                     //
00324                                     coupleFacePoints.append
00325                                         (P + alpha*curMasterEdge.vec(points_));
00326                                 }
00327                             }
00328                         }
00329                     }
00330                     else
00331                     {
00332                         // Add special cases, for intersection of two
00333                         // parallel line Warning. Here, typically, no new
00334                         // points will be created. Either one or two of
00335                         // the slave edge points need to be added to the
00336                         // master edge and vice versa. The problem is that
00337                         // no symmetry exists, i.e. both operations needs
00338                         // to be done separately for both master and slave
00339                         // side.
00340 
00341                         // Master side
00342                         // check if the first or second point of slave edge is
00343                         // on the master edge
00344                         vector ps = S - P;
00345 
00346                         bool colinear = false;
00347 
00348                         if (mag(ps) < SMALL)
00349                         {
00350                             // colinear because P and S are the same point
00351                             colinear = true;
00352                         }
00353                         else if
00354                         (
00355                             (ps & d)/(mag(ps)*mag(d)) > 1.0 - smallMergeTol_
00356                         )
00357                         {
00358                             // colinear because ps and d are parallel
00359                             colinear = true;
00360                         }
00361 
00362                         if (colinear)
00363                         {
00364                             scalar alpha1 = (ps & d)/magSqr(d);
00365 
00366                             if
00367                             (
00368                                 alpha1 > -smallMergeTol_
00369                              && alpha1 < 1 + smallMergeTol_
00370                             )
00371                             {
00372 #                               ifdef DEBUG_COUPLE_INTERSECTION
00373                                 Info<< "adding irregular master "
00374                                     << "intersection1: "
00375                                     << points_[slaveEdges[slaveEdgeI].start()]
00376                                     << endl;
00377 #                               endif
00378 
00379                                 masterEdgePoints[masterEdgeI].append
00380                                 (
00381                                     slaveEdges[slaveEdgeI].start()
00382                                 );
00383                             }
00384 
00385                              scalar alpha2 = ((ps + e) & d)/magSqr(d);
00386 
00387                             if
00388                             (
00389                                 alpha2 > -smallMergeTol_
00390                              && alpha2 < 1 + smallMergeTol_
00391                             )
00392                             {
00393 #                               ifdef DEBUG_COUPLE_INTERSECTION
00394                                 Info<< "adding irregular master "
00395                                     << "intersection2: "
00396                                     << points_[slaveEdges[slaveEdgeI].end()]
00397                                     << endl;
00398 #                               endif
00399 
00400                                 masterEdgePoints[masterEdgeI].append
00401                                 (
00402                                     slaveEdges[slaveEdgeI].end()
00403                                 );
00404                             }
00405 
00406                             // Slave side
00407                             // check if the first or second point of
00408                             // master edge is on the slave edge
00409 
00410                             vector sp = P - S;
00411 
00412                             scalar beta1 = (sp & e)/magSqr(e);
00413 
00414 #                           ifdef DEBUG_COUPLE_INTERSECTION
00415                             Info << "P: " << P << " S: " << S << " d: " << d
00416                                 << " e: " << e << " sp: " << sp
00417                                 << " beta1: " << beta1 << endl;
00418 #                           endif
00419 
00420                             if
00421                             (
00422                                 beta1 > -smallMergeTol_
00423                              && beta1 < 1 + smallMergeTol_
00424                             )
00425                             {
00426 #                               ifdef DEBUG_COUPLE_INTERSECTION
00427                                 Info<< "adding irregular slave "
00428                                     << "intersection1: "
00429                                     << points_[masterEdges[masterEdgeI].start()]
00430                                     << endl;
00431 #                               endif
00432 
00433                                 slaveEdgePoints[slaveEdgeI].append
00434                                 (
00435                                     masterEdges[masterEdgeI].start()
00436                                 );
00437                             }
00438 
00439                             scalar beta2 = ((sp + d) & e)/magSqr(e);
00440 
00441                             if
00442                             (
00443                                 beta2 > -smallMergeTol_
00444                              && beta2 < 1 + smallMergeTol_
00445                             )
00446                             {
00447 #                               ifdef DEBUG_COUPLE_INTERSECTION
00448                                 Info << "adding irregular slave "
00449                                     << "intersection2: "
00450                                     << points_[masterEdges[masterEdgeI].end()]
00451                                     << endl;
00452 #                               endif
00453 
00454                                 slaveEdgePoints[slaveEdgeI].append
00455                                 (
00456                                     masterEdges[masterEdgeI].end()
00457                                 );
00458                             }
00459                         } // end of colinear
00460                     } // end of singular intersection
00461                 } // end of slave edges
00462             } // end of master edges
00463 
00464 #           ifdef DEBUG_COUPLE_INTERSECTION
00465             Info << "additional slave edge points: " << endl;
00466             forAll (slaveEdgePoints, edgeI)
00467             {
00468                 Info << "edge: " << edgeI << ": " << slaveEdgePoints[edgeI]
00469                     << endl;
00470             }
00471 #           endif
00472 
00473             // Add new points
00474             if (nLivePoints + coupleFacePoints.size() >= points_.size())
00475             {
00476                 // increase the size of the points list
00477                 Info << "Resizing points list" << endl;
00478                 points_.setSize(points_.size() + couples_.size());
00479             }
00480 
00481             for
00482             (
00483                 SLList<point>::iterator coupleFacePointsIter =
00484                     coupleFacePoints.begin();
00485                 coupleFacePointsIter != coupleFacePoints.end();
00486                 ++coupleFacePointsIter
00487             )
00488             {
00489                 points_[nLivePoints] = coupleFacePointsIter();
00490                 nLivePoints++;
00491             }
00492 
00493             // edge intersection finished
00494 
00495             // Creating new master side
00496 
00497             // count the number of additional points for face
00498             label nAdditionalMasterPoints = 0;
00499 
00500             forAll (masterEdgePoints, edgeI)
00501             {
00502                 nAdditionalMasterPoints += masterEdgePoints[edgeI].size();
00503             }
00504 
00505             face tmpMasterFace
00506             (
00507                 masterFace.size()
00508               + nAdditionalMasterPoints
00509             );
00510             label nTmpMasterLabels = 0;
00511 
00512 #           ifdef DEBUG_COUPLE_INTERSECTION
00513             Info << "masterFace: " << masterFace << endl
00514                 << "nAdditionalMasterPoints: " << nAdditionalMasterPoints
00515                 << endl;
00516 #           endif
00517 
00518             forAll (masterEdges, masterEdgeI)
00519             {
00520                 // Insert the starting point of the edge
00521                 tmpMasterFace[nTmpMasterLabels] =
00522                     masterEdges[masterEdgeI].start();
00523                 nTmpMasterLabels++;
00524 
00525                 // get reference to added points of current edge
00526                 const SLList<label>& curMEdgePoints =
00527                     masterEdgePoints[masterEdgeI];
00528 
00529                 // create a markup list of points that have been used
00530                 boolList usedMasterPoint(curMEdgePoints.size(), false);
00531 
00532                 vector edgeVector = masterEdges[masterEdgeI].vec(points_);
00533 
00534 #               ifdef DEBUG_FACE_ORDERING
00535                 Info<< "edgeVector: " << edgeVector << endl
00536                     << "curMEdgePoints.size(): " << curMEdgePoints.size()
00537                     << endl;
00538 #               endif
00539 
00540                 // renormalise
00541                 edgeVector /= magSqr(edgeVector);
00542 
00543                 point edgeStartPoint =
00544                     points_[masterEdges[masterEdgeI].start()];
00545 
00546                 // loop until the next label to add is -1
00547                 for(;;)
00548                 {
00549                     label nextPointLabel = -1;
00550                     label usedI = -1;
00551                     scalar minAlpha = GREAT;
00552 
00553                     label i = 0;
00554 
00555                     for
00556                     (
00557                         SLList<label>::const_iterator curMEdgePointsIter =
00558                             curMEdgePoints.begin();
00559                         curMEdgePointsIter != curMEdgePoints.end();
00560                         ++curMEdgePointsIter
00561                     )
00562                     {
00563                         if (!usedMasterPoint[i])
00564                         {
00565                             scalar alpha =
00566                                 edgeVector
00567                               & (
00568                                     points_[curMEdgePointsIter()]
00569                                   - edgeStartPoint
00570                                 );
00571 
00572 #                           ifdef DEBUG_FACE_ORDERING
00573                             Info<< " edgeStartPoint: " << edgeStartPoint
00574                                 << " edgeEndPoint: "
00575                                 << points_[masterEdges[masterEdgeI].end()]
00576                                 << " other point: "
00577                                 << points_[curMEdgePointsIter()]
00578                                 << " alpha: " << alpha << endl;
00579 #                           endif
00580 
00581                             if (alpha < minAlpha)
00582                             {
00583                                 minAlpha = alpha;
00584                                 usedI = i;
00585                                 nextPointLabel = curMEdgePointsIter();
00586                             }
00587                         }
00588 
00589 #                       ifdef DEBUG_FACE_ORDERING
00590                         Info << "nextPointLabel: " << nextPointLabel << endl;
00591 #                       endif
00592 
00593                         i++;
00594                     }
00595 
00596                     if (nextPointLabel > -1)
00597                     {
00598 #                       ifdef DEBUG_FACE_ORDERING
00599                         Info<< "added nextPointLabel: " << nextPointLabel
00600                             << " nTmpMasterLabels: " << nTmpMasterLabels
00601                             << " to place " << nTmpMasterLabels << endl;
00602 #                       endif
00603 
00604                         usedMasterPoint[usedI] = true;
00605                         // add the next point
00606                         tmpMasterFace[nTmpMasterLabels] =
00607                             nextPointLabel;
00608                         nTmpMasterLabels++;
00609                     }
00610                     else
00611                     {
00612                         break;
00613                     }
00614                 }
00615             }
00616 
00617             // reset the size of master
00618             tmpMasterFace.setSize(nTmpMasterLabels);
00619 
00620 #           ifdef DEBUG_FACE_ORDERING
00621             Info << "tmpMasterFace: " << tmpMasterFace << endl;
00622 #           endif
00623 
00624             // Eliminate all zero-length edges
00625             face newMasterFace(labelList(tmpMasterFace.size(), labelMax));
00626 
00627             // insert first point by hand. Careful: the first one is
00628             // used for comparison to allow the edge collapse across
00629             // point zero.
00630             //
00631             newMasterFace[0] = tmpMasterFace[0];
00632             label nMaster = 0;
00633 
00634             edgeList mstEdgesToCollapse = tmpMasterFace.edges();
00635 
00636             scalar masterTol =
00637                 cpMergePointTol_*boundBox(tmpMasterFace.points(points_)).mag();
00638 
00639             forAll (mstEdgesToCollapse, edgeI)
00640             {
00641 #               ifdef DEBUG_FACE_ORDERING
00642                 Info<< "edgeI: " << edgeI << " curEdge: "
00643                     << mstEdgesToCollapse[edgeI] << endl
00644                     << "master edge " << edgeI << ", "
00645                     << mstEdgesToCollapse[edgeI].mag(points_) << endl;
00646 #               endif
00647 
00648                 // Edge merge tolerance = masterTol
00649                 if (mstEdgesToCollapse[edgeI].mag(points_) < masterTol)
00650                 {
00651                     newMasterFace[nMaster] =
00652                         min
00653                         (
00654                             newMasterFace[nMaster],
00655                             mstEdgesToCollapse[edgeI].end()
00656                         );
00657 
00658 #                   ifdef DEBUG_FACE_ORDERING
00659                     Info << "Collapsed: nMaster: " << nMaster
00660                         << " label: " << newMasterFace[nMaster] << endl;
00661 #                   endif
00662 
00663                 }
00664                 else
00665                 {
00666                     nMaster++;
00667 
00668                     if (edgeI < mstEdgesToCollapse.size() - 1)
00669                     {
00670                         // last edge does not add the point
00671 #                   ifdef DEBUG_FACE_ORDERING
00672                         Info<< "Added: nMaster: " << nMaster
00673                             << " label: " << mstEdgesToCollapse[edgeI].end()
00674                             << endl;
00675 #                   endif
00676 
00677                         newMasterFace[nMaster] =
00678                             mstEdgesToCollapse[edgeI].end();
00679                     }
00680                 }
00681             }
00682 
00683             newMasterFace.setSize(nMaster);
00684 
00685 #           ifdef DEBUG_COUPLE
00686             Info<< "newMasterFace: " << newMasterFace << endl
00687                 << "points: " << newMasterFace.points(points_) << endl;
00688 #           endif
00689 
00690             // Creating new slave side
00691 
00692             // count the number of additional points for face
00693             label nAdditionalSlavePoints = 0;
00694 
00695             forAll (slaveEdgePoints, edgeI)
00696             {
00697                 nAdditionalSlavePoints += slaveEdgePoints[edgeI].size();
00698             }
00699 
00700             face tmpSlaveFace
00701             (
00702                 slaveFace.size()
00703               + nAdditionalSlavePoints
00704             );
00705             label nTmpSlaveLabels = 0;
00706 
00707 #           ifdef DEBUG_COUPLE_INTERSECTION
00708             Info<< "slaveFace: " << slaveFace << endl
00709                 << "nAdditionalSlavePoints: " << nAdditionalSlavePoints << endl;
00710 #           endif
00711 
00712             forAll (slaveEdges, slaveEdgeI)
00713             {
00714                 // Insert the starting point of the edge
00715                 tmpSlaveFace[nTmpSlaveLabels] =
00716                     slaveEdges[slaveEdgeI].start();
00717                 nTmpSlaveLabels++;
00718 
00719                 // get reference to added points of current edge
00720                 const SLList<label>& curSEdgePoints =
00721                     slaveEdgePoints[slaveEdgeI];
00722 
00723                 // create a markup list of points that have been used
00724                 boolList usedSlavePoint(curSEdgePoints.size(), false);
00725 
00726                 vector edgeVector = slaveEdges[slaveEdgeI].vec(points_);
00727 
00728 #               ifdef DEBUG_FACE_ORDERING
00729                 Info << "curSEdgePoints.size(): "
00730                     << curSEdgePoints.size() << endl
00731                     << "edgeVector: " << edgeVector << endl;
00732 #               endif
00733 
00734                 // renormalise
00735                 edgeVector /= magSqr(edgeVector);
00736 
00737                 point edgeStartPoint =
00738                     points_[slaveEdges[slaveEdgeI].start()];
00739 
00740                 // loop until the next label to add is -1
00741                 for(;;)
00742                 {
00743                     label nextPointLabel = -1;
00744                     label usedI = -1;
00745                     scalar minAlpha = GREAT;
00746 
00747                     label i = 0;
00748 
00749                     for
00750                     (
00751                         SLList<label>::const_iterator curSEdgePointsIter =
00752                             curSEdgePoints.begin();
00753                         curSEdgePointsIter != curSEdgePoints.end();
00754                         ++curSEdgePointsIter
00755                     )
00756                     {
00757                         if (!usedSlavePoint[i])
00758                         {
00759                             scalar alpha =
00760                                 edgeVector
00761                               & (
00762                                     points_[curSEdgePointsIter()]
00763                                   - edgeStartPoint
00764                                 );
00765 
00766 #                           ifdef DEBUG_FACE_ORDERING
00767                             Info<< " edgeStartPoint: " << edgeStartPoint
00768                                 << " edgeEndPoint: "
00769                                 << points_[slaveEdges[slaveEdgeI].end()]
00770                                 << " other point: "
00771                                 << points_[curSEdgePointsIter()]
00772                                 << " alpha: " << alpha << endl;
00773 #                           endif
00774 
00775                             if (alpha < minAlpha)
00776                             {
00777                                 minAlpha = alpha;
00778                                 usedI = i;
00779                                 nextPointLabel = curSEdgePointsIter();
00780                             }
00781                         }
00782 
00783 #                       ifdef DEBUG_FACE_ORDERING
00784                         Info << "nextPointLabel: " << nextPointLabel << endl;
00785 #                       endif
00786 
00787                         i++;
00788                     }
00789 
00790                     if (nextPointLabel > -1)
00791                     {
00792 #                       ifdef DEBUG_FACE_ORDERING
00793                         Info<< "added nextPointLabel: " << nextPointLabel
00794                             << " nTmpSlaveLabels: " << nTmpSlaveLabels
00795                             << " to place " << nTmpSlaveLabels << endl;
00796 #                       endif
00797 
00798                         usedSlavePoint[usedI] = true;
00799                         // add the next point
00800                         tmpSlaveFace[nTmpSlaveLabels] =
00801                             nextPointLabel;
00802                         nTmpSlaveLabels++;
00803                     }
00804                     else
00805                     {
00806                         break;
00807                     }
00808                 }
00809             }
00810 
00811             // reset the size of slave
00812             tmpSlaveFace.setSize(nTmpSlaveLabels);
00813 
00814 #           ifdef DEBUG_FACE_ORDERING
00815             Info << "tmpSlaveFace: " << tmpSlaveFace << endl;
00816 #           endif
00817 
00818             // Eliminate all zero-length edges
00819             face newSlaveFace(labelList(tmpSlaveFace.size(), labelMax));
00820 
00821             // insert first point by hand. Careful: the first one is
00822             // used for comparison to allow the edge collapse across
00823             // point zero.
00824             //
00825             newSlaveFace[0] = tmpSlaveFace[0];
00826             label nSlave = 0;
00827 
00828             edgeList slvEdgesToCollapse = tmpSlaveFace.edges();
00829 
00830             scalar slaveTol =
00831                 cpMergePointTol_*boundBox(tmpSlaveFace.points(points_)).mag();
00832 
00833             forAll(slvEdgesToCollapse, edgeI)
00834             {
00835 #               ifdef DEBUG_FACE_ORDERING
00836                 Info << "slave edge length: " << edgeI << ", "
00837                     << slvEdgesToCollapse[edgeI].mag(points_)<< endl;
00838 #               endif
00839 
00840                  // edge merge tolerance = slaveTol
00841                 if (slvEdgesToCollapse[edgeI].mag(points_) < slaveTol)
00842                 {
00843                     newSlaveFace[nSlave] =
00844                         min
00845                         (
00846                             newSlaveFace[nSlave],
00847                             slvEdgesToCollapse[edgeI].end()
00848                         );
00849                 }
00850                 else
00851                 {
00852                     nSlave++;
00853                     if (edgeI < slvEdgesToCollapse.size() - 1)
00854                     {
00855                         // last edge does not add the point
00856                         newSlaveFace[nSlave] = slvEdgesToCollapse[edgeI].end();
00857                     }
00858                 }
00859             }
00860 
00861             newSlaveFace.setSize(nSlave);
00862 
00863 #           ifdef DEBUG_COUPLE
00864             Info<< "newSlaveFace: " << newSlaveFace << endl
00865                 << "points: " << newSlaveFace.points(points_) << endl << endl;
00866 #           endif
00867 
00868             // Create the intersection face
00869 
00870             // Algorithm:
00871             // Loop through
00872             // points of the master and try to find one which falls
00873             // within the slave.  If not found, look through all
00874             // edges of the slave and find one which falls within the
00875             // master.  This point will be the starting location for
00876             // the cut face.
00877 
00878             edgeList newMasterEdges = newMasterFace.edges();
00879             edgeList newSlaveEdges = newSlaveFace.edges();
00880 
00881 #           ifdef DEBUG_RIGHT_HAND_WALK
00882             Info << "newMasterEdges: " << newMasterEdges << endl
00883                 << "newSlaveEdges: " << newSlaveEdges << endl;
00884 #           endif
00885 
00886             edge startEdge(-1, -1);
00887 
00888             // Remember where the start edge was found:
00889             // 0 for not found
00890             // 1 for master
00891             // 2 for slave
00892             label startEdgeFound = 0;
00893 
00894             vector masterProjDir = -newMasterFace.normal(points_);
00895 
00896             forAll (newSlaveEdges, edgeI)
00897             {
00898                 // Take the slave edge points and project into the master.
00899                 // In order to create a good intersection, move the
00900                 // point away from the master in the direction of its
00901                 // normal.
00902                 point pointStart = points_[newSlaveEdges[edgeI].start()];
00903 
00904                 point pointEnd = points_[newSlaveEdges[edgeI].end()];
00905 
00906                 if
00907                 (
00908                     newMasterFace.ray
00909                     (
00910                         pointStart,
00911                         masterProjDir,
00912                         points_,
00913                         intersection::FULL_RAY
00914                     ).hit()
00915                  && newMasterFace.ray
00916                     (
00917                         pointEnd,
00918                         masterProjDir,
00919                         points_,
00920                         intersection::FULL_RAY
00921                     ).hit()
00922                 )
00923                 {
00924                     startEdge = newSlaveEdges[edgeI];
00925                     startEdgeFound = 2;
00926 
00927 #                   ifdef DEBUG_RIGHT_HAND_WALK
00928                     Info << "slave edge found" << endl;
00929 #                   endif
00930 
00931                     break;
00932                 }
00933             }
00934 
00935             if (startEdgeFound == 0)
00936             {
00937                 vector slaveProjDir = -newSlaveFace.normal(points_);
00938 
00939                 forAll (newMasterEdges, edgeI)
00940                 {
00941                     // Take the edge master points and project into the slave.
00942                     // In order to create a good intersection, move the
00943                     // point away from the slave in the direction of its
00944                     // normal.
00945                     point pointStart = points_[newMasterEdges[edgeI].start()];
00946 
00947                     point pointEnd = points_[newMasterEdges[edgeI].end()];
00948 
00949                     if
00950                     (
00951                         newSlaveFace.ray
00952                         (
00953                             pointStart,
00954                             slaveProjDir,
00955                             points_,
00956                         intersection::FULL_RAY
00957                         ).hit()
00958                      && newSlaveFace.ray
00959                         (
00960                             pointEnd,
00961                             slaveProjDir,
00962                             points_,
00963                             intersection::FULL_RAY
00964                         ).hit()
00965                     )
00966                     {
00967                         startEdge = newMasterEdges[edgeI];
00968                         startEdgeFound = 1;
00969 
00970 #                       ifdef DEBUG_RIGHT_HAND_WALK
00971                         Info << "master edge found" << endl;
00972 #                       endif
00973 
00974                         break;
00975                     }
00976                 }
00977             }
00978 
00979             // create the intersected face using right-hand walk rule
00980             face intersectedFace
00981             (
00982                 labelList(newMasterFace.size() + newSlaveFace.size(), -1)
00983             );
00984 
00985             if (startEdgeFound > 0)
00986             {
00987 #               ifdef DEBUG_RIGHT_HAND_WALK
00988                 Info << "start edge: " << startEdge << endl;
00989 #               endif
00990 
00991                 // Loop through both faces and add all edges
00992                 // containing the current point and add them to the
00993                 // list of edges to consider.  Make sure all edges are
00994                 // added such that the current point is their start.
00995                 // Loop through all edges to consider and find the one
00996                 // which produces the buggest right-hand-turn.  This
00997                 // is the next edge to be added to the face.  If its
00998                 // end is the same as the starting point, the face is
00999                 // complete; resize it to the number of active points
01000                 // and exit.
01001 
01002                 vector planeNormal = newMasterFace.normal(points_);
01003                 planeNormal /= mag(planeNormal) + VSMALL;
01004 
01005 #               ifdef DEBUG_RIGHT_HAND_WALK
01006                 Info << "planeNormal: " << planeNormal << endl;
01007 #               endif
01008 
01009                 // Do a check to control the right-hand turn.  This is
01010                 // based on the triple product of the edge start
01011                 // vector to face centre, the edge vector and the
01012                 // plane normal.  If the triple product is negative,
01013                 // the edge needs to be reversed to allow the
01014                 // right-hand-turn rule to work.
01015 
01016                 vector faceCentre;
01017 
01018                 if (startEdgeFound == 1)
01019                 {
01020                     faceCentre = newMasterFace.centre(points_);
01021                 }
01022                 else
01023                 {
01024                     faceCentre = newSlaveFace.centre(points_);
01025                 }
01026 
01027                 scalar tripleProduct =
01028                     (
01029                         (faceCentre - points_[startEdge.start()])
01030                       ^ startEdge.vec(points_)
01031                     ) & planeNormal;
01032 
01033                 if (tripleProduct < 0)
01034                 {
01035 #                   ifdef DEBUG_RIGHT_HAND_WALK
01036                     Info << "Turning edge for right-hand turn rule" << endl;
01037 #                   endif
01038                     startEdge = startEdge.reverseEdge();
01039                 }
01040 
01041                 // prepare the loop for the right-hand walk
01042                 intersectedFace[0] = startEdge.start();
01043                 intersectedFace[1] = startEdge.end();
01044                 label nIntFacePoints = 2;
01045 
01046                 edge curEdge = startEdge;
01047 
01048                 bool completedFace = false;
01049 
01050                 do
01051                 {
01052                     SLList<edge> edgesToConsider;
01053 
01054                     // collect master edges
01055                     forAll (newMasterEdges, edgeI)
01056                     {
01057                         const edge& cme = newMasterEdges[edgeI];
01058 
01059                         if (cme != curEdge)
01060                         {
01061                             if (cme.start() == curEdge.end())
01062                             {
01063                                 edgesToConsider.append(cme);
01064                             }
01065                             else if (cme.end() == curEdge.end())
01066                             {
01067                                 edgesToConsider.append(cme.reverseEdge());
01068                             }
01069                             // otherwise, it does not have the current point
01070                         }
01071                     }
01072 
01073                     // collect slave edges
01074                     forAll (newSlaveEdges, edgeI)
01075                     {
01076                         const edge& cse = newSlaveEdges[edgeI];
01077 
01078                         if (cse != curEdge)
01079                         {
01080                             if (cse.start() == curEdge.end())
01081                             {
01082                                 edgesToConsider.append(cse);
01083                             }
01084                             else if (cse.end() == curEdge.end())
01085                             {
01086                                 edgesToConsider.append(cse.reverseEdge());
01087                             }
01088                             // otherwise, it does not have the current point
01089                         }
01090                     }
01091 
01092 #                   ifdef DEBUG_RIGHT_HAND_WALK
01093                     Info<< "number of edges to consider: "
01094                         << edgesToConsider.size() << endl
01095                         << "edges to consider: " << edgesToConsider << endl;
01096 #                   endif
01097 
01098                     if (edgesToConsider.empty())
01099                     {
01100                         FatalErrorIn("void starMesh::createCoupleMatches()")
01101                             << setprecision(12)
01102                             << "void starMesh::createCoupleMatches() : "
01103                             << endl << "error in face intersection: "
01104                             << "no edges to consider for closing the loop"
01105                             << coupleI << ". STAR couple ID: "
01106                             << couples_[coupleI].coupleID() << endl
01107                             << "Cut Master Face: " << newMasterFace << endl
01108                             << "points: " << newMasterFace.points(points_)
01109                             << endl
01110                             << "Cut Slave Face: " << newSlaveFace << endl
01111                             << "points: " << newSlaveFace.points(points_)
01112                             << endl << "intersected face: "
01113                             << intersectedFace
01114                             << abort(FatalError);
01115                     }
01116 
01117                     // vector along the edge
01118                     vector ahead = curEdge.vec(points_);
01119                     ahead -= planeNormal*(planeNormal & ahead);
01120                     ahead /= mag(ahead) + VSMALL;
01121 
01122                     // vector pointing right
01123                     vector right = ahead ^ planeNormal;
01124                     right /= mag(right) + VSMALL;
01125 
01126                     // first edge taken for reference
01127                     edge nextEdge = edgesToConsider.first();
01128                     vector nextEdgeVec = nextEdge.vec(points_);
01129                     nextEdgeVec -= planeNormal*(planeNormal & nextEdgeVec);
01130                     nextEdgeVec /= mag(nextEdgeVec) + VSMALL;
01131 
01132                     scalar rightTurn = nextEdgeVec & right;
01133                     scalar goStraight = nextEdgeVec & ahead;
01134 
01135 #                   ifdef DEBUG_RIGHT_HAND_WALK
01136                     Info<< "rightTurn: " << rightTurn
01137                         << " goStraight: " << goStraight << endl;
01138 #                   endif
01139 
01140                     for
01141                     (
01142                         SLList<edge>::iterator etcIter =
01143                             edgesToConsider.begin();
01144                         etcIter != edgesToConsider.end();
01145                         ++etcIter
01146                     )
01147                     {
01148                         // right-hand walk rule
01149                         vector newDir = etcIter().vec(points_);
01150                         newDir -= planeNormal*(planeNormal & newDir);
01151                         newDir /= mag(newDir) + VSMALL;
01152 
01153                         scalar curRightTurn = newDir & right;
01154                         scalar curGoStraight = newDir & ahead;
01155 
01156 #                       ifdef DEBUG_RIGHT_HAND_WALK
01157                         Info << "curRightTurn: " << curRightTurn
01158                             << " curGoStraight: " << curGoStraight << endl;
01159 #                       endif
01160 
01161                         if (rightTurn < 0) // old edge turning left
01162                         {
01163                             if (curRightTurn < 0) // new edge turning left
01164                             {
01165                                 // both go left. Grab one with greater ahead
01166                                 if (curGoStraight > goStraight)
01167                                 {
01168 #                                   ifdef DEBUG_RIGHT_HAND_WALK
01169                                     Info << "a" << endl;
01170 #                                   endif
01171 
01172                                     // Good edge, turning left less than before
01173                                     nextEdge = etcIter();
01174                                     rightTurn = curRightTurn;
01175                                     goStraight = curGoStraight;
01176                                 }
01177                             }
01178                             else // new edge turning right
01179                             {
01180 #                               ifdef DEBUG_RIGHT_HAND_WALK
01181                                 Info << "b" << endl;
01182 #                               endif
01183 
01184                                 // good edge, turning right
01185                                 nextEdge = etcIter();
01186                                 rightTurn = curRightTurn;
01187                                 goStraight = curGoStraight;
01188                             }
01189                         }
01190                         else // old edge turning right
01191                         {
01192                             // new edge turning left rejected
01193                             if (curRightTurn >= 0) // new edge turning right
01194                             {
01195                                 // grab one with smaller ahead
01196                                 if (curGoStraight < goStraight)
01197                                 {
01198 #                                   ifdef DEBUG_RIGHT_HAND_WALK
01199                                     Info << "c" << endl;
01200 #                                   endif
01201 
01202                                     // good edge, turning right more than before
01203                                     nextEdge = etcIter();
01204                                     rightTurn = curRightTurn;
01205                                     goStraight = curGoStraight;
01206                                 }
01207                             }
01208                         }
01209                     }
01210 
01211                     // check if the loop is completed
01212                     if (nextEdge.end() == intersectedFace[0])
01213                     {
01214                         // loop is completed. No point to add
01215                         completedFace = true;
01216                     }
01217                     else
01218                     {
01219                         // Check if there is room for more points
01220                         if (nIntFacePoints >= intersectedFace.size())
01221                         {
01222                             FatalErrorIn("void starMesh::createCoupleMatches()")
01223                                 << setprecision(12)
01224                                 << "void starMesh::createCoupleMatches() : "
01225                                 << endl << "error in intersected face: "
01226                                 << "lost thread for intersection of couple "
01227                                 << coupleI << ". STAR couple ID: "
01228                                 << couples_[coupleI].coupleID() << endl
01229                                 << "Cut Master Face: " << newMasterFace << endl
01230                                 << "points: " << newMasterFace.points(points_)
01231                                 << endl
01232                                 << "Cut Slave Face: " << newSlaveFace << endl
01233                                 << "points: " << newSlaveFace.points(points_)
01234                                 << endl << "intersected face: "
01235                                 << intersectedFace
01236                                 << abort(FatalError);
01237                         }
01238 
01239                         // insert the point
01240                         intersectedFace[nIntFacePoints] = nextEdge.end();
01241                         nIntFacePoints++;
01242 
01243                         // grab the current point and the current edge
01244                         curEdge = nextEdge;
01245 
01246 #                       ifdef DEBUG_RIGHT_HAND_WALK
01247                         Info<< "inserted point " << nextEdge.end() << endl
01248                             << "curEdge: " << curEdge << endl;
01249 #                       endif
01250                     }
01251                 }
01252                 while (!completedFace);
01253 
01254                 // resize the face
01255                 intersectedFace.setSize(nIntFacePoints);
01256 
01257 #               ifdef DEBUG_COUPLE
01258                 Info << "intersectedFace: " << intersectedFace << endl;
01259 #               endif
01260 
01261                 // check the intersection face for duplicate points
01262                 forAll (intersectedFace, checkI)
01263                 {
01264                     for
01265                     (
01266                         label checkJ = checkI + 1;
01267                         checkJ < intersectedFace.size();
01268                         checkJ++
01269                     )
01270                     {
01271                         if (intersectedFace[checkI] == intersectedFace[checkJ])
01272                         {
01273                             FatalErrorIn("void starMesh::createCoupleMatches()")
01274                                 << setprecision(12)
01275                                 << "void starMesh::createCoupleMatches() : "
01276                                 << endl << "error in intersected face: "
01277                                 << "duplicate point in intersected face "
01278                                 << "for couple no " << coupleI
01279                                 << ". STAR couple ID: "
01280                                 << couples_[coupleI].coupleID() << endl
01281                                 << "Duplicate point label: "
01282                                 << intersectedFace[checkI] << endl
01283                                 << "Cut Master Face: " << newMasterFace << endl
01284                                 << "points: " << newMasterFace.points(points_)
01285                                 << endl
01286                                 << "Cut Slave Face: " << newSlaveFace << endl
01287                                 << "points: " << newSlaveFace.points(points_)
01288                                 << endl << "intersected face: "
01289                                 << intersectedFace
01290                                 << abort(FatalError);
01291                         }
01292                     }
01293                 }
01294             }
01295             else
01296             {
01297                 FatalErrorIn("void starMesh::createCoupleMatches()")
01298                     << setprecision(12)
01299                     << "void starMesh::createCoupleMatches() : " << endl
01300                     << "could not find start edge for intersection of couple "
01301                     << coupleI << ". STAR couple ID: "
01302                     << couples_[coupleI].coupleID() << endl
01303                     << "Cut Master Face: " << newMasterFace << endl
01304                     << "points: " << newMasterFace.points(points_) << endl
01305                     << "Cut Slave Face: " << newSlaveFace << endl
01306                     << "points: " << newSlaveFace.points(points_)
01307                     << abort(FatalError);
01308             }
01309 
01310             // Project all points of the intersected face
01311             // onto the master face to ensure closedness
01312             vector pointProjectionNormal = -masterFace.normal(points_);
01313 
01314             forAll (intersectedFace, intPointI)
01315             {
01316 #               ifdef DEBUG_COUPLE_PROJECTION
01317                 Info << "Proj: old point: "
01318                     << points_[intersectedFace[intPointI]] << endl;
01319 #               endif
01320 
01321                 pointHit projHit =
01322                     masterFace.ray
01323                     (
01324                         points_[intersectedFace[intPointI]],
01325                         pointProjectionNormal,
01326                         points_,
01327                         intersection::FULL_RAY
01328                     );
01329 
01330                 if (projHit.hit())
01331                 {
01332                     points_[intersectedFace[intPointI]] =
01333                         projHit.hitPoint();
01334 
01335 #                   ifdef DEBUG_COUPLE_PROJECTION
01336                     Info<< "      new point: "
01337                         << points_[intersectedFace[intPointI]] << endl;
01338 #                   endif
01339                 }
01340             }
01341 
01342             // Check the direction of the intersection face
01343             if
01344             (
01345                 (
01346                     masterFace.normal(points_)
01347                   & intersectedFace.normal(points_)
01348                 ) < VSMALL
01349             )
01350             {
01351                 intersectedFace = intersectedFace.reverseFace();
01352             }
01353 
01354             // Add the new face to both master and slave
01355 
01356             // Master face is replaced by a set of slave faces
01357             Map<SLList<label> >::iterator crfMasterIter =
01358                 cellRemovedFaces.find(fp.masterCell());
01359 
01360             if (crfMasterIter == cellRemovedFaces.end())
01361             {
01362                 cellRemovedFaces.insert
01363                 (
01364                     fp.masterCell(),
01365                     fp.masterFace()
01366                 );
01367             }
01368             else
01369             {
01370                 crfMasterIter().append(fp.masterFace());
01371             }
01372 
01373             Map<SLList<label> >::iterator crfSlaveIter =
01374                 cellRemovedFaces.find(fp.slaveCell());
01375 
01376             if (crfSlaveIter == cellRemovedFaces.end())
01377             {
01378                 cellRemovedFaces.insert
01379                 (
01380                     fp.slaveCell(),
01381                     fp.slaveFace()
01382                 );
01383             }
01384             else
01385             {
01386                 crfSlaveIter().append(fp.slaveFace());
01387             }
01388 
01389             Map<SLList<face> >::iterator cafMasterIter =
01390                 cellAddedFaces.find(fp.masterCell());
01391             if (cafMasterIter == cellAddedFaces.end())
01392             {
01393                 cellAddedFaces.insert
01394                 (
01395                     fp.masterCell(),
01396                     SLList<face>(intersectedFace)
01397                 );
01398             }
01399             else
01400             {
01401                 cafMasterIter().append(intersectedFace);
01402             }
01403 
01404             Map<SLList<face> >::iterator cafSlaveIter =
01405                 cellAddedFaces.find(fp.slaveCell());
01406             if (cafSlaveIter == cellAddedFaces.end())
01407             {
01408                 cellAddedFaces.insert
01409                 (
01410                     fp.slaveCell(),
01411                     SLList<face>(intersectedFace.reverseFace())
01412                 );
01413             }
01414             else
01415             {
01416                 cafSlaveIter().append(intersectedFace.reverseFace());
01417             }
01418         } // end of arbitrary match
01419     }
01420 
01421     if (couples_.size())
01422     {
01423         // Loop through all cells and reset faces for removal to zero size
01424         const labelList crfToc = cellRemovedFaces.toc();
01425 
01426         forAll (crfToc, cellI)
01427         {
01428             const label curCell = crfToc[cellI];
01429 
01430             const SLList<label>& curRemovedFaces = cellRemovedFaces[curCell];
01431 
01432             for
01433             (
01434                 SLList<label>::const_iterator curRemovedFacesIter =
01435                     curRemovedFaces.begin();
01436                 curRemovedFacesIter != curRemovedFaces.end();
01437                 ++curRemovedFacesIter
01438             )
01439             {
01440                 cellFaces_[curCell][curRemovedFacesIter()].setSize(0);
01441             }
01442 
01443             if (curRemovedFaces.size())
01444             {
01445                 // reset the shape pointer to unknown
01446                 cellShapes_[curCell] = cellShape(*unknownPtr_, labelList(0));
01447             }
01448         }
01449 
01450         const labelList cafToc = cellAddedFaces.toc();
01451 
01452         // Insert the new faces into the list
01453         forAll (cafToc, cellI)
01454         {
01455             const label curCell = cafToc[cellI];
01456 
01457             const SLList<face>& curAddedFaces = cellAddedFaces[curCell];
01458 
01459             faceList oldFaces = cellFaces_[curCell];
01460 
01461             faceList& newFaces = cellFaces_[curCell];
01462 
01463             newFaces.setSize(oldFaces.size() + curAddedFaces.size());
01464             label nNewFaces = 0;
01465 
01466             // copy original faces that have not been removed
01467             forAll (oldFaces, faceI)
01468             {
01469                 if (oldFaces[faceI].size())
01470                 {
01471                     newFaces[nNewFaces] = oldFaces[faceI];
01472                     nNewFaces++;
01473                 }
01474             }
01475 
01476             // add new faces
01477             for
01478             (
01479                 SLList<face>::const_iterator curAddedFacesIter =
01480                     curAddedFaces.begin();
01481                 curAddedFacesIter != curAddedFaces.end();
01482                 ++curAddedFacesIter
01483             )
01484             {
01485                 newFaces[nNewFaces] = curAddedFacesIter();
01486                 nNewFaces++;
01487             }
01488 
01489             // reset the size of the face list
01490             newFaces.setSize(nNewFaces);
01491 
01492             if (curAddedFaces.size())
01493             {
01494                 // reset the shape pointer to unknown
01495                 cellShapes_[curCell] = cellShape(*unknownPtr_, labelList(0));
01496             }
01497         }
01498 
01499         // Resize the point list to the number of created points
01500         points_.setSize(nLivePoints);
01501 
01502         // Finished
01503         Info << "Finished doing couples" << endl;
01504     }
01505 }
01506 
01507 
01508 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines