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

STARCDMeshWriter.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include <conversion/STARCDMeshWriter.H>
00027 
00028 #include <OpenFOAM/Time.H>
00029 #include <OpenFOAM/SortableList.H>
00030 #include <OpenFOAM/OFstream.H>
00031 
00032 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00033 
00034 const char* Foam::meshWriters::STARCD::defaultBoundaryName =
00035     "Default_Boundary_Region";
00036 
00037 const Foam::label Foam::meshWriters::STARCD::foamToStarFaceAddr[4][6] =
00038 {
00039     { 4, 5, 2, 3, 0, 1 },     // 11 = pro-STAR hex
00040     { 0, 1, 4, 5, 2, -1 },    // 12 = pro-STAR prism
00041     { 5, 4, 2, 0, -1, -1 },   // 13 = pro-STAR tetra
00042     { 0, 4, 3, 5, 2, -1 }     // 14 = pro-STAR pyramid
00043 };
00044 
00045 
00046 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00047 
00048 Foam::label Foam::meshWriters::STARCD::findDefaultBoundary() const
00049 {
00050     const polyBoundaryMesh& patches = mesh_.boundaryMesh();
00051 
00052     label id = -1;
00053 
00054     // find Default_Boundary_Region if it exists
00055     forAll(patches, patchI)
00056     {
00057         if (defaultBoundaryName == patches[patchI].name())
00058         {
00059             id = patchI;
00060             break;
00061         }
00062     }
00063     return id;
00064 }
00065 
00066 
00067 void Foam::meshWriters::STARCD::getCellTable()
00068 {
00069     // read constant/polyMesh/propertyName
00070     IOList<label> ioList
00071     (
00072         IOobject
00073         (
00074             "cellTableId",
00075             "constant",
00076             polyMesh::meshSubDir,
00077             mesh_,
00078             IOobject::READ_IF_PRESENT,
00079             IOobject::NO_WRITE,
00080             false
00081         )
00082     );
00083 
00084     bool useCellZones = false;
00085     cellTableId_.setSize(mesh_.nCells(), -1);
00086 
00087     // get information from constant/polyMesh/cellTableId if possible
00088     if (ioList.headerOk())
00089     {
00090         if (ioList.size() == mesh_.nCells())
00091         {
00092             cellTableId_.transfer(ioList);
00093 
00094             if (cellTable_.empty())
00095             {
00096                 Info<< "no cellTable information available" << endl;
00097             }
00098         }
00099         else
00100         {
00101             WarningIn("STARCD::getCellTable()")
00102                 << ioList.objectPath() << " has incorrect number of cells "
00103                 << " - use cellZone information"
00104                 << endl;
00105 
00106             ioList.clear();
00107             useCellZones = true;
00108         }
00109     }
00110     else
00111     {
00112         useCellZones = true;
00113     }
00114 
00115 
00116     if (useCellZones)
00117     {
00118         if (cellTable_.empty())
00119         {
00120             Info<< "created cellTable from cellZones" << endl;
00121             cellTable_ = mesh_;
00122         }
00123 
00124         // track if there are unzoned cells
00125         label nUnzoned = mesh_.nCells();
00126 
00127         // get the cellZone <-> cellTable correspondence
00128         Info<< "matching cellZones to cellTable" << endl;
00129 
00130         forAll (mesh_.cellZones(), zoneI)
00131         {
00132             const cellZone& cZone = mesh_.cellZones()[zoneI];
00133             if (cZone.size())
00134             {
00135                 nUnzoned -= cZone.size();
00136 
00137                 label tableId = cellTable_.findIndex(cZone.name());
00138                 if (tableId < 0)
00139                 {
00140                     dictionary dict;
00141 
00142                     dict.add("Label", cZone.name());
00143                     dict.add("MaterialType", "fluid");
00144                     tableId = cellTable_.append(dict);
00145                 }
00146 
00147                 forAll (cZone, i)
00148                 {
00149                     cellTableId_[cZone[i]] = tableId;
00150                 }
00151             }
00152         }
00153 
00154         if (nUnzoned)
00155         {
00156             dictionary dict;
00157 
00158             dict.add("Label", "__unZonedCells__");
00159             dict.add("MaterialType", "fluid");
00160             label tableId = cellTable_.append(dict);
00161 
00162             forAll (cellTableId_, i)
00163             {
00164                 if (cellTableId_[i] < 0)
00165                 {
00166                     cellTableId_[i] = tableId;
00167                 }
00168             }
00169         }
00170     }
00171 }
00172 
00173 
00174 void Foam::meshWriters::STARCD::writeHeader(Ostream& os, const char* filetype)
00175 {
00176     os  << "PROSTAR_" << filetype << nl
00177         << 4000
00178         << " " << 0
00179         << " " << 0
00180         << " " << 0
00181         << " " << 0
00182         << " " << 0
00183         << " " << 0
00184         << " " << 0
00185         << endl;
00186 }
00187 
00188 
00189 void Foam::meshWriters::STARCD::writePoints(const fileName& prefix) const
00190 {
00191     OFstream os(prefix + ".vrt");
00192     writeHeader(os, "VERTEX");
00193 
00194     // Set the precision of the points data to 10
00195     os.precision(10);
00196 
00197     // force decimal point for Fortran input
00198     os.setf(std::ios::showpoint);
00199 
00200     const pointField& points = mesh_.points();
00201 
00202     Info<< "Writing " << os.name() << " : "
00203         << points.size() << " points" << endl;
00204 
00205     forAll(points, ptI)
00206     {
00207         // convert [m] -> [mm]
00208         os
00209             << ptI + 1 << " "
00210             << scaleFactor_ * points[ptI].x() << " "
00211             << scaleFactor_ * points[ptI].y() << " "
00212             << scaleFactor_ * points[ptI].z() << nl;
00213     }
00214     os.flush();
00215 
00216 }
00217 
00218 
00219 void Foam::meshWriters::STARCD::writeCells(const fileName& prefix) const
00220 {
00221     OFstream os(prefix + ".cel");
00222     writeHeader(os, "CELL");
00223 
00224     // this is what we seem to need
00225     // map foam cellModeller index -> star shape
00226     Map<label> shapeLookupIndex;
00227     shapeLookupIndex.insert(hexModel->index(), 11);
00228     shapeLookupIndex.insert(prismModel->index(), 12);
00229     shapeLookupIndex.insert(tetModel->index(), 13);
00230     shapeLookupIndex.insert(pyrModel->index(), 14);
00231 
00232     const cellShapeList& shapes = mesh_.cellShapes();
00233     const cellList& cells  = mesh_.cells();
00234     const faceList& faces  = mesh_.faces();
00235     const labelList& owner = mesh_.faceOwner();
00236 
00237     Info<< "Writing " << os.name() << " : "
00238         << cells.size() << " cells" << endl;
00239 
00240     forAll(cells, cellId)
00241     {
00242         label tableId = cellTableId_[cellId];
00243         label materialType  = 1;        // 1(fluid)
00244         if (cellTable_.found(tableId))
00245         {
00246             const dictionary& dict = cellTable_[tableId];
00247             if (dict.found("MaterialType"))
00248             {
00249                 word matType;
00250                 dict.lookup("MaterialType") >> matType;
00251                 if (matType == "solid")
00252                 {
00253                     materialType = 2;
00254                 }
00255 
00256             }
00257         }
00258 
00259         const cellShape& shape = shapes[cellId];
00260         label mapIndex = shape.model().index();
00261 
00262         // a registered primitive type
00263         if (shapeLookupIndex.found(mapIndex))
00264         {
00265             label shapeId = shapeLookupIndex[mapIndex];
00266             const labelList& vrtList = shapes[cellId];
00267 
00268             os  << cellId + 1
00269                 << " " << shapeId
00270                 << " " << vrtList.size()
00271                 << " " << tableId
00272                 << " " << materialType;
00273 
00274             // primitives have <= 8 vertices, but prevent overrun anyhow
00275             // indent following lines for ease of reading
00276             label count = 0;
00277             forAll(vrtList, i)
00278             {
00279                 if ((count % 8) == 0)
00280                 {
00281                     os  << nl
00282                         << "  " << cellId + 1;
00283                 }
00284                 os << " " << vrtList[i] + 1;
00285                 count++;
00286             }
00287             os << endl;
00288 
00289         }
00290         else
00291         {
00292             label shapeId = 255;        // treat as general polyhedral
00293             const labelList& cFaces  = cells[cellId];
00294 
00295             // create (beg,end) indices
00296             List<label> indices(cFaces.size() + 1);
00297             indices[0] = indices.size();
00298 
00299             label count = indices.size();
00300             // determine the total number of vertices
00301             forAll(cFaces, faceI)
00302             {
00303                 count += faces[cFaces[faceI]].size();
00304                 indices[faceI+1] = count;
00305             }
00306 
00307             os  << cellId + 1
00308                 << " " << shapeId
00309                 << " " << count
00310                 << " " << tableId
00311                 << " " << materialType;
00312 
00313             // write indices - max 8 per line
00314             // indent following lines for ease of reading
00315             count = 0;
00316             forAll(indices, i)
00317             {
00318                 if ((count % 8) == 0)
00319                 {
00320                     os  << nl
00321                         << "  " << cellId + 1;
00322                 }
00323                 os << " " << indices[i];
00324                 count++;
00325             }
00326 
00327             // write faces - max 8 per line
00328             forAll(cFaces, faceI)
00329             {
00330                 label meshFace = cFaces[faceI];
00331                 face f;
00332 
00333                 if (owner[meshFace] == cellId)
00334                 {
00335                     f = faces[meshFace];
00336                 }
00337                 else
00338                 {
00339                     f = faces[meshFace].reverseFace();
00340                 }
00341 
00342                 forAll(f, i)
00343                 {
00344                     if ((count % 8) == 0)
00345                     {
00346                         os  << nl
00347                             << "  " << cellId + 1;
00348                     }
00349 
00350                     os << " " << f[i] + 1;
00351                     count++;
00352                 }
00353             }
00354 
00355             os << endl;
00356         }
00357     }
00358 }
00359 
00360 
00361 void Foam::meshWriters::STARCD::writeBoundary(const fileName& prefix) const
00362 {
00363     OFstream os(prefix + ".bnd");
00364     writeHeader(os, "BOUNDARY");
00365 
00366     const cellShapeList& shapes = mesh_.cellShapes();
00367     const cellList& cells  = mesh_.cells();
00368     const faceList& faces  = mesh_.faces();
00369     const labelList& owner = mesh_.faceOwner();
00370     const polyBoundaryMesh& patches = mesh_.boundaryMesh();
00371 
00372     // this is what we seem to need
00373     // these MUST correspond to foamToStarFaceAddr
00374     //
00375     Map<label> faceLookupIndex;
00376     faceLookupIndex.insert(hexModel->index(), 0);
00377     faceLookupIndex.insert(prismModel->index(), 1);
00378     faceLookupIndex.insert(tetModel->index(), 2);
00379     faceLookupIndex.insert(pyrModel->index(), 3);
00380 
00381     Info<< "Writing " << os.name() << " : "
00382         << (mesh_.nFaces() - patches[0].start()) << " boundaries" << endl;
00383 
00384 
00385     label defaultId = findDefaultBoundary();
00386 
00387     //
00388     // write boundary faces - skip Default_Boundary_Region entirely
00389     //
00390     label boundId = 0;
00391     forAll(patches, patchI)
00392     {
00393         label regionId = patchI;
00394         if (regionId == defaultId)
00395         {
00396             continue;   // skip - already written
00397         }
00398         else if (defaultId == -1 || regionId < defaultId)
00399         {
00400             regionId++;
00401         }
00402 
00403         label patchStart = patches[patchI].start();
00404         label patchSize  = patches[patchI].size();
00405         word  bndType = boundaryRegion_.boundaryType(patches[patchI].name());
00406 
00407         for
00408         (
00409             label faceI = patchStart;
00410             faceI < (patchStart + patchSize);
00411             ++faceI
00412         )
00413         {
00414             label cellId = owner[faceI];
00415             const labelList& cFaces  = cells[cellId];
00416             const cellShape& shape = shapes[cellId];
00417             label cellFaceId = findIndex(cFaces, faceI);
00418 
00419             //      Info<< "cell " << cellId + 1 << " face " << faceI
00420             //          << " == " << faces[faceI]
00421             //          << " is index " << cellFaceId << " from " << cFaces;
00422 
00423             // Unfortunately, the order of faces returned by
00424             //   primitiveMesh::cells() is not necessarily the same
00425             //   as defined by primitiveMesh::cellShapes()
00426             // Thus, for registered primitive types, do the lookup ourselves.
00427             // Finally, the cellModel face number is re-mapped to the
00428             // STAR-CD local face number
00429 
00430             label mapIndex = shape.model().index();
00431 
00432             // a registered primitive type
00433             if (faceLookupIndex.found(mapIndex))
00434             {
00435                 const faceList sFaces = shape.faces();
00436                 forAll(sFaces, sFaceI)
00437                 {
00438                     if (faces[faceI] == sFaces[sFaceI])
00439                     {
00440                         cellFaceId = sFaceI;
00441                         break;
00442                     }
00443                 }
00444 
00445                 mapIndex = faceLookupIndex[mapIndex];
00446                 cellFaceId = foamToStarFaceAddr[mapIndex][cellFaceId];
00447             }
00448             // Info<< endl;
00449 
00450             boundId++;
00451 
00452             os
00453                 << boundId
00454                 << " " << cellId + 1
00455                 << " " << cellFaceId + 1
00456                 << " " << regionId
00457                 << " " << 0
00458                 << " " << bndType.c_str()
00459                 << endl;
00460         }
00461     }
00462 }
00463 
00464 
00465 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00466 
00467 Foam::meshWriters::STARCD::STARCD
00468 (
00469     const polyMesh& mesh,
00470     const scalar scaleFactor
00471 )
00472 :
00473     meshWriter(mesh, scaleFactor)
00474 {
00475     boundaryRegion_.readDict(mesh_);
00476     cellTable_.readDict(mesh_);
00477     getCellTable();
00478 }
00479 
00480 
00481 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00482 
00483 Foam::meshWriters::STARCD::~STARCD()
00484 {}
00485 
00486 
00487 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00488 
00489 void Foam::meshWriters::STARCD::rmFiles(const fileName& baseName) const
00490 {
00491     rm(baseName + ".vrt");
00492     rm(baseName + ".cel");
00493     rm(baseName + ".bnd");
00494     rm(baseName + ".inp");
00495 }
00496 
00497 
00498 bool Foam::meshWriters::STARCD::write(const fileName& meshName) const
00499 {
00500     fileName baseName(meshName);
00501 
00502     if (baseName.empty())
00503     {
00504         baseName = meshWriter::defaultMeshName;
00505 
00506         if
00507         (
00508             mesh_.time().timeName() != "0"
00509          && mesh_.time().timeName() != "constant"
00510         )
00511         {
00512             baseName += "_" + mesh_.time().timeName();
00513         }
00514     }
00515 
00516     rmFiles(baseName);
00517     writePoints(baseName);
00518     writeCells(baseName);
00519 
00520     if (writeBoundary_)
00521     {
00522         writeBoundary(baseName);
00523     }
00524 
00525     return true;
00526 }
00527 
00528 
00529 bool Foam::meshWriters::STARCD::writeSurface
00530 (
00531     const fileName& meshName,
00532     const bool& triangulate
00533 ) const
00534 {
00535     fileName baseName(meshName);
00536 
00537     if (baseName.empty())
00538     {
00539         baseName = meshWriter::defaultSurfaceName;
00540 
00541         if
00542         (
00543             mesh_.time().timeName() != "0"
00544          && mesh_.time().timeName() != "constant"
00545         )
00546         {
00547             baseName += "_" + mesh_.time().timeName();
00548         }
00549     }
00550 
00551     rmFiles(baseName);
00552 
00553     OFstream celFile(baseName + ".cel");
00554     writeHeader(celFile, "CELL");
00555 
00556     Info << "Writing " << celFile.name() << endl;
00557 
00558     // mesh and patch info
00559     const pointField& points = mesh_.points();
00560     const labelList& owner = mesh_.faceOwner();
00561     const faceList& meshFaces = mesh_.faces();
00562     const polyBoundaryMesh& patches = mesh_.boundaryMesh();
00563 
00564     label shapeId = 3;  // shell/baffle element
00565     label typeId  = 4;  // 4(shell)
00566 
00567     // remember which points need to be written
00568     labelHashSet pointHash;
00569 
00570     // write boundary faces as normal STAR-CD mesh
00571     if (triangulate)
00572     {
00573         // cell Id has no particular meaning - just increment
00574         // use the cellTable id from the patch Number
00575         label cellId = 0;
00576 
00577         forAll(patches, patchI)
00578         {
00579             label patchStart = patches[patchI].start();
00580             label patchSize  = patches[patchI].size();
00581 
00582             label ctableId = patchI + 1;
00583 
00584             for
00585             (
00586                 label faceI = patchStart;
00587                 faceI < (patchStart + patchSize);
00588                 ++faceI
00589             )
00590             {
00591                 const face& f = meshFaces[faceI];
00592 
00593                 label nTri = f.nTriangles(points);
00594                 faceList triFaces;
00595 
00596                 // triangulate polygons, but not quads
00597                 if (nTri <= 2)
00598                 {
00599                     triFaces.setSize(1);
00600                     triFaces[0] = f;
00601                 }
00602                 else
00603                 {
00604                     triFaces.setSize(nTri);
00605                     nTri = 0;
00606                     f.triangles(points, nTri, triFaces);
00607                 }
00608 
00609                 forAll(triFaces, faceI)
00610                 {
00611                     const labelList& vrtList = triFaces[faceI];
00612 
00613                     celFile
00614                         << cellId + 1 << " "
00615                         << shapeId << " "
00616                         << vrtList.size() << " "
00617                         << ctableId << " "
00618                         << typeId;
00619 
00620                     // must be 3 (triangle) but could be quad
00621                     label count = 0;
00622                     forAll(vrtList, i)
00623                     {
00624                         if ((count % 8) == 0)
00625                         {
00626                             celFile
00627                                 << nl
00628                                 << "  " << cellId + 1;
00629                         }
00630                         // remember which points we'll need to write
00631                         pointHash.insert(vrtList[i]);
00632                         celFile << " " << vrtList[i] + 1;
00633                         count++;
00634                     }
00635                     celFile << endl;
00636 
00637                     cellId++;
00638                 }
00639             }
00640         }
00641     }
00642     else
00643     {
00644         // cell Id is the OpenFOAM face Id
00645         // use the cellTable id from the face owner
00646         // - allows separation of parts
00647         forAll(patches, patchI)
00648         {
00649             label patchStart = patches[patchI].start();
00650             label patchSize  = patches[patchI].size();
00651 
00652             for
00653             (
00654                 label faceI = patchStart;
00655                 faceI < (patchStart + patchSize);
00656                 ++faceI
00657             )
00658             {
00659                 const labelList& vrtList = meshFaces[faceI];
00660                 label cellId = faceI;
00661 
00662                 celFile
00663                     << cellId + 1 << " "
00664                     << shapeId << " "
00665                     << vrtList.size() << " "
00666                     << cellTableId_[owner[faceI]] << " "
00667                     << typeId;
00668 
00669                 // likely <= 8 vertices, but prevent overrun anyhow
00670                 label count = 0;
00671                 forAll(vrtList, i)
00672                 {
00673                     if ((count % 8) == 0)
00674                     {
00675                         celFile
00676                             << nl
00677                             << "  " << cellId + 1;
00678                     }
00679                     // remember which points we'll need to write
00680                     pointHash.insert(vrtList[i]);
00681                     celFile << " " << vrtList[i] + 1;
00682                     count++;
00683                 }
00684                 celFile << endl;
00685             }
00686         }
00687     }
00688 
00689     OFstream vrtFile(baseName + ".vrt");
00690     writeHeader(vrtFile, "VERTEX");
00691 
00692     vrtFile.precision(10);
00693     vrtFile.setf(std::ios::showpoint);  // force decimal point for Fortran
00694 
00695     Info << "Writing " << vrtFile.name() << endl;
00696 
00697     // build sorted table of contents
00698     SortableList<label> toc(pointHash.size());
00699     {
00700         label i = 0;
00701         forAllConstIter(labelHashSet, pointHash, iter)
00702         {
00703             toc[i++] = iter.key();
00704         }
00705     }
00706     toc.sort();
00707     toc.shrink();
00708     pointHash.clear();
00709 
00710     // write points in sorted order
00711     forAll(toc, i)
00712     {
00713         label vrtId = toc[i];
00714         vrtFile
00715             << vrtId + 1
00716             << " " << scaleFactor_ * points[vrtId].x()
00717             << " " << scaleFactor_ * points[vrtId].y()
00718             << " " << scaleFactor_ * points[vrtId].z()
00719             << endl;
00720     }
00721 
00722     return true;
00723 }
00724 
00725 
00726 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines