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 
00030 
00031 
00032 
00033 
00034 
00035 
00036 
00037 
00038 
00039 
00040 
00041 
00042 
00043 
00044 
00045 
00046 
00047 
00048 
00049 
00050 
00051 
00052 
00053 
00054 
00055 
00056 
00057 
00058 
00059 
00060 
00061 
00062 
00063 
00064 
00065 
00066 
00067 
00068 
00069 
00070 
00071 
00072 
00073 
00074 
00075 
00076 #include <OpenFOAM/argList.H>
00077 #include <OpenFOAM/Time.H>
00078 #include <OpenFOAM/timeSelector.H>
00079 #include <finiteVolume/fvMesh.H>
00080 #include <OpenFOAM/mathematicalConstants.H>
00081 #include <dynamicMesh/polyTopoChange.H>
00082 #include <OpenFOAM/mapPolyMesh.H>
00083 #include <OpenFOAM/PackedBoolList.H>
00084 #include <meshTools/meshTools.H>
00085 #include <OpenFOAM/OFstream.H>
00086 #include "meshDualiser.H"
00087 #include <OpenFOAM/ReadFields.H>
00088 #include <finiteVolume/volFields.H>
00089 #include <finiteVolume/surfaceFields.H>
00090 
00091 using namespace Foam;
00092 
00093 
00094 
00095 
00096 
00097 
00098 void simpleMarkFeatures
00099 (
00100     const polyMesh& mesh,
00101     const PackedBoolList& isBoundaryEdge,
00102     const scalar featureAngle,
00103     const bool concaveMultiCells,
00104     const bool doNotPreserveFaceZones,
00105 
00106     labelList& featureFaces,
00107     labelList& featureEdges,
00108     labelList& singleCellFeaturePoints,
00109     labelList& multiCellFeaturePoints
00110 )
00111 {
00112     scalar minCos = Foam::cos(featureAngle * mathematicalConstant::pi/180.0);
00113 
00114     const polyBoundaryMesh& patches = mesh.boundaryMesh();
00115 
00116     
00117     labelHashSet featureEdgeSet;
00118     labelHashSet singleCellFeaturePointSet;
00119     labelHashSet multiCellFeaturePointSet;
00120 
00121 
00122     
00123     
00124 
00125     forAll(patches, patchI)
00126     {
00127         const polyPatch& pp = patches[patchI];
00128         const labelList& meshEdges = pp.meshEdges();
00129 
00130         
00131         for (label edgeI = pp.nInternalEdges(); edgeI < pp.nEdges(); edgeI++)
00132         {
00133             label meshEdgeI = meshEdges[edgeI];
00134             featureEdgeSet.insert(meshEdgeI);
00135             singleCellFeaturePointSet.insert(mesh.edges()[meshEdgeI][0]);
00136             singleCellFeaturePointSet.insert(mesh.edges()[meshEdgeI][1]);
00137         }
00138     }
00139 
00140 
00141 
00142     
00143     
00144     
00145     
00146     
00147 
00148     
00149     primitivePatch allBoundary
00150     (
00151         SubList<face>
00152         (
00153             mesh.faces(),
00154             mesh.nFaces()-mesh.nInternalFaces(),
00155             mesh.nInternalFaces()
00156         ),
00157         mesh.points()
00158     );
00159 
00160     
00161     allBoundary.checkPointManifold(false, &singleCellFeaturePointSet);
00162 
00163     
00164     const labelListList& edgeFaces = allBoundary.edgeFaces();
00165     const labelList& meshPoints = allBoundary.meshPoints();
00166 
00167     forAll(edgeFaces, edgeI)
00168     {
00169         const labelList& eFaces = edgeFaces[edgeI];
00170 
00171         if (eFaces.size() > 2)
00172         {
00173             const edge& e = allBoundary.edges()[edgeI];
00174 
00175             
00176             
00177             
00178             
00179 
00180             singleCellFeaturePointSet.insert(meshPoints[e[0]]);
00181             singleCellFeaturePointSet.insert(meshPoints[e[1]]);
00182         }
00183     }
00184 
00185     
00186     forAll(edgeFaces, edgeI)
00187     {
00188         const labelList& eFaces = edgeFaces[edgeI];
00189 
00190         if (eFaces.size() == 2)
00191         {
00192             label f0 = eFaces[0];
00193             label f1 = eFaces[1];
00194 
00195             
00196             const vector& n0 = allBoundary.faceNormals()[f0];
00197             const vector& n1 = allBoundary.faceNormals()[f1];
00198 
00199             if ((n0 & n1) < minCos)
00200             {
00201                 const edge& e = allBoundary.edges()[edgeI];
00202                 label v0 = meshPoints[e[0]];
00203                 label v1 = meshPoints[e[1]];
00204 
00205                 label meshEdgeI = meshTools::findEdge(mesh, v0, v1);
00206                 featureEdgeSet.insert(meshEdgeI);
00207 
00208                 
00209                 
00210                 vector c1c0
00211                 (
00212                     allBoundary[f1].centre(allBoundary.points())
00213                   - allBoundary[f0].centre(allBoundary.points())
00214                 );
00215 
00216                 if (concaveMultiCells && (c1c0 & n0) > SMALL)
00217                 {
00218                     
00219                     Info<< "Detected concave feature edge:" << edgeI
00220                         << " cos:" << (c1c0 & n0)
00221                         << " coords:"
00222                         << allBoundary.points()[v0]
00223                         << allBoundary.points()[v1]
00224                         << endl;
00225 
00226                     singleCellFeaturePointSet.erase(v0);
00227                     multiCellFeaturePointSet.insert(v0);
00228                     singleCellFeaturePointSet.erase(v1);
00229                     multiCellFeaturePointSet.insert(v1);
00230                 }
00231                 else
00232                 {
00233                     
00234                     if (!multiCellFeaturePointSet.found(v0))
00235                     {
00236                         singleCellFeaturePointSet.insert(v0);
00237                     }
00238                     if (!multiCellFeaturePointSet.found(v1))
00239                     {
00240                         singleCellFeaturePointSet.insert(v1);
00241                     }
00242                 }
00243             }
00244         }
00245     }
00246 
00247 
00248     
00249     
00250 
00251     
00252     labelHashSet featureFaceSet(mesh.nFaces()-mesh.nInternalFaces());
00253     
00254     for (label faceI = mesh.nInternalFaces(); faceI < mesh.nFaces(); faceI++)
00255     {
00256         featureFaceSet.insert(faceI);
00257     }
00258 
00259     
00260     const faceZoneMesh& faceZones = mesh.faceZones();
00261 
00262     if (doNotPreserveFaceZones)
00263     {
00264         if (faceZones.size() > 0)
00265         {
00266             WarningIn("simpleMarkFeatures(..)")
00267                 << "Detected " << faceZones.size()
00268                 << " faceZones. These will not be preserved."
00269                 << endl;
00270         }
00271     }
00272     else
00273     {
00274         if (faceZones.size() > 0)
00275         {
00276             Info<< "Detected " << faceZones.size()
00277                 << " faceZones. Preserving these by marking their"
00278                 << " points, edges and faces as features." << endl;
00279         }
00280 
00281         forAll(faceZones, zoneI)
00282         {
00283             const faceZone& fz = faceZones[zoneI];
00284 
00285             Info<< "Inserting all faces in faceZone " << fz.name()
00286                 << " as features." << endl;
00287 
00288             forAll(fz, i)
00289             {
00290                 label faceI = fz[i];
00291                 const face& f = mesh.faces()[faceI];
00292                 const labelList& fEdges = mesh.faceEdges()[faceI];
00293 
00294                 featureFaceSet.insert(faceI);
00295                 forAll(f, fp)
00296                 {
00297                     
00298                     
00299                     singleCellFeaturePointSet.erase(f[fp]);
00300                     multiCellFeaturePointSet.insert(f[fp]);
00301 
00302                     
00303                     featureEdgeSet.insert(fEdges[fp]);
00304                 }
00305             }
00306         }
00307     }
00308 
00309     
00310     featureFaces = featureFaceSet.toc();
00311     featureEdges = featureEdgeSet.toc();
00312     singleCellFeaturePoints = singleCellFeaturePointSet.toc();
00313     multiCellFeaturePoints = multiCellFeaturePointSet.toc();
00314 }
00315 
00316 
00317 
00318 void dumpFeatures
00319 (
00320     const polyMesh& mesh,
00321     const labelList& featureFaces,
00322     const labelList& featureEdges,
00323     const labelList& singleCellFeaturePoints,
00324     const labelList& multiCellFeaturePoints
00325 )
00326 {
00327     {
00328         OFstream str("featureFaces.obj");
00329         Info<< "Dumping centres of featureFaces to obj file " << str.name()
00330             << endl;
00331         forAll(featureFaces, i)
00332         {
00333             meshTools::writeOBJ(str, mesh.faceCentres()[featureFaces[i]]);
00334         }
00335     }
00336     {
00337         OFstream str("featureEdges.obj");
00338         Info<< "Dumping featureEdges to obj file " << str.name() << endl;
00339         label vertI = 0;
00340 
00341         forAll(featureEdges, i)
00342         {
00343             const edge& e = mesh.edges()[featureEdges[i]];
00344             meshTools::writeOBJ(str, mesh.points()[e[0]]);
00345             vertI++;
00346             meshTools::writeOBJ(str, mesh.points()[e[1]]);
00347             vertI++;
00348             str<< "l " << vertI-1 << ' ' << vertI << nl;
00349         }
00350     }
00351     {
00352         OFstream str("singleCellFeaturePoints.obj");
00353         Info<< "Dumping featurePoints that become a single cell to obj file "
00354             << str.name() << endl;
00355         forAll(singleCellFeaturePoints, i)
00356         {
00357             meshTools::writeOBJ(str, mesh.points()[singleCellFeaturePoints[i]]);
00358         }
00359     }
00360     {
00361         OFstream str("multiCellFeaturePoints.obj");
00362         Info<< "Dumping featurePoints that become multiple cells to obj file "
00363             << str.name() << endl;
00364         forAll(multiCellFeaturePoints, i)
00365         {
00366             meshTools::writeOBJ(str, mesh.points()[multiCellFeaturePoints[i]]);
00367         }
00368     }
00369 }
00370 
00371 
00372 int main(int argc, char *argv[])
00373 {
00374     argList::noParallel();
00375     timeSelector::addOptions(true, false);
00376 
00377     argList::validArgs.append("feature angle[0-180]");
00378     argList::validOptions.insert("splitAllFaces", "");
00379     argList::validOptions.insert("concaveMultiCells", "");
00380     argList::validOptions.insert("doNotPreserveFaceZones", "");
00381     argList::validOptions.insert("overwrite", "");
00382 
00383 #   include <OpenFOAM/setRootCase.H>
00384 #   include <OpenFOAM/createTime.H>
00385 
00386     instantList timeDirs = timeSelector::select0(runTime, args);
00387 
00388 #   include <OpenFOAM/createMesh.H>
00389 
00390     const word oldInstance = mesh.pointsInstance();
00391 
00392     
00393     
00394     
00395     PackedBoolList isBoundaryEdge(mesh.nEdges());
00396     for (label faceI = mesh.nInternalFaces(); faceI < mesh.nFaces(); faceI++)
00397     {
00398         const labelList& fEdges = mesh.faceEdges()[faceI];
00399 
00400         forAll(fEdges, i)
00401         {
00402             isBoundaryEdge.set(fEdges[i], 1);
00403         }
00404     }
00405 
00406     scalar featureAngle(readScalar(IStringStream(args.additionalArgs()[0])()));
00407 
00408     scalar minCos = Foam::cos(featureAngle * mathematicalConstant::pi/180.0);
00409 
00410     Info<< "Feature:" << featureAngle << endl
00411         << "minCos :" << minCos << endl
00412         << endl;
00413 
00414 
00415     const bool splitAllFaces = args.optionFound("splitAllFaces");
00416     if (splitAllFaces)
00417     {
00418         Info<< "Splitting all internal faces to create multiple faces"
00419             << " between two cells." << nl
00420             << endl;
00421     }
00422 
00423     const bool overwrite = args.optionFound("overwrite");
00424     const bool doNotPreserveFaceZones = args.optionFound
00425     (
00426         "doNotPreserveFaceZones"
00427     );
00428     const bool concaveMultiCells = args.optionFound("concaveMultiCells");
00429     if (concaveMultiCells)
00430     {
00431         Info<< "Generating multiple cells for points on concave feature edges."
00432             << nl << endl;
00433     }
00434 
00435 
00436     
00437     labelList featureFaces;
00438     
00439     labelList featureEdges;
00440     
00441     labelList singleCellFeaturePoints;
00442     
00443     labelList multiCellFeaturePoints;
00444 
00445     
00446     simpleMarkFeatures
00447     (
00448         mesh,
00449         isBoundaryEdge,
00450         featureAngle,
00451         concaveMultiCells,
00452         doNotPreserveFaceZones,
00453 
00454         featureFaces,
00455         featureEdges,
00456         singleCellFeaturePoints,
00457         multiCellFeaturePoints
00458     );
00459 
00460     
00461     
00462     
00463     
00464     if (splitAllFaces)
00465     {
00466         featureEdges = identity(mesh.nEdges());
00467         featureFaces = identity(mesh.nFaces());
00468     }
00469 
00470     
00471     dumpFeatures
00472     (
00473         mesh,
00474         featureFaces,
00475         featureEdges,
00476         singleCellFeaturePoints,
00477         multiCellFeaturePoints
00478     );
00479 
00480 
00481 
00482     
00483     IOobjectList objects(mesh, runTime.timeName());
00484 
00485     
00486     PtrList<volScalarField> vsFlds;
00487     ReadFields(mesh, objects, vsFlds);
00488 
00489     PtrList<volVectorField> vvFlds;
00490     ReadFields(mesh, objects, vvFlds);
00491 
00492     PtrList<volSphericalTensorField> vstFlds;
00493     ReadFields(mesh, objects, vstFlds);
00494 
00495     PtrList<volSymmTensorField> vsymtFlds;
00496     ReadFields(mesh, objects, vsymtFlds);
00497 
00498     PtrList<volTensorField> vtFlds;
00499     ReadFields(mesh, objects, vtFlds);
00500 
00501     
00502     PtrList<surfaceScalarField> ssFlds;
00503     ReadFields(mesh, objects, ssFlds);
00504 
00505     PtrList<surfaceVectorField> svFlds;
00506     ReadFields(mesh, objects, svFlds);
00507 
00508     PtrList<surfaceSphericalTensorField> sstFlds;
00509     ReadFields(mesh, objects, sstFlds);
00510 
00511     PtrList<surfaceSymmTensorField> ssymtFlds;
00512     ReadFields(mesh, objects, ssymtFlds);
00513 
00514     PtrList<surfaceTensorField> stFlds;
00515     ReadFields(mesh, objects, stFlds);
00516 
00517 
00518     
00519     polyTopoChange meshMod(mesh.boundaryMesh().size());
00520 
00521     
00522     meshDualiser dualMaker(mesh);
00523 
00524     
00525     
00526     dualMaker.setRefinement
00527     (
00528         splitAllFaces,
00529         featureFaces,
00530         featureEdges,
00531         singleCellFeaturePoints,
00532         multiCellFeaturePoints,
00533         meshMod
00534     );
00535 
00536     
00537     autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
00538 
00539     
00540     mesh.updateMesh(map);
00541 
00542     
00543     if (map().hasMotionPoints())
00544     {
00545         mesh.movePoints(map().preMotionPoints());
00546     }
00547 
00548     if (!overwrite)
00549     {
00550         runTime++;
00551     }
00552     else
00553     {
00554         mesh.setInstance(oldInstance);
00555     }
00556 
00557     Info<< "Writing dual mesh to " << runTime.timeName() << endl;
00558 
00559     mesh.write();
00560 
00561     Info<< "End\n" << endl;
00562 
00563     return 0;
00564 }
00565 
00566 
00567