00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
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
00042
00043
00044
00045
00046
00047 const label cellMapSize = min
00048 (
00049 cellShapes_.size()/10,
00050 couples_.size()*2
00051 );
00052
00053
00054 Map<SLList<face> > cellAddedFaces(cellMapSize);
00055
00056 Map<SLList<label> > cellRemovedFaces(cellMapSize);
00057
00058
00059
00060
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
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
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
00114 if (fp.integralMatch())
00115 {
00116
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
00152
00153
00154 SLList<point> coupleFacePoints;
00155
00156
00157 edgeList masterEdges = masterFace.edges();
00158 List<SLList<label> > masterEdgePoints(masterEdges.size());
00159
00160
00161 edgeList slaveEdges = slaveFace.edges();
00162 List<SLList<label> > slaveEdgePoints(slaveEdges.size());
00163
00164
00165 vector n = masterFace.normal(points_);
00166 n /= mag(n) + VSMALL;
00167
00168
00169
00170 forAll (masterEdges, masterEdgeI)
00171 {
00172 const edge& curMasterEdge = masterEdges[masterEdgeI];
00173
00174 point P = points_[curMasterEdge.start()];
00175
00176
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
00186
00187
00188 forAll (slaveEdges, slaveEdgeI)
00189 {
00190 const edge& curSlaveEdge = slaveEdges[slaveEdgeI];
00191
00192 point S = points_[curSlaveEdge.start()];
00193
00194
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
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
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
00231 # ifdef DEBUG_COUPLE_INTERSECTION
00232 Info<< "intersection of non-parallel edges"
00233 << endl;
00234 # endif
00235
00236
00237
00238
00239
00240 if (alpha < smallMergeTol_)
00241 {
00242
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
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
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
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
00322
00323
00324 coupleFacePoints.append
00325 (P + alpha*curMasterEdge.vec(points_));
00326 }
00327 }
00328 }
00329 }
00330 else
00331 {
00332
00333
00334
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344 vector ps = S - P;
00345
00346 bool colinear = false;
00347
00348 if (mag(ps) < SMALL)
00349 {
00350
00351 colinear = true;
00352 }
00353 else if
00354 (
00355 (ps & d)/(mag(ps)*mag(d)) > 1.0 - smallMergeTol_
00356 )
00357 {
00358
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
00407
00408
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 }
00460 }
00461 }
00462 }
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
00474 if (nLivePoints + coupleFacePoints.size() >= points_.size())
00475 {
00476
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
00494
00495
00496
00497
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
00521 tmpMasterFace[nTmpMasterLabels] =
00522 masterEdges[masterEdgeI].start();
00523 nTmpMasterLabels++;
00524
00525
00526 const SLList<label>& curMEdgePoints =
00527 masterEdgePoints[masterEdgeI];
00528
00529
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
00541 edgeVector /= magSqr(edgeVector);
00542
00543 point edgeStartPoint =
00544 points_[masterEdges[masterEdgeI].start()];
00545
00546
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
00606 tmpMasterFace[nTmpMasterLabels] =
00607 nextPointLabel;
00608 nTmpMasterLabels++;
00609 }
00610 else
00611 {
00612 break;
00613 }
00614 }
00615 }
00616
00617
00618 tmpMasterFace.setSize(nTmpMasterLabels);
00619
00620 # ifdef DEBUG_FACE_ORDERING
00621 Info << "tmpMasterFace: " << tmpMasterFace << endl;
00622 # endif
00623
00624
00625 face newMasterFace(labelList(tmpMasterFace.size(), labelMax));
00626
00627
00628
00629
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
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
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
00691
00692
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
00715 tmpSlaveFace[nTmpSlaveLabels] =
00716 slaveEdges[slaveEdgeI].start();
00717 nTmpSlaveLabels++;
00718
00719
00720 const SLList<label>& curSEdgePoints =
00721 slaveEdgePoints[slaveEdgeI];
00722
00723
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
00735 edgeVector /= magSqr(edgeVector);
00736
00737 point edgeStartPoint =
00738 points_[slaveEdges[slaveEdgeI].start()];
00739
00740
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
00800 tmpSlaveFace[nTmpSlaveLabels] =
00801 nextPointLabel;
00802 nTmpSlaveLabels++;
00803 }
00804 else
00805 {
00806 break;
00807 }
00808 }
00809 }
00810
00811
00812 tmpSlaveFace.setSize(nTmpSlaveLabels);
00813
00814 # ifdef DEBUG_FACE_ORDERING
00815 Info << "tmpSlaveFace: " << tmpSlaveFace << endl;
00816 # endif
00817
00818
00819 face newSlaveFace(labelList(tmpSlaveFace.size(), labelMax));
00820
00821
00822
00823
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
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
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
00869
00870
00871
00872
00873
00874
00875
00876
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
00889
00890
00891
00892 label startEdgeFound = 0;
00893
00894 vector masterProjDir = -newMasterFace.normal(points_);
00895
00896 forAll (newSlaveEdges, edgeI)
00897 {
00898
00899
00900
00901
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
00942
00943
00944
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
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
00992
00993
00994
00995
00996
00997
00998
00999
01000
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
01010
01011
01012
01013
01014
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
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
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
01070 }
01071 }
01072
01073
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
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
01118 vector ahead = curEdge.vec(points_);
01119 ahead -= planeNormal*(planeNormal & ahead);
01120 ahead /= mag(ahead) + VSMALL;
01121
01122
01123 vector right = ahead ^ planeNormal;
01124 right /= mag(right) + VSMALL;
01125
01126
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
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)
01162 {
01163 if (curRightTurn < 0)
01164 {
01165
01166 if (curGoStraight > goStraight)
01167 {
01168 # ifdef DEBUG_RIGHT_HAND_WALK
01169 Info << "a" << endl;
01170 # endif
01171
01172
01173 nextEdge = etcIter();
01174 rightTurn = curRightTurn;
01175 goStraight = curGoStraight;
01176 }
01177 }
01178 else
01179 {
01180 # ifdef DEBUG_RIGHT_HAND_WALK
01181 Info << "b" << endl;
01182 # endif
01183
01184
01185 nextEdge = etcIter();
01186 rightTurn = curRightTurn;
01187 goStraight = curGoStraight;
01188 }
01189 }
01190 else
01191 {
01192
01193 if (curRightTurn >= 0)
01194 {
01195
01196 if (curGoStraight < goStraight)
01197 {
01198 # ifdef DEBUG_RIGHT_HAND_WALK
01199 Info << "c" << endl;
01200 # endif
01201
01202
01203 nextEdge = etcIter();
01204 rightTurn = curRightTurn;
01205 goStraight = curGoStraight;
01206 }
01207 }
01208 }
01209 }
01210
01211
01212 if (nextEdge.end() == intersectedFace[0])
01213 {
01214
01215 completedFace = true;
01216 }
01217 else
01218 {
01219
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
01240 intersectedFace[nIntFacePoints] = nextEdge.end();
01241 nIntFacePoints++;
01242
01243
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
01255 intersectedFace.setSize(nIntFacePoints);
01256
01257 # ifdef DEBUG_COUPLE
01258 Info << "intersectedFace: " << intersectedFace << endl;
01259 # endif
01260
01261
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
01311
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
01343 if
01344 (
01345 (
01346 masterFace.normal(points_)
01347 & intersectedFace.normal(points_)
01348 ) < VSMALL
01349 )
01350 {
01351 intersectedFace = intersectedFace.reverseFace();
01352 }
01353
01354
01355
01356
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 }
01419 }
01420
01421 if (couples_.size())
01422 {
01423
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
01446 cellShapes_[curCell] = cellShape(*unknownPtr_, labelList(0));
01447 }
01448 }
01449
01450 const labelList cafToc = cellAddedFaces.toc();
01451
01452
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
01467 forAll (oldFaces, faceI)
01468 {
01469 if (oldFaces[faceI].size())
01470 {
01471 newFaces[nNewFaces] = oldFaces[faceI];
01472 nNewFaces++;
01473 }
01474 }
01475
01476
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
01490 newFaces.setSize(nNewFaces);
01491
01492 if (curAddedFaces.size())
01493 {
01494
01495 cellShapes_[curCell] = cellShape(*unknownPtr_, labelList(0));
01496 }
01497 }
01498
01499
01500 points_.setSize(nLivePoints);
01501
01502
01503 Info << "Finished doing couples" << endl;
01504 }
01505 }
01506
01507
01508