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 #include "boundaryMesh.H"
00027 #include <OpenFOAM/Time.H>
00028 #include <OpenFOAM/polyMesh.H>
00029 #include <dynamicMesh/repatchPolyTopoChanger.H>
00030 #include <OpenFOAM/faceList.H>
00031 #include <meshTools/indexedOctree.H>
00032 #include "treeDataPrimitivePatch.H"
00033 #include <triSurface/triSurface.H>
00034 #include <OpenFOAM/SortableList.H>
00035 #include <OpenFOAM/OFstream.H>
00036 #include "uindirectPrimitivePatch.H"
00037
00038
00039
00040 defineTypeNameAndDebug(Foam::boundaryMesh, 0);
00041
00042
00043 const Foam::vector Foam::boundaryMesh::splitNormal_(3, 2, 1);
00044
00045
00046 const Foam::scalar Foam::boundaryMesh::distanceTol_ = 1E-2;
00047
00048
00049
00050
00051
00052 Foam::label Foam::boundaryMesh::nFeatureEdges(label pointI) const
00053 {
00054 label nFeats = 0;
00055
00056 const labelList& pEdges = mesh().pointEdges()[pointI];
00057
00058 forAll(pEdges, pEdgeI)
00059 {
00060 label edgeI = pEdges[pEdgeI];
00061
00062 if (edgeToFeature_[edgeI] != -1)
00063 {
00064 nFeats++;
00065 }
00066 }
00067 return nFeats;
00068 }
00069
00070
00071
00072 Foam::label Foam::boundaryMesh::nextFeatureEdge
00073 (
00074 const label edgeI,
00075 const label vertI
00076 ) const
00077 {
00078 const labelList& pEdges = mesh().pointEdges()[vertI];
00079
00080 forAll(pEdges, pEdgeI)
00081 {
00082 label nbrEdgeI = pEdges[pEdgeI];
00083
00084 if (nbrEdgeI != edgeI)
00085 {
00086 label featI = edgeToFeature_[nbrEdgeI];
00087
00088 if (featI != -1)
00089 {
00090 return nbrEdgeI;
00091 }
00092 }
00093 }
00094
00095 return -1;
00096 }
00097
00098
00099
00100
00101
00102 Foam::labelList Foam::boundaryMesh::collectSegment
00103 (
00104 const boolList& isFeaturePoint,
00105 const label startEdgeI,
00106 boolList& featVisited
00107 ) const
00108 {
00109
00110
00111 label edgeI = startEdgeI;
00112
00113 const edge& e = mesh().edges()[edgeI];
00114
00115 label vertI = e.start();
00116
00117 while (!isFeaturePoint[vertI])
00118 {
00119
00120
00121 edgeI = nextFeatureEdge(edgeI, vertI);
00122
00123 if ((edgeI == -1) || (edgeI == startEdgeI))
00124 {
00125 break;
00126 }
00127
00128
00129
00130 const edge& e = mesh().edges()[edgeI];
00131
00132 vertI = e.otherVertex(vertI);
00133 }
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144 labelList featLabels(featureEdges_.size());
00145
00146 label featLabelI = 0;
00147
00148 label initEdgeI = edgeI;
00149
00150 do
00151 {
00152
00153 label featI = edgeToFeature_[edgeI];
00154
00155 if (featI == -1)
00156 {
00157 FatalErrorIn("boundaryMesh::collectSegment")
00158 << "Problem" << abort(FatalError);
00159 }
00160 featLabels[featLabelI++] = featI;
00161
00162 featVisited[featI] = true;
00163
00164
00165
00166 const edge& e = mesh().edges()[edgeI];
00167
00168 vertI = e.otherVertex(vertI);
00169
00170
00171
00172 edgeI = nextFeatureEdge(edgeI, vertI);
00173
00174 if ((edgeI == -1) || (edgeI == initEdgeI))
00175 {
00176 break;
00177 }
00178 }
00179 while (!isFeaturePoint[vertI]);
00180
00181
00182
00183 featLabels.setSize(featLabelI);
00184
00185 return featLabels;
00186 }
00187
00188
00189 void Foam::boundaryMesh::markEdges
00190 (
00191 const label maxDistance,
00192 const label edgeI,
00193 const label distance,
00194 labelList& minDistance,
00195 DynamicList<label>& visited
00196 ) const
00197 {
00198 if (distance < maxDistance)
00199 {
00200
00201
00202 if (minDistance[edgeI] == -1)
00203 {
00204
00205 visited.append(edgeI);
00206 }
00207 else if (minDistance[edgeI] <= distance)
00208 {
00209
00210 return;
00211 }
00212
00213 minDistance[edgeI] = distance;
00214
00215 const edge& e = mesh().edges()[edgeI];
00216
00217
00218 const labelList& startEdges = mesh().pointEdges()[e.start()];
00219
00220 forAll(startEdges, pEdgeI)
00221 {
00222 markEdges
00223 (
00224 maxDistance,
00225 startEdges[pEdgeI],
00226 distance+1,
00227 minDistance,
00228 visited
00229 );
00230 }
00231
00232
00233 const labelList& endEdges = mesh().pointEdges()[e.end()];
00234
00235 forAll(endEdges, pEdgeI)
00236 {
00237 markEdges
00238 (
00239 maxDistance,
00240 endEdges[pEdgeI],
00241 distance+1,
00242 minDistance,
00243 visited
00244 );
00245 }
00246 }
00247 }
00248
00249
00250 Foam::label Foam::boundaryMesh::findPatchID
00251 (
00252 const polyPatchList& patches,
00253 const word& patchName
00254 ) const
00255 {
00256 forAll(patches, patchI)
00257 {
00258 if (patches[patchI].name() == patchName)
00259 {
00260 return patchI;
00261 }
00262 }
00263
00264 return -1;
00265 }
00266
00267
00268 Foam::wordList Foam::boundaryMesh::patchNames() const
00269 {
00270 wordList names(patches_.size());
00271
00272 forAll(patches_, patchI)
00273 {
00274 names[patchI] = patches_[patchI].name();
00275 }
00276 return names;
00277 }
00278
00279
00280 Foam::label Foam::boundaryMesh::whichPatch
00281 (
00282 const polyPatchList& patches,
00283 const label faceI
00284 ) const
00285 {
00286 forAll(patches, patchI)
00287 {
00288 const polyPatch& pp = patches[patchI];
00289
00290 if ((faceI >= pp.start()) && (faceI < (pp.start() + pp.size())))
00291 {
00292 return patchI;
00293 }
00294 }
00295 return -1;
00296 }
00297
00298
00299
00300
00301 Foam::labelList Foam::boundaryMesh::faceToEdge
00302 (
00303 const boolList& regionEdge,
00304 const label region,
00305 const labelList& changedFaces,
00306 labelList& edgeRegion
00307 ) const
00308 {
00309 labelList changedEdges(mesh().nEdges(), -1);
00310 label changedI = 0;
00311
00312 forAll(changedFaces, i)
00313 {
00314 label faceI = changedFaces[i];
00315
00316 const labelList& fEdges = mesh().faceEdges()[faceI];
00317
00318 forAll(fEdges, fEdgeI)
00319 {
00320 label edgeI = fEdges[fEdgeI];
00321
00322 if (!regionEdge[edgeI] && (edgeRegion[edgeI] == -1))
00323 {
00324 edgeRegion[edgeI] = region;
00325
00326 changedEdges[changedI++] = edgeI;
00327 }
00328 }
00329 }
00330
00331 changedEdges.setSize(changedI);
00332
00333 return changedEdges;
00334 }
00335
00336
00337
00338 Foam::labelList Foam::boundaryMesh::edgeToFace
00339 (
00340 const label region,
00341 const labelList& changedEdges,
00342 labelList& faceRegion
00343 ) const
00344 {
00345 labelList changedFaces(mesh().size(), -1);
00346 label changedI = 0;
00347
00348 forAll(changedEdges, i)
00349 {
00350 label edgeI = changedEdges[i];
00351
00352 const labelList& eFaces = mesh().edgeFaces()[edgeI];
00353
00354 forAll(eFaces, eFaceI)
00355 {
00356 label faceI = eFaces[eFaceI];
00357
00358 if (faceRegion[faceI] == -1)
00359 {
00360 faceRegion[faceI] = region;
00361
00362 changedFaces[changedI++] = faceI;
00363 }
00364 }
00365 }
00366
00367 changedFaces.setSize(changedI);
00368
00369 return changedFaces;
00370 }
00371
00372
00373
00374 void Foam::boundaryMesh::markZone
00375 (
00376 const boolList& borderEdge,
00377 label faceI,
00378 label currentZone,
00379 labelList& faceZone
00380 ) const
00381 {
00382 faceZone[faceI] = currentZone;
00383
00384
00385 labelList changedFaces(1, faceI);
00386
00387 labelList changedEdges;
00388
00389
00390 labelList edgeZone(mesh().nEdges(), -1);
00391
00392 while(1)
00393 {
00394 changedEdges = faceToEdge
00395 (
00396 borderEdge,
00397 currentZone,
00398 changedFaces,
00399 edgeZone
00400 );
00401
00402 if (debug)
00403 {
00404 Pout<< "From changedFaces:" << changedFaces.size()
00405 << " to changedEdges:" << changedEdges.size()
00406 << endl;
00407 }
00408
00409 if (changedEdges.empty())
00410 {
00411 break;
00412 }
00413
00414 changedFaces = edgeToFace(currentZone, changedEdges, faceZone);
00415
00416 if (debug)
00417 {
00418 Pout<< "From changedEdges:" << changedEdges.size()
00419 << " to changedFaces:" << changedFaces.size()
00420 << endl;
00421 }
00422
00423 if (changedFaces.empty())
00424 {
00425 break;
00426 }
00427 }
00428 }
00429
00430
00431
00432
00433
00434 Foam::boundaryMesh::boundaryMesh()
00435 :
00436 meshPtr_(NULL),
00437 patches_(),
00438 meshFace_(),
00439 featurePoints_(),
00440 featureEdges_(),
00441 featureToEdge_(),
00442 edgeToFeature_(),
00443 featureSegments_(),
00444 extraEdges_()
00445 {}
00446
00447
00448
00449
00450 Foam::boundaryMesh::~boundaryMesh()
00451 {
00452 clearOut();
00453 }
00454
00455
00456 void Foam::boundaryMesh::clearOut()
00457 {
00458 if (meshPtr_)
00459 {
00460 delete meshPtr_;
00461
00462 meshPtr_ = NULL;
00463 }
00464 }
00465
00466
00467
00468
00469 void Foam::boundaryMesh::read(const polyMesh& mesh)
00470 {
00471 patches_.clear();
00472
00473 patches_.setSize(mesh.boundaryMesh().size());
00474
00475
00476 label nBFaces = mesh.nFaces() - mesh.nInternalFaces();
00477
00478 faceList bFaces(nBFaces);
00479
00480 meshFace_.setSize(nBFaces);
00481
00482 label bFaceI = 0;
00483
00484
00485 forAll(mesh.boundaryMesh(), patchI)
00486 {
00487 const polyPatch& pp = mesh.boundaryMesh()[patchI];
00488
00489 patches_.set
00490 (
00491 patchI,
00492 new boundaryPatch
00493 (
00494 pp.name(),
00495 patchI,
00496 pp.size(),
00497 bFaceI,
00498 pp.type()
00499 )
00500 );
00501
00502
00503 forAll(pp, patchFaceI)
00504 {
00505 meshFace_[bFaceI] = pp.start() + patchFaceI;
00506
00507 bFaces[bFaceI] = pp[patchFaceI];
00508
00509 bFaceI++;
00510 }
00511 }
00512
00513
00514 if (debug)
00515 {
00516 Pout<< "read : patches now:" << endl;
00517
00518 forAll(patches_, patchI)
00519 {
00520 const boundaryPatch& bp = patches_[patchI];
00521
00522 Pout<< " name : " << bp.name() << endl
00523 << " size : " << bp.size() << endl
00524 << " start : " << bp.start() << endl
00525 << " type : " << bp.physicalType() << endl
00526 << endl;
00527 }
00528 }
00529
00530
00531
00532
00533
00534
00535 PrimitivePatch<face, List, const pointField&> globalPatch
00536 (
00537 bFaces,
00538 mesh.points()
00539 );
00540
00541
00542 clearOut();
00543
00544 meshPtr_ = new bMesh(globalPatch.localFaces(), globalPatch.localPoints());
00545
00546
00547 if (debug & 2)
00548 {
00549 const bMesh& msh = *meshPtr_;
00550
00551 Pout<< "** Start of Faces **" << endl;
00552
00553 forAll(msh, faceI)
00554 {
00555 const face& f = msh[faceI];
00556
00557 point ctr(vector::zero);
00558
00559 forAll(f, fp)
00560 {
00561 ctr += msh.points()[f[fp]];
00562 }
00563 ctr /= f.size();
00564
00565 Pout<< " " << faceI
00566 << " ctr:" << ctr
00567 << " verts:" << f
00568 << endl;
00569 }
00570
00571 Pout<< "** End of Faces **" << endl;
00572
00573 Pout<< "** Start of Points **" << endl;
00574
00575 forAll(msh.points(), pointI)
00576 {
00577 Pout<< " " << pointI
00578 << " coord:" << msh.points()[pointI]
00579 << endl;
00580 }
00581
00582 Pout<< "** End of Points **" << endl;
00583 }
00584
00585
00586 featurePoints_.setSize(0);
00587 featureEdges_.setSize(0);
00588
00589 featureToEdge_.setSize(0);
00590 edgeToFeature_.setSize(meshPtr_->nEdges());
00591 edgeToFeature_ = -1;
00592
00593 featureSegments_.setSize(0);
00594
00595 extraEdges_.setSize(0);
00596 }
00597
00598
00599 void Foam::boundaryMesh::readTriSurface(const fileName& fName)
00600 {
00601 triSurface surf(fName);
00602
00603 if (surf.empty())
00604 {
00605 return;
00606 }
00607
00608
00609 SortableList<label> regions(surf.size());
00610
00611 forAll(surf, triI)
00612 {
00613 regions[triI] = surf[triI].region();
00614 }
00615 regions.sort();
00616
00617
00618 Map<label> regionToBoundaryPatch;
00619
00620 label oldRegion = -1111;
00621 label boundPatch = 0;
00622
00623 forAll(regions, i)
00624 {
00625 if (regions[i] != oldRegion)
00626 {
00627 regionToBoundaryPatch.insert(regions[i], boundPatch);
00628
00629 oldRegion = regions[i];
00630 boundPatch++;
00631 }
00632 }
00633
00634 const geometricSurfacePatchList& surfPatches = surf.patches();
00635
00636 patches_.clear();
00637
00638 if (surfPatches.size() == regionToBoundaryPatch.size())
00639 {
00640
00641
00642
00643 patches_.setSize(surfPatches.size());
00644
00645
00646 forAll(surfPatches, patchI)
00647 {
00648 const geometricSurfacePatch& surfPatch = surfPatches[patchI];
00649
00650 patches_.set
00651 (
00652 patchI,
00653 new boundaryPatch
00654 (
00655 surfPatch.name(),
00656 patchI,
00657 0,
00658 0,
00659 surfPatch.geometricType()
00660 )
00661 );
00662 }
00663 }
00664 else
00665 {
00666
00667
00668 patches_.setSize(regionToBoundaryPatch.size());
00669
00670 forAll(patches_, patchI)
00671 {
00672 patches_.set
00673 (
00674 patchI,
00675 new boundaryPatch
00676 (
00677 "patch" + name(patchI),
00678 patchI,
00679 0,
00680 0,
00681 "empty"
00682 )
00683 );
00684 }
00685 }
00686
00687
00688
00689
00690
00691 const labelList& indices = regions.indices();
00692
00693 faceList bFaces(surf.size());
00694
00695 meshFace_.setSize(surf.size());
00696
00697 label bFaceI = 0;
00698
00699
00700 label surfRegion = regions[0];
00701 label foamRegion = regionToBoundaryPatch[surfRegion];
00702
00703 Pout<< "Surface region " << surfRegion << " becomes boundary patch "
00704 << foamRegion << " with name " << patches_[foamRegion].name() << endl;
00705
00706
00707
00708 label startFaceI = 0;
00709
00710 forAll(indices, indexI)
00711 {
00712 label triI = indices[indexI];
00713
00714 const labelledTri& tri = surf.localFaces()[triI];
00715
00716 if (tri.region() != surfRegion)
00717 {
00718
00719 boundaryPatch& bp = patches_[foamRegion];
00720
00721 bp.size() = bFaceI - startFaceI;
00722 bp.start() = startFaceI;
00723
00724 surfRegion = tri.region();
00725 foamRegion = regionToBoundaryPatch[surfRegion];
00726
00727 Pout<< "Surface region " << surfRegion << " becomes boundary patch "
00728 << foamRegion << " with name " << patches_[foamRegion].name()
00729 << endl;
00730
00731 startFaceI = bFaceI;
00732 }
00733
00734 meshFace_[bFaceI] = triI;
00735
00736 bFaces[bFaceI++] = face(tri);
00737 }
00738
00739
00740 boundaryPatch& bp = patches_[foamRegion];
00741
00742 bp.size() = bFaceI - startFaceI;
00743 bp.start() = startFaceI;
00744
00745
00746
00747
00748
00749 clearOut();
00750
00751
00752 meshPtr_ = new bMesh(bFaces, surf.localPoints());
00753
00754
00755 featurePoints_.setSize(0);
00756 featureEdges_.setSize(0);
00757
00758 featureToEdge_.setSize(0);
00759 edgeToFeature_.setSize(meshPtr_->nEdges());
00760 edgeToFeature_ = -1;
00761
00762 featureSegments_.setSize(0);
00763
00764 extraEdges_.setSize(0);
00765 }
00766
00767
00768 void Foam::boundaryMesh::writeTriSurface(const fileName& fName) const
00769 {
00770 geometricSurfacePatchList surfPatches(patches_.size());
00771
00772 forAll(patches_, patchI)
00773 {
00774 const boundaryPatch& bp = patches_[patchI];
00775
00776 surfPatches[patchI] =
00777 geometricSurfacePatch
00778 (
00779 bp.physicalType(),
00780 bp.name(),
00781 patchI
00782 );
00783 }
00784
00785
00786
00787
00788
00789
00790 labelList nTris(mesh().size());
00791
00792 label totalNTris = getNTris(0, mesh().size(), nTris);
00793
00794
00795 labelList startTri(mesh().size());
00796
00797 label triI = 0;
00798
00799 forAll(mesh(), faceI)
00800 {
00801 startTri[faceI] = triI;
00802
00803 triI += nTris[faceI];
00804 }
00805
00806
00807 labelList triVerts(3*totalNTris);
00808
00809 triangulate(0, mesh().size(), totalNTris, triVerts);
00810
00811
00812
00813 List<labelledTri> tris(totalNTris);
00814
00815 triI = 0;
00816
00817 forAll(patches_, patchI)
00818 {
00819 const boundaryPatch& bp = patches_[patchI];
00820
00821 forAll(bp, patchFaceI)
00822 {
00823 label faceI = bp.start() + patchFaceI;
00824
00825 label triVertI = 3*startTri[faceI];
00826
00827 for (label faceTriI = 0; faceTriI < nTris[faceI]; faceTriI++)
00828 {
00829 label v0 = triVerts[triVertI++];
00830 label v1 = triVerts[triVertI++];
00831 label v2 = triVerts[triVertI++];
00832
00833 tris[triI++] = labelledTri(v0, v1, v2, patchI);
00834 }
00835 }
00836 }
00837
00838 triSurface surf(tris, surfPatches, mesh().points());
00839
00840 OFstream surfStream(fName);
00841
00842 surf.write(surfStream);
00843 }
00844
00845
00846
00847
00848
00849
00850
00851
00852
00853
00854
00855
00856
00857 Foam::labelList Foam::boundaryMesh::getNearest
00858 (
00859 const primitiveMesh& pMesh,
00860 const vector& searchSpan
00861 ) const
00862 {
00863
00864
00865
00866
00867 DynamicList<label> leftFaces(mesh().size()/2);
00868 DynamicList<label> rightFaces(mesh().size()/2);
00869
00870 forAll(mesh(), bFaceI)
00871 {
00872 scalar sign = mesh().faceNormals()[bFaceI] & splitNormal_;
00873
00874 if (sign > -1E-5)
00875 {
00876 rightFaces.append(bFaceI);
00877 }
00878 if (sign < 1E-5)
00879 {
00880 leftFaces.append(bFaceI);
00881 }
00882 }
00883
00884 leftFaces.shrink();
00885 rightFaces.shrink();
00886
00887 if (debug)
00888 {
00889 Pout<< "getNearest :"
00890 << " rightBin:" << rightFaces.size()
00891 << " leftBin:" << leftFaces.size()
00892 << endl;
00893 }
00894
00895 uindirectPrimitivePatch leftPatch
00896 (
00897 UIndirectList<face>(mesh(), leftFaces),
00898 mesh().points()
00899 );
00900 uindirectPrimitivePatch rightPatch
00901 (
00902 UIndirectList<face>(mesh(), rightFaces),
00903 mesh().points()
00904 );
00905
00906
00907
00908 treeBoundBox overallBb(mesh().localPoints());
00909
00910
00911
00912 scalar tol = 1E-6 * overallBb.avgDim();
00913
00914 point& bbMin = overallBb.min();
00915 bbMin.x() -= tol;
00916 bbMin.y() -= tol;
00917 bbMin.z() -= tol;
00918
00919 point& bbMax = overallBb.max();
00920 bbMax.x() += 2*tol;
00921 bbMax.y() += 2*tol;
00922 bbMax.z() += 2*tol;
00923
00924
00925 indexedOctree
00926 <
00927 treeDataPrimitivePatch<face, UIndirectList, const pointField&>
00928 > leftTree
00929 (
00930 treeDataPrimitivePatch<face, UIndirectList, const pointField&>
00931 (
00932 false,
00933 leftPatch
00934 ),
00935 overallBb,
00936 10,
00937 10,
00938 3.0
00939 );
00940 indexedOctree
00941 <
00942 treeDataPrimitivePatch<face, UIndirectList, const pointField&>
00943 > rightTree
00944 (
00945 treeDataPrimitivePatch<face, UIndirectList, const pointField&>
00946 (
00947 false,
00948 rightPatch
00949 ),
00950 overallBb,
00951 10,
00952 10,
00953 3.0
00954 );
00955
00956 if (debug)
00957 {
00958 Pout<< "getNearest : built trees" << endl;
00959 }
00960
00961
00962 const vectorField& ns = mesh().faceNormals();
00963
00964
00965
00966
00967
00968
00969 labelList nearestBFaceI(pMesh.nFaces() - pMesh.nInternalFaces());
00970
00971 treeBoundBox tightest;
00972
00973 const scalar searchDimSqr = magSqr(searchSpan);
00974
00975 forAll(nearestBFaceI, patchFaceI)
00976 {
00977 label meshFaceI = pMesh.nInternalFaces() + patchFaceI;
00978
00979 const point& ctr = pMesh.faceCentres()[meshFaceI];
00980
00981 if (debug && (patchFaceI % 1000) == 0)
00982 {
00983 Pout<< "getNearest : patchFace:" << patchFaceI
00984 << " meshFaceI:" << meshFaceI << " ctr:" << ctr << endl;
00985 }
00986
00987
00988
00989 vector n = pMesh.faceAreas()[meshFaceI];
00990 scalar area = mag(n);
00991 n /= area;
00992
00993 scalar typDim = -GREAT;
00994 const face& f = pMesh.faces()[meshFaceI];
00995
00996 forAll(f, fp)
00997 {
00998 typDim = max(typDim, mag(pMesh.points()[f[fp]] - ctr));
00999 }
01000
01001
01002 pointIndexHit rightInfo = rightTree.findNearest(ctr, searchDimSqr);
01003
01004
01005
01006 pointIndexHit leftInfo = leftTree.findNearest(ctr, searchDimSqr);
01007
01008 if (rightInfo.hit())
01009 {
01010 if (leftInfo.hit())
01011 {
01012
01013 label rightFaceI = rightFaces[rightInfo.index()];
01014 label leftFaceI = leftFaces[leftInfo.index()];
01015
01016 label rightDist = mag(rightInfo.hitPoint()-ctr);
01017 label leftDist = mag(leftInfo.hitPoint()-ctr);
01018
01019 scalar rightSign = n & ns[rightFaceI];
01020 scalar leftSign = n & ns[leftFaceI];
01021
01022 if
01023 (
01024 (rightSign > 0 && leftSign > 0)
01025 || (rightSign < 0 && leftSign < 0)
01026 )
01027 {
01028
01029 if (rightDist < leftDist)
01030 {
01031 nearestBFaceI[patchFaceI] = rightFaceI;
01032 }
01033 else
01034 {
01035 nearestBFaceI[patchFaceI] = leftFaceI;
01036 }
01037 }
01038 else
01039 {
01040
01041
01042
01043
01044
01045
01046
01047 typDim *= distanceTol_;
01048
01049 if (rightDist < typDim && leftDist < typDim)
01050 {
01051
01052 if (rightSign > 0)
01053 {
01054 nearestBFaceI[patchFaceI] = rightFaceI;
01055 }
01056 else
01057 {
01058 nearestBFaceI[patchFaceI] = leftFaceI;
01059 }
01060 }
01061 else
01062 {
01063
01064 if (rightDist < leftDist)
01065 {
01066 nearestBFaceI[patchFaceI] = rightFaceI;
01067 }
01068 else
01069 {
01070 nearestBFaceI[patchFaceI] = leftFaceI;
01071 }
01072 }
01073 }
01074 }
01075 else
01076 {
01077
01078
01079 label rightFaceI = rightFaces[rightInfo.index()];
01080 nearestBFaceI[patchFaceI] = rightFaceI;
01081 }
01082 }
01083 else
01084 {
01085
01086
01087 if (leftInfo.hit())
01088 {
01089
01090
01091 nearestBFaceI[patchFaceI] = leftFaces[leftInfo.index()];
01092 }
01093 else
01094 {
01095
01096 nearestBFaceI[patchFaceI] = -1;
01097 }
01098 }
01099 }
01100
01101 return nearestBFaceI;
01102 }
01103
01104
01105 void Foam::boundaryMesh::patchify
01106 (
01107 const labelList& nearest,
01108 const polyBoundaryMesh& oldPatches,
01109 polyMesh& newMesh
01110 ) const
01111 {
01112
01113
01114
01115
01116
01117
01118 HashTable<label> nameToIndex(2*patches_.size());
01119
01120 Map<word> indexToName(2*patches_.size());
01121
01122
01123 label nNewPatches = patches_.size();
01124
01125 forAll(oldPatches, oldPatchI)
01126 {
01127 const polyPatch& patch = oldPatches[oldPatchI];
01128
01129 label newPatchI = findPatchID(patch.name());
01130
01131 if (newPatchI != -1)
01132 {
01133 nameToIndex.insert(patch.name(), newPatchI);
01134
01135 indexToName.insert(newPatchI, patch.name());
01136 }
01137 }
01138
01139
01140
01141 forAll(patches_, bPatchI)
01142 {
01143 const boundaryPatch& bp = patches_[bPatchI];
01144
01145 if (!nameToIndex.found(bp.name()))
01146 {
01147 nameToIndex.insert(bp.name(), bPatchI);
01148
01149 indexToName.insert(bPatchI, bp.name());
01150 }
01151 }
01152
01153
01154
01155
01156
01157
01158
01159 List<polyPatch*> newPatchPtrList(nNewPatches);
01160
01161 label meshFaceI = newMesh.nInternalFaces();
01162
01163
01164 label facesToBeDone = newMesh.nFaces() - newMesh.nInternalFaces();
01165
01166 forAll(patches_, bPatchI)
01167 {
01168 const boundaryPatch& bp = patches_[bPatchI];
01169
01170 label newPatchI = nameToIndex[bp.name()];
01171
01172
01173 label oldPatchI = findPatchID(oldPatches, bp.name());
01174
01175 if (oldPatchI == -1)
01176 {
01177
01178 if (debug)
01179 {
01180 Pout<< "patchify : Creating new polyPatch:" << bp.name()
01181 << " type:" << bp.physicalType() << endl;
01182 }
01183
01184 newPatchPtrList[newPatchI] = polyPatch::New
01185 (
01186 bp.physicalType(),
01187 bp.name(),
01188 facesToBeDone,
01189 meshFaceI,
01190 newPatchI,
01191 newMesh.boundaryMesh()
01192 ).ptr();
01193
01194 meshFaceI += facesToBeDone;
01195
01196
01197 facesToBeDone = 0;
01198 }
01199 else
01200 {
01201
01202 const polyPatch& oldPatch = oldPatches[oldPatchI];
01203
01204 if (debug)
01205 {
01206 Pout<< "patchify : Cloning existing polyPatch:"
01207 << oldPatch.name() << endl;
01208 }
01209
01210 newPatchPtrList[newPatchI] = oldPatch.clone
01211 (
01212 newMesh.boundaryMesh(),
01213 newPatchI,
01214 facesToBeDone,
01215 meshFaceI
01216 ).ptr();
01217
01218 meshFaceI += facesToBeDone;
01219
01220
01221 facesToBeDone = 0;
01222 }
01223 }
01224
01225
01226 if (debug)
01227 {
01228 Pout<< "Patchify : new polyPatch list:" << endl;
01229
01230 forAll(newPatchPtrList, patchI)
01231 {
01232 const polyPatch& newPatch = *newPatchPtrList[patchI];
01233
01234 if (debug)
01235 {
01236 Pout<< "polyPatch:" << newPatch.name() << endl
01237 << " type :" << newPatch.typeName << endl
01238 << " size :" << newPatch.size() << endl
01239 << " start:" << newPatch.start() << endl
01240 << " index:" << patchI << endl;
01241 }
01242 }
01243 }
01244
01245
01246 repatchPolyTopoChanger polyMeshRepatcher(newMesh);
01247 polyMeshRepatcher.changePatches(newPatchPtrList);
01248
01249
01250
01251
01252
01253 if (newPatchPtrList.size())
01254 {
01255 List<DynamicList<label> > patchFaces(nNewPatches);
01256
01257
01258 label nAvgFaces =
01259 (newMesh.nFaces() - newMesh.nInternalFaces())
01260 / nNewPatches;
01261
01262 forAll(patchFaces, newPatchI)
01263 {
01264 patchFaces[newPatchI].setCapacity(nAvgFaces);
01265 }
01266
01267
01268
01269
01270
01271
01272 forAll(oldPatches, oldPatchI)
01273 {
01274 const polyPatch& patch = oldPatches[oldPatchI];
01275
01276 forAll(patch, patchFaceI)
01277 {
01278
01279
01280 label meshFaceI = patch.start() + patchFaceI;
01281
01282 label bFaceI = meshFaceI - newMesh.nInternalFaces();
01283
01284 patchFaces[whichPatch(nearest[bFaceI])].append(meshFaceI);
01285 }
01286 }
01287
01288 forAll(patchFaces, newPatchI)
01289 {
01290 patchFaces[newPatchI].shrink();
01291 }
01292
01293
01294
01295
01296
01297 for (label newPatchI = 1; newPatchI < patchFaces.size(); newPatchI++)
01298 {
01299 const labelList& pFaces = patchFaces[newPatchI];
01300
01301 forAll(pFaces, pFaceI)
01302 {
01303 polyMeshRepatcher.changePatchID(pFaces[pFaceI], newPatchI);
01304 }
01305 }
01306
01307 polyMeshRepatcher.repatch();
01308 }
01309 }
01310
01311
01312 void Foam::boundaryMesh::setFeatureEdges(const scalar minCos)
01313 {
01314 edgeToFeature_.setSize(mesh().nEdges());
01315
01316 edgeToFeature_ = -1;
01317
01318
01319
01320
01321 featureToEdge_.setSize(mesh().nEdges());
01322
01323 label featureI = 0;
01324
01325 if (minCos >= 0.9999)
01326 {
01327
01328 forAll(mesh().edges(), edgeI)
01329 {
01330 edgeToFeature_[edgeI] = featureI;
01331 featureToEdge_[featureI++] = edgeI;
01332 }
01333 }
01334 else
01335 {
01336 forAll(mesh().edges(), edgeI)
01337 {
01338 const labelList& eFaces = mesh().edgeFaces()[edgeI];
01339
01340 if (eFaces.size() == 2)
01341 {
01342 label face0I = eFaces[0];
01343
01344 label face1I = eFaces[1];
01345
01348
01349
01350
01351
01352
01353
01354 {
01355 const vector& n0 = mesh().faceNormals()[face0I];
01356
01357 const vector& n1 = mesh().faceNormals()[face1I];
01358
01359 float cosAng = n0 & n1;
01360
01361 if (cosAng < minCos)
01362 {
01363 edgeToFeature_[edgeI] = featureI;
01364 featureToEdge_[featureI++] = edgeI;
01365 }
01366 }
01367 }
01368 else
01369 {
01370
01371 edgeToFeature_[edgeI] = featureI;
01372 featureToEdge_[featureI++] = edgeI;
01373 }
01374 }
01375 }
01376
01377
01378 featureToEdge_.setSize(featureI);
01379
01380
01381
01382
01383
01384 featureEdges_.setSize(featureI);
01385 featurePoints_.setSize(mesh().nPoints());
01386
01387 labelList featToMeshPoint(mesh().nPoints(), -1);
01388
01389 label featPtI = 0;
01390
01391 forAll(featureToEdge_, fEdgeI)
01392 {
01393 label edgeI = featureToEdge_[fEdgeI];
01394
01395 const edge& e = mesh().edges()[edgeI];
01396
01397 label start = featToMeshPoint[e.start()];
01398
01399 if (start == -1)
01400 {
01401 featToMeshPoint[e.start()] = featPtI;
01402
01403 featurePoints_[featPtI] = mesh().points()[e.start()];
01404
01405 start = featPtI;
01406
01407 featPtI++;
01408 }
01409
01410 label end = featToMeshPoint[e.end()];
01411
01412 if (end == -1)
01413 {
01414 featToMeshPoint[e.end()] = featPtI;
01415
01416 featurePoints_[featPtI] = mesh().points()[e.end()];
01417
01418 end = featPtI;
01419
01420 featPtI++;
01421 }
01422
01423
01424 featureEdges_[fEdgeI] = edge(start, end);
01425 }
01426
01427
01428 featurePoints_.setSize(featPtI);
01429
01430
01431
01432
01433
01434
01435
01436
01437
01438 boolList isFeaturePoint(mesh().nPoints(), false);
01439
01440 forAll(featureToEdge_, featI)
01441 {
01442 label edgeI = featureToEdge_[featI];
01443
01444 const edge& e = mesh().edges()[edgeI];
01445
01446 if (nFeatureEdges(e.start()) != 2)
01447 {
01448 isFeaturePoint[e.start()] = true;
01449 }
01450
01451 if (nFeatureEdges(e.end()) != 2)
01452 {
01453 isFeaturePoint[e.end()] = true;
01454 }
01455 }
01456
01457
01458
01459
01460
01461
01462
01463 DynamicList<labelList> segments;
01464
01465
01466 boolList featVisited(featureToEdge_.size(), false);
01467
01468 do
01469 {
01470 label startFeatI = -1;
01471
01472 forAll(featVisited, featI)
01473 {
01474 if (!featVisited[featI])
01475 {
01476 startFeatI = featI;
01477
01478 break;
01479 }
01480 }
01481
01482 if (startFeatI == -1)
01483 {
01484
01485 break;
01486 }
01487
01488 segments.append
01489 (
01490 collectSegment
01491 (
01492 isFeaturePoint,
01493 featureToEdge_[startFeatI],
01494 featVisited
01495 )
01496 );
01497 }
01498 while (true);
01499
01500
01501
01502
01503
01504 featureSegments_.setSize(segments.size());
01505
01506 forAll(featureSegments_, segmentI)
01507 {
01508 featureSegments_[segmentI] = segments[segmentI];
01509 }
01510 }
01511
01512
01513 void Foam::boundaryMesh::setExtraEdges(const label edgeI)
01514 {
01515 labelList minDistance(mesh().nEdges(), -1);
01516
01517
01518 DynamicList<label> visitedEdges;
01519
01520
01521 markEdges(8, edgeI, 0, minDistance, visitedEdges);
01522
01523
01524 extraEdges_.transfer(visitedEdges);
01525 }
01526
01527
01528 Foam::label Foam::boundaryMesh::whichPatch(const label faceI) const
01529 {
01530 forAll(patches_, patchI)
01531 {
01532 const boundaryPatch& bp = patches_[patchI];
01533
01534 if ((faceI >= bp.start()) && (faceI < (bp.start() + bp.size())))
01535 {
01536 return patchI;
01537 }
01538 }
01539
01540 FatalErrorIn("boundaryMesh::whichPatch(const label)")
01541 << "Cannot find face " << faceI << " in list of boundaryPatches "
01542 << patches_
01543 << abort(FatalError);
01544
01545 return -1;
01546 }
01547
01548
01549 Foam::label Foam::boundaryMesh::findPatchID(const word& patchName) const
01550 {
01551 forAll(patches_, patchI)
01552 {
01553 if (patches_[patchI].name() == patchName)
01554 {
01555 return patchI;
01556 }
01557 }
01558
01559 return -1;
01560 }
01561
01562
01563 void Foam::boundaryMesh::addPatch(const word& patchName)
01564 {
01565 patches_.setSize(patches_.size() + 1);
01566
01567
01568
01569 label patchI = patches_.size()-1;
01570
01571 boundaryPatch* bpPtr = new boundaryPatch
01572 (
01573 patchName,
01574 patchI,
01575 0,
01576 mesh().size(),
01577 "empty"
01578 );
01579
01580 patches_.set(patchI, bpPtr);
01581
01582 if (debug)
01583 {
01584 Pout<< "addPatch : patches now:" << endl;
01585
01586 forAll(patches_, patchI)
01587 {
01588 const boundaryPatch& bp = patches_[patchI];
01589
01590 Pout<< " name : " << bp.name() << endl
01591 << " size : " << bp.size() << endl
01592 << " start : " << bp.start() << endl
01593 << " type : " << bp.physicalType() << endl
01594 << endl;
01595 }
01596 }
01597 }
01598
01599
01600 void Foam::boundaryMesh::deletePatch(const word& patchName)
01601 {
01602 label delPatchI = findPatchID(patchName);
01603
01604 if (delPatchI == -1)
01605 {
01606 FatalErrorIn("boundaryMesh::deletePatch(const word&")
01607 << "Can't find patch named " << patchName
01608 << abort(FatalError);
01609 }
01610
01611 if (patches_[delPatchI].size())
01612 {
01613 FatalErrorIn("boundaryMesh::deletePatch(const word&")
01614 << "Trying to delete non-empty patch " << patchName
01615 << endl << "Current size:" << patches_[delPatchI].size()
01616 << abort(FatalError);
01617 }
01618
01619 PtrList<boundaryPatch> newPatches(patches_.size() - 1);
01620
01621 for (label patchI = 0; patchI < delPatchI; patchI++)
01622 {
01623 newPatches.set(patchI, patches_[patchI].clone());
01624 }
01625
01626
01627
01628 for (label patchI = delPatchI + 1; patchI < patches_.size(); patchI++)
01629 {
01630 newPatches.set(patchI - 1, patches_[patchI].clone());
01631 }
01632
01633 patches_.clear();
01634
01635 patches_ = newPatches;
01636
01637 if (debug)
01638 {
01639 Pout<< "deletePatch : patches now:" << endl;
01640
01641 forAll(patches_, patchI)
01642 {
01643 const boundaryPatch& bp = patches_[patchI];
01644
01645 Pout<< " name : " << bp.name() << endl
01646 << " size : " << bp.size() << endl
01647 << " start : " << bp.start() << endl
01648 << " type : " << bp.physicalType() << endl
01649 << endl;
01650 }
01651 }
01652 }
01653
01654
01655 void Foam::boundaryMesh::changePatchType
01656 (
01657 const word& patchName,
01658 const word& patchType
01659 )
01660 {
01661 label changeI = findPatchID(patchName);
01662
01663 if (changeI == -1)
01664 {
01665 FatalErrorIn("boundaryMesh::changePatchType(const word&, const word&)")
01666 << "Can't find patch named " << patchName
01667 << abort(FatalError);
01668 }
01669
01670
01671
01672
01673
01674
01675 PtrList<boundaryPatch> newPatches(patches_.size());
01676
01677 forAll(patches_, patchI)
01678 {
01679 if (patchI == changeI)
01680 {
01681
01682 const boundaryPatch& oldBp = patches_[patchI];
01683
01684 boundaryPatch* bpPtr = new boundaryPatch
01685 (
01686 oldBp.name(),
01687 oldBp.index(),
01688 oldBp.size(),
01689 oldBp.start(),
01690 patchType
01691 );
01692
01693 newPatches.set(patchI, bpPtr);
01694 }
01695 else
01696 {
01697
01698 newPatches.set(patchI, patches_[patchI].clone());
01699 }
01700 }
01701
01702 patches_ = newPatches;
01703 }
01704
01705
01706 void Foam::boundaryMesh::changeFaces
01707 (
01708 const labelList& patchIDs,
01709 labelList& oldToNew
01710 )
01711 {
01712 if (patchIDs.size() != mesh().size())
01713 {
01714 FatalErrorIn("boundaryMesh::changeFaces(const labelList& patchIDs)")
01715 << "List of patchIDs not equal to number of faces." << endl
01716 << "PatchIDs size:" << patchIDs.size()
01717 << " nFaces::" << mesh().size()
01718 << abort(FatalError);
01719 }
01720
01721
01722
01723 labelList nFaces(patches_.size(), 0);
01724
01725 forAll(patchIDs, faceI)
01726 {
01727 label patchID = patchIDs[faceI];
01728
01729 if (patchID < 0 || patchID >= patches_.size())
01730 {
01731 FatalErrorIn("boundaryMesh::changeFaces(const labelList&)")
01732 << "PatchID " << patchID << " out of range"
01733 << abort(FatalError);
01734 }
01735 nFaces[patchID]++;
01736 }
01737
01738
01739
01740
01741 labelList startFace(patches_.size());
01742
01743 startFace[0] = 0;
01744
01745 for (label patchI = 1; patchI < patches_.size(); patchI++)
01746 {
01747 startFace[patchI] = startFace[patchI-1] + nFaces[patchI-1];
01748 }
01749
01750
01751 PtrList<boundaryPatch> newPatches(patches_.size());
01752
01753 forAll(patches_, patchI)
01754 {
01755 const boundaryPatch& bp = patches_[patchI];
01756
01757 newPatches.set
01758 (
01759 patchI,
01760 new boundaryPatch
01761 (
01762 bp.name(),
01763 patchI,
01764 nFaces[patchI],
01765 startFace[patchI],
01766 bp.physicalType()
01767 )
01768 );
01769 }
01770 patches_ = newPatches;
01771
01772 if (debug)
01773 {
01774 Pout<< "changeFaces : patches now:" << endl;
01775
01776 forAll(patches_, patchI)
01777 {
01778 const boundaryPatch& bp = patches_[patchI];
01779
01780 Pout<< " name : " << bp.name() << endl
01781 << " size : " << bp.size() << endl
01782 << " start : " << bp.start() << endl
01783 << " type : " << bp.physicalType() << endl
01784 << endl;
01785 }
01786 }
01787
01788
01789
01790 oldToNew.setSize(patchIDs.size());
01791
01792 forAll(patchIDs, faceI)
01793 {
01794 int patchID = patchIDs[faceI];
01795
01796 oldToNew[faceI] = startFace[patchID]++;
01797 }
01798
01799
01800 faceList newFaces(mesh().size());
01801
01802 labelList newMeshFace(mesh().size());
01803
01804 forAll(oldToNew, faceI)
01805 {
01806 newFaces[oldToNew[faceI]] = mesh()[faceI];
01807 newMeshFace[oldToNew[faceI]] = meshFace_[faceI];
01808 }
01809
01810
01811 bMesh* newMeshPtr_ = new bMesh(newFaces, mesh().points());
01812
01813
01814 meshFace_.transfer(newMeshFace);
01815
01816
01817
01818 clearOut();
01819
01820
01821 meshPtr_ = newMeshPtr_;
01822 }
01823
01824
01825 Foam::label Foam::boundaryMesh::getNTris(const label faceI) const
01826 {
01827 const face& f = mesh()[faceI];
01828
01829 return f.nTriangles(mesh().points());
01830 }
01831
01832
01833 Foam::label Foam::boundaryMesh::getNTris
01834 (
01835 const label startFaceI,
01836 const label nFaces,
01837 labelList& nTris
01838 ) const
01839 {
01840 label totalNTris = 0;
01841
01842 nTris.setSize(nFaces);
01843
01844 for (label i = 0; i < nFaces; i++)
01845 {
01846 label faceNTris = getNTris(startFaceI + i);
01847
01848 nTris[i] = faceNTris;
01849
01850 totalNTris += faceNTris;
01851 }
01852 return totalNTris;
01853 }
01854
01855
01856
01857
01858 void Foam::boundaryMesh::triangulate
01859 (
01860 const label startFaceI,
01861 const label nFaces,
01862 const label totalNTris,
01863 labelList& triVerts
01864 ) const
01865 {
01866
01867 triVerts.setSize(3*totalNTris);
01868
01869 label vertI = 0;
01870
01871 for (label i = 0; i < nFaces; i++)
01872 {
01873 label faceI = startFaceI + i;
01874
01875 const face& f = mesh()[faceI];
01876
01877
01878 faceList triFaces(f.nTriangles(mesh().points()));
01879
01880 label nTri = 0;
01881
01882 f.triangles(mesh().points(), nTri, triFaces);
01883
01884
01885
01886 forAll(triFaces, triFaceI)
01887 {
01888 const face& triF = triFaces[triFaceI];
01889
01890 triVerts[vertI++] = triF[0];
01891 triVerts[vertI++] = triF[1];
01892 triVerts[vertI++] = triF[2];
01893 }
01894 }
01895 }
01896
01897
01898
01899 Foam::label Foam::boundaryMesh::getNPoints
01900 (
01901 const label startFaceI,
01902 const label nFaces
01903 ) const
01904 {
01905 SubList<face> patchFaces(mesh(), nFaces, startFaceI);
01906
01907 primitivePatch patch(patchFaces, mesh().points());
01908
01909 return patch.nPoints();
01910 }
01911
01912
01913
01914 void Foam::boundaryMesh::triangulateLocal
01915 (
01916 const label startFaceI,
01917 const label nFaces,
01918 const label totalNTris,
01919 labelList& triVerts,
01920 labelList& localToGlobal
01921 ) const
01922 {
01923 SubList<face> patchFaces(mesh(), nFaces, startFaceI);
01924
01925 primitivePatch patch(patchFaces, mesh().points());
01926
01927
01928 localToGlobal = patch.meshPoints();
01929
01930
01931 triVerts.setSize(3*totalNTris);
01932
01933 label vertI = 0;
01934
01935 for (label i = 0; i < nFaces; i++)
01936 {
01937
01938 const face& f = patch.localFaces()[i];
01939
01940
01941 faceList triFaces(f.nTriangles(patch.localPoints()));
01942
01943 label nTri = 0;
01944
01945 f.triangles(patch.localPoints(), nTri, triFaces);
01946
01947
01948
01949 forAll(triFaces, triFaceI)
01950 {
01951 const face& triF = triFaces[triFaceI];
01952
01953 triVerts[vertI++] = triF[0];
01954 triVerts[vertI++] = triF[1];
01955 triVerts[vertI++] = triF[2];
01956 }
01957 }
01958 }
01959
01960
01961 void Foam::boundaryMesh::markFaces
01962 (
01963 const labelList& protectedEdges,
01964 const label seedFaceI,
01965 boolList& visited
01966 ) const
01967 {
01968 boolList protectedEdge(mesh().nEdges(), false);
01969
01970 forAll(protectedEdges, i)
01971 {
01972 protectedEdge[protectedEdges[i]] = true;
01973 }
01974
01975
01976
01977 labelList currentZone(mesh().size(), -1);
01978
01979
01980 markZone(protectedEdge, seedFaceI, 0, currentZone);
01981
01982
01983 visited.setSize(mesh().size());
01984
01985 forAll(currentZone, faceI)
01986 {
01987 if (currentZone[faceI] == 0)
01988 {
01989 visited[faceI] = true;
01990 }
01991 else
01992 {
01993 visited[faceI] = false;
01994 }
01995 }
01996 }
01997
01998
01999