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
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091 #include <OpenFOAM/argList.H>
00092 #include <OpenFOAM/timeSelector.H>
00093 #include <OpenFOAM/Time.H>
00094 #include <OpenFOAM/polyMesh.H>
00095 #include <OpenFOAM/OFstream.H>
00096 #include <meshTools/meshTools.H>
00097 #include <meshTools/cellSet.H>
00098 #include <meshTools/faceSet.H>
00099
00100 using namespace Foam;
00101
00102
00103
00104 void writeOBJ(const point& pt, Ostream& os)
00105 {
00106 os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
00107 }
00108
00109
00110 void writePoints(const polyMesh& mesh, const fileName& timeName)
00111 {
00112 label vertI = 0;
00113
00114 fileName pointFile(mesh.time().path()/"meshPoints_" + timeName + ".obj");
00115
00116 Info << "Writing mesh points and edges to " << pointFile << endl;
00117
00118 OFstream pointStream(pointFile);
00119
00120 forAll(mesh.points(), pointI)
00121 {
00122 writeOBJ(mesh.points()[pointI], pointStream);
00123 vertI++;
00124 }
00125
00126 forAll(mesh.edges(), edgeI)
00127 {
00128 const edge& e = mesh.edges()[edgeI];
00129
00130 pointStream << "l " << e.start() + 1 << ' ' << e.end() + 1 << nl;
00131 }
00132 }
00133
00134
00135
00136 void writePoints
00137 (
00138 const polyMesh& mesh,
00139 const labelList& cellLabels,
00140 const fileName& timeName
00141 )
00142 {
00143 fileName fName(mesh.time().path()/"meshPoints_" + timeName + ".obj");
00144
00145 Info << "Writing mesh points and edges to " << fName << endl;
00146
00147 OFstream str(fName);
00148
00149
00150 label vertI = 0;
00151
00152
00153 Map<label> pointToObj(6*cellLabels.size());
00154
00155 forAll(cellLabels, i)
00156 {
00157 label cellI = cellLabels[i];
00158
00159 const labelList& cEdges = mesh.cellEdges()[cellI];
00160
00161 forAll(cEdges, cEdgeI)
00162 {
00163 const edge& e = mesh.edges()[cEdges[cEdgeI]];
00164
00165 label v0;
00166
00167 Map<label>::iterator e0Fnd = pointToObj.find(e[0]);
00168
00169 if (e0Fnd == pointToObj.end())
00170 {
00171 meshTools::writeOBJ(str, mesh.points()[e[0]]);
00172 v0 = vertI++;
00173 pointToObj.insert(e[0], v0);
00174 }
00175 else
00176 {
00177 v0 = e0Fnd();
00178 }
00179
00180 label v1;
00181
00182 Map<label>::iterator e1Fnd = pointToObj.find(e[1]);
00183
00184 if (e1Fnd == pointToObj.end())
00185 {
00186 meshTools::writeOBJ(str, mesh.points()[e[1]]);
00187 v1 = vertI++;
00188 pointToObj.insert(e[1], v1);
00189 }
00190 else
00191 {
00192 v1 = e1Fnd();
00193 }
00194
00195
00196 str << "l " << v0+1 << ' ' << v1+1 << nl;
00197 }
00198 }
00199 }
00200
00201
00202
00203 void writePoints
00204 (
00205 const polyMesh& mesh,
00206 const label cellI,
00207 const fileName& timeName
00208 )
00209 {
00210 fileName fName
00211 (
00212 mesh.time().path()
00213 / "meshPoints_" + timeName + '_' + name(cellI) + ".obj"
00214 );
00215
00216 Info << "Writing mesh points and edges to " << fName << endl;
00217
00218 OFstream pointStream(fName);
00219
00220 const cell& cFaces = mesh.cells()[cellI];
00221
00222 meshTools::writeOBJ(pointStream, mesh.faces(), mesh.points(), cFaces);
00223 }
00224
00225
00226
00227
00228 void writeFaceCentres(const polyMesh& mesh,const fileName& timeName)
00229 {
00230 fileName faceFile
00231 (
00232 mesh.time().path()
00233 / "meshFaceCentres_" + timeName + ".obj"
00234 );
00235
00236 Info << "Writing mesh face centres to " << faceFile << endl;
00237
00238 OFstream faceStream(faceFile);
00239
00240 forAll(mesh.faceCentres(), faceI)
00241 {
00242 writeOBJ(mesh.faceCentres()[faceI], faceStream);
00243 }
00244 }
00245
00246
00247 void writeCellCentres(const polyMesh& mesh, const fileName& timeName)
00248 {
00249 fileName cellFile
00250 (
00251 mesh.time().path()/"meshCellCentres_" + timeName + ".obj"
00252 );
00253
00254 Info << "Writing mesh cell centres to " << cellFile << endl;
00255
00256 OFstream cellStream(cellFile);
00257
00258 forAll(mesh.cellCentres(), cellI)
00259 {
00260 writeOBJ(mesh.cellCentres()[cellI], cellStream);
00261 }
00262 }
00263
00264
00265 void writePatchCentres
00266 (
00267 const polyMesh& mesh,
00268 const fileName& timeName
00269 )
00270 {
00271 const polyBoundaryMesh& patches = mesh.boundaryMesh();
00272
00273 forAll(patches, patchI)
00274 {
00275 const polyPatch& pp = patches[patchI];
00276
00277 fileName faceFile
00278 (
00279 mesh.time().path()/"patch_" + pp.name() + '_' + timeName + ".obj"
00280 );
00281
00282 Info << "Writing patch face centres to " << faceFile << endl;
00283
00284 OFstream patchFaceStream(faceFile);
00285
00286 forAll(pp.faceCentres(), faceI)
00287 {
00288 writeOBJ(pp.faceCentres()[faceI], patchFaceStream);
00289 }
00290 }
00291 }
00292
00293
00294 void writePatchFaces
00295 (
00296 const polyMesh& mesh,
00297 const fileName& timeName
00298 )
00299 {
00300 const polyBoundaryMesh& patches = mesh.boundaryMesh();
00301
00302 forAll(patches, patchI)
00303 {
00304 const polyPatch& pp = patches[patchI];
00305
00306 fileName faceFile
00307 (
00308 mesh.time().path()
00309 / "patchFaces_" + pp.name() + '_' + timeName + ".obj"
00310 );
00311
00312 Info << "Writing patch faces to " << faceFile << endl;
00313
00314 OFstream patchFaceStream(faceFile);
00315
00316 forAll(pp.localPoints(), pointI)
00317 {
00318 writeOBJ(pp.localPoints()[pointI], patchFaceStream);
00319 }
00320
00321 forAll(pp.localFaces(), faceI)
00322 {
00323 const face& f = pp.localFaces()[faceI];
00324
00325 patchFaceStream<< 'f';
00326
00327 forAll(f, fp)
00328 {
00329 patchFaceStream << ' ' << f[fp]+1;
00330 }
00331 patchFaceStream << nl;
00332 }
00333 }
00334 }
00335
00336
00337 void writePatchBoundaryEdges
00338 (
00339 const polyMesh& mesh,
00340 const fileName& timeName
00341 )
00342 {
00343 const polyBoundaryMesh& patches = mesh.boundaryMesh();
00344
00345 forAll(patches, patchI)
00346 {
00347 const polyPatch& pp = patches[patchI];
00348
00349 fileName edgeFile
00350 (
00351 mesh.time().path()
00352 / "patchEdges_" + pp.name() + '_' + timeName + ".obj"
00353 );
00354
00355 Info << "Writing patch edges to " << edgeFile << endl;
00356
00357 OFstream patchEdgeStream(edgeFile);
00358
00359 forAll(pp.localPoints(), pointI)
00360 {
00361 writeOBJ(pp.localPoints()[pointI], patchEdgeStream);
00362 }
00363
00364 for (label edgeI = pp.nInternalEdges(); edgeI < pp.nEdges(); edgeI++)
00365 {
00366 if (pp.edgeFaces()[edgeI].size() == 1)
00367 {
00368 const edge& e = pp.edges()[edgeI];
00369
00370 patchEdgeStream<< "l " << e[0]+1 << ' ' << e[1]+1 << nl;
00371 }
00372 }
00373 }
00374 }
00375
00376
00377 void writePointCells
00378 (
00379 const polyMesh& mesh,
00380 const label pointI,
00381 const fileName& timeName
00382 )
00383 {
00384 const labelList& pCells = mesh.pointCells()[pointI];
00385
00386 labelHashSet allEdges(6*pCells.size());
00387
00388 forAll(pCells, i)
00389 {
00390 const labelList& cEdges = mesh.cellEdges()[pCells[i]];
00391
00392 forAll(cEdges, i)
00393 {
00394 allEdges.insert(cEdges[i]);
00395 }
00396 }
00397
00398
00399 fileName pFile
00400 (
00401 mesh.time().path()
00402 / "pointEdges_" + timeName + '_' + name(pointI) + ".obj"
00403 );
00404
00405 Info << "Writing pointEdges to " << pFile << endl;
00406
00407 OFstream pointStream(pFile);
00408
00409 label vertI = 0;
00410
00411 for
00412 (
00413 labelHashSet::const_iterator iter = allEdges.begin();
00414 iter != allEdges.end();
00415 ++iter
00416 )
00417 {
00418 const edge& e = mesh.edges()[iter.key()];
00419
00420 meshTools::writeOBJ(pointStream, mesh.points()[e[0]]); vertI++;
00421 meshTools::writeOBJ(pointStream, mesh.points()[e[1]]); vertI++;
00422 pointStream<< "l " << vertI-1 << ' ' << vertI << nl;
00423 }
00424 }
00425
00426
00427
00428
00429 int main(int argc, char *argv[])
00430 {
00431 timeSelector::addOptions();
00432 argList::validOptions.insert("patchFaces", "");
00433 argList::validOptions.insert("patchEdges", "");
00434 argList::validOptions.insert("cell", "cellI");
00435 argList::validOptions.insert("face", "faceI");
00436 argList::validOptions.insert("point", "pointI");
00437 argList::validOptions.insert("cellSet", "setName");
00438 argList::validOptions.insert("faceSet", "setName");
00439 # include <OpenFOAM/addRegionOption.H>
00440
00441 # include <OpenFOAM/setRootCase.H>
00442 # include <OpenFOAM/createTime.H>
00443 runTime.functionObjects().off();
00444
00445 bool patchFaces = args.optionFound("patchFaces");
00446 bool patchEdges = args.optionFound("patchEdges");
00447 bool doCell = args.optionFound("cell");
00448 bool doPoint = args.optionFound("point");
00449 bool doFace = args.optionFound("face");
00450 bool doCellSet = args.optionFound("cellSet");
00451 bool doFaceSet = args.optionFound("faceSet");
00452
00453
00454 Info<< "Writing mesh objects as .obj files such that the object"
00455 << " numbering" << endl
00456 << "(for points, faces, cells) is consistent with"
00457 << " Foam numbering (starting from 0)." << endl << endl;
00458
00459 instantList timeDirs = timeSelector::select0(runTime, args);
00460
00461 # include <OpenFOAM/createNamedPolyMesh.H>
00462
00463 forAll(timeDirs, timeI)
00464 {
00465 runTime.setTime(timeDirs[timeI], timeI);
00466
00467 Info<< "Time = " << runTime.timeName() << endl;
00468
00469 polyMesh::readUpdateState state = mesh.readUpdate();
00470
00471 if (!timeI || state != polyMesh::UNCHANGED)
00472 {
00473 if (patchFaces)
00474 {
00475 writePatchFaces(mesh, runTime.timeName());
00476 }
00477 if (patchEdges)
00478 {
00479 writePatchBoundaryEdges(mesh, runTime.timeName());
00480 }
00481 if (doCell)
00482 {
00483 label cellI = args.optionRead<label>("cell");
00484
00485 writePoints(mesh, cellI, runTime.timeName());
00486 }
00487 if (doPoint)
00488 {
00489 label pointI = args.optionRead<label>("point");
00490
00491 writePointCells(mesh, pointI, runTime.timeName());
00492 }
00493 if (doFace)
00494 {
00495 label faceI = args.optionRead<label>("face");
00496
00497 fileName fName
00498 (
00499 mesh.time().path()
00500 / "meshPoints_"
00501 + runTime.timeName()
00502 + '_'
00503 + name(faceI)
00504 + ".obj"
00505 );
00506
00507 Info<< "Writing mesh points and edges to " << fName << endl;
00508
00509 OFstream str(fName);
00510
00511 const face& f = mesh.faces()[faceI];
00512
00513 meshTools::writeOBJ(str, faceList(1, f), mesh.points());
00514 }
00515 if (doCellSet)
00516 {
00517 word setName(args.option("cellSet"));
00518
00519 cellSet cells(mesh, setName);
00520
00521 Info<< "Read " << cells.size() << " cells from set " << setName
00522 << endl;
00523
00524 writePoints(mesh, cells.toc(), runTime.timeName());
00525
00526 }
00527 if (doFaceSet)
00528 {
00529 word setName(args.option("faceSet"));
00530
00531 faceSet faces(mesh, setName);
00532
00533 Info<< "Read " << faces.size() << " faces from set " << setName
00534 << endl;
00535
00536 fileName fName
00537 (
00538 mesh.time().path()
00539 / "meshPoints_"
00540 + runTime.timeName()
00541 + '_'
00542 + setName
00543 + ".obj"
00544 );
00545
00546 Info << "Writing mesh points and edges to " << fName << endl;
00547
00548 OFstream str(fName);
00549
00550 meshTools::writeOBJ
00551 (
00552 str,
00553 mesh.faces(),
00554 mesh.points(),
00555 faces.toc()
00556 );
00557 }
00558 else if
00559 (
00560 !patchFaces
00561 && !patchEdges
00562 && !doCell
00563 && !doPoint
00564 && !doFace
00565 && !doCellSet
00566 && !doFaceSet
00567 )
00568 {
00569
00570 writePoints(mesh, runTime.timeName());
00571
00572
00573 writeFaceCentres(mesh, runTime.timeName());
00574
00575
00576 writeCellCentres(mesh, runTime.timeName());
00577
00578
00579 writePatchCentres(mesh, runTime.timeName());
00580 }
00581 }
00582 else
00583 {
00584 Info << "No mesh." << endl;
00585 }
00586
00587 Info << nl << endl;
00588 }
00589
00590
00591 Info << "End\n" << endl;
00592
00593 return 0;
00594 }
00595
00596
00597