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

readKivaGrid.H

Go to the documentation of this file.
00001 ifstream kivaFile(kivaFileName.c_str());
00002 
00003 if (!kivaFile.good())
00004 {
00005     FatalErrorIn(args.executable())
00006         << "Cannot open file " << kivaFileName
00007         << exit(FatalError);
00008 }
00009 
00010 Info << "Reading kiva grid from file " << kivaFileName << endl;
00011 
00012 
00013 char title[120];
00014 kivaFile.getline(title, 120, '\n');
00015 
00016 label nPoints, nCells, nRegs;
00017 kivaFile >> nCells >> nPoints >> nRegs;
00018 
00019 pointField points(nPoints);
00020 label i4;
00021 labelList idface(nPoints), fv(nPoints);
00022 
00023 for(label i=0; i<nPoints; i++)
00024 {
00025     scalar ffv;
00026     kivaFile
00027         >> i4
00028         >> points[i].x() >> points[i].y() >> points[i].z()
00029         >> ffv;
00030 
00031     if (kivaVersion == kiva3v)
00032     {
00033         kivaFile >> idface[i];
00034     }
00035 
00036     // Convert from KIVA cgs units to SI
00037     points[i] *= 0.01;
00038 
00039     fv[i] = label(ffv);
00040 }
00041 
00042 labelList i1tab(nPoints), i3tab(nPoints), i8tab(nPoints), idreg(nPoints),
00043       f(nPoints), bcl(nPoints), bcf(nPoints), bcb(nPoints);
00044 
00045 label nBfaces = 0;
00046 
00047 for(label i=0; i<nPoints; i++)
00048 {
00049     label i1, i3, i8;
00050     scalar ff, fbcl, fbcf, fbcb;
00051 
00052     kivaFile
00053         >> i1 >> i3 >> i4 >> i8
00054         >> ff >> fbcl >> fbcf >> fbcb
00055         >> idreg[i];
00056 
00057     // Correct indices to start from 0
00058     i1tab[i] = i1 - 1;
00059     i3tab[i] = i3 - 1;
00060     i8tab[i] = i8 - 1;
00061 
00062     // Convert scalar indices into hexLabels
00063     f[i] = label(ff);
00064     bcl[i] = label(fbcl);
00065     bcf[i] = label(fbcf);
00066     bcb[i] = label(fbcb);
00067 
00068     if (bcl[i] > 0 && bcl[i] != 4) nBfaces++;
00069     if (bcf[i] > 0 && bcf[i] != 4) nBfaces++;
00070     if (bcb[i] > 0 && bcb[i] != 4) nBfaces++;
00071 }
00072 
00073 
00074 label mTable;
00075 kivaFile >> mTable;
00076 
00077 if (mTable == 0)
00078 {
00079     FatalErrorIn(args.executable())
00080         << "mTable == 0. This is not supported."
00081            " kivaToFoam needs complete neighbour information"
00082         << exit(FatalError);
00083 }
00084 
00085 
00086 labelList imtab(nPoints), jmtab(nPoints), kmtab(nPoints);
00087 
00088 for(label i=0; i<nPoints; i++)
00089 {
00090     label im, jm, km;
00091     kivaFile >> i4 >> im >> jm >> km;
00092 
00093     // Correct indices to start from 0
00094     imtab[i] = im - 1;
00095     jmtab[i] = jm - 1;
00096     kmtab[i] = km - 1;
00097 }
00098 
00099 Info << "Finished reading KIVA file" << endl;
00100 
00101 cellShapeList cellShapes(nPoints);
00102 
00103 labelList cellZoning(nPoints, -1);
00104 
00105 const cellModel& hex = *(cellModeller::lookup("hex"));
00106 labelList hexLabels(8);
00107 label activeCells = 0;
00108 
00109 // Create and set the collocated point collapse map
00110 labelList pointMap(nPoints);
00111 
00112 forAll (pointMap, i)
00113 {
00114     pointMap[i] = i;
00115 }
00116 
00117 // Initialise all cells to hex and search and set map for collocated points
00118 for(label i=0; i<nPoints; i++)
00119 {
00120     if (f[i] > 0.0)
00121     {
00122         hexLabels[0] = i;
00123         hexLabels[1] = i1tab[i];
00124         hexLabels[2] = i3tab[i1tab[i]];
00125         hexLabels[3] = i3tab[i];
00126         hexLabels[4] = i8tab[i];
00127         hexLabels[5] = i1tab[i8tab[i]];
00128         hexLabels[6] = i3tab[i1tab[i8tab[i]]];
00129         hexLabels[7] = i3tab[i8tab[i]];
00130 
00131         cellShapes[activeCells] = cellShape(hex, hexLabels);
00132 
00133         edgeList edges = cellShapes[activeCells].edges();
00134 
00135         forAll (edges, ei)
00136         {
00137             if (edges[ei].mag(points) < SMALL)
00138             {
00139                 label start = pointMap[edges[ei].start()];
00140                 while (start != pointMap[start])
00141                 {
00142                     start = pointMap[start];
00143                 }
00144 
00145                 label end = pointMap[edges[ei].end()];
00146                 while (end != pointMap[end])
00147                 {
00148                     end = pointMap[end];
00149                 }
00150 
00151                 label minLabel = min(start, end);
00152 
00153                 pointMap[start] = pointMap[end] = minLabel;
00154             }
00155         }
00156 
00157         // Grab the cell zoning info
00158         cellZoning[activeCells] = idreg[i];
00159 
00160         activeCells++;
00161     }
00162 }
00163 
00164 // Contract list of cells to the active ones
00165 cellShapes.setSize(activeCells);
00166 cellZoning.setSize(activeCells);
00167 
00168 // Map collocated points to refer to the same point and collapse cell shape
00169 // to the corresponding hex-degenerate.
00170 forAll (cellShapes, celli)
00171 {
00172     cellShape& cs = cellShapes[celli];
00173 
00174     forAll (cs, i)
00175     {
00176         cs[i] = pointMap[cs[i]];
00177     }
00178 
00179     cs.collapse();
00180 }
00181 
00182 label bcIDs[11] = {-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};
00183 
00184 const label nBCs = 12;
00185 
00186 const word* kivaPatchTypes[nBCs] =
00187 {
00188     &wallPolyPatch::typeName,
00189     &wallPolyPatch::typeName,
00190     &wallPolyPatch::typeName,
00191     &wallPolyPatch::typeName,
00192     &symmetryPolyPatch::typeName,
00193     &wedgePolyPatch::typeName,
00194     &polyPatch::typeName,
00195     &polyPatch::typeName,
00196     &polyPatch::typeName,
00197     &polyPatch::typeName,
00198     &symmetryPolyPatch::typeName,
00199     &cyclicPolyPatch::typeName
00200 };
00201 
00202 enum patchTypeNames
00203 {
00204     PISTON,
00205     VALVE,
00206     LINER,
00207     CYLINDERHEAD,
00208     AXIS,
00209     WEDGE,
00210     INFLOW,
00211     OUTFLOW,
00212     PRESIN,
00213     PRESOUT,
00214     SYMMETRYPLANE,
00215     CYCLIC
00216 };
00217 
00218 const char* kivaPatchNames[nBCs] =
00219 {
00220     "piston",
00221     "valve",
00222     "liner",
00223     "cylinderHead",
00224     "axis",
00225     "wedge",
00226     "inflow", 
00227     "outflow",
00228     "presin",
00229     "presout",
00230     "symmetryPlane",
00231     "cyclic"
00232 };
00233 
00234 
00235 List<SLList<face> > pFaces[nBCs];
00236 
00237 face quadFace(4);
00238 face triFace(3);
00239 
00240 for(label i=0; i<nPoints; i++)
00241 {
00242     if (f[i] > 0)
00243     {
00244         // left
00245         label bci = bcl[i];
00246         if (bci != 4)
00247         {
00248             quadFace[0] = i3tab[i];
00249             quadFace[1] = i;
00250             quadFace[2] = i8tab[i];
00251             quadFace[3] = i3tab[i8tab[i]];
00252 
00253 #           include "checkPatch.H"
00254         }
00255 
00256         // right
00257         bci = bcl[i1tab[i]];
00258         if (bci != 4)
00259         {
00260             quadFace[0] = i1tab[i];
00261             quadFace[1] = i3tab[i1tab[i]];
00262             quadFace[2] = i8tab[i3tab[i1tab[i]]];
00263             quadFace[3] = i8tab[i1tab[i]];
00264 
00265 #           include "checkPatch.H"
00266         }
00267 
00268         // front
00269         bci = bcf[i];
00270         if (bci != 4)
00271         {
00272             quadFace[0] = i;
00273             quadFace[1] = i1tab[i];
00274             quadFace[2] = i8tab[i1tab[i]];
00275             quadFace[3] = i8tab[i];
00276 
00277 #           include "checkPatch.H"
00278         }
00279 
00280         // derriere (back)
00281         bci = bcf[i3tab[i]];
00282         if (bci != 4)
00283         {
00284             quadFace[0] = i3tab[i];
00285             quadFace[1] = i8tab[i3tab[i]];
00286             quadFace[2] = i8tab[i1tab[i3tab[i]]];
00287             quadFace[3] = i1tab[i3tab[i]];
00288 
00289 #           include "checkPatch.H"
00290         }
00291 
00292         // bottom
00293         bci = bcb[i];
00294         if (bci != 4)
00295         {
00296             quadFace[0] = i;
00297             quadFace[1] = i3tab[i];
00298             quadFace[2] = i1tab[i3tab[i]];
00299             quadFace[3] = i1tab[i];
00300 
00301 #           include "checkPatch.H"
00302         }
00303 
00304         // top
00305         bci = bcb[i8tab[i]];
00306         if (bci != 4)
00307         {
00308             quadFace[0] = i8tab[i];
00309             quadFace[1] = i1tab[i8tab[i]];
00310             quadFace[2] = i3tab[i1tab[i8tab[i]]];
00311             quadFace[3] = i3tab[i8tab[i]];
00312 
00313 #           include "checkPatch.H"
00314         }
00315     }
00316 }
00317 
00318 // Tranfer liner faces that are above the minimum cylinder-head z height
00319 // into the cylinder-head region
00320 if
00321 (
00322     pFaces[LINER].size()
00323  && pFaces[LINER][0].size()
00324  && pFaces[CYLINDERHEAD].size()
00325  && pFaces[CYLINDERHEAD][0].size()
00326 )
00327 {
00328     scalar minz = GREAT;
00329 
00330     for
00331     (
00332         SLList<face>::iterator iter = pFaces[CYLINDERHEAD][0].begin();
00333         iter != pFaces[CYLINDERHEAD][0].end();
00334         ++iter
00335     )
00336     {
00337         const face& pf = iter();
00338 
00339         forAll (pf, pfi)
00340         {
00341             minz = min(minz, points[pf[pfi]].z());
00342         }
00343     }
00344 
00345     minz += SMALL;
00346 
00347     SLList<face> newLinerFaces;
00348 
00349     for
00350     (
00351         SLList<face>::iterator iter = pFaces[LINER][0].begin();
00352         iter != pFaces[LINER][0].end();
00353         ++iter
00354     )
00355     {
00356         const face& pf = iter();
00357 
00358         scalar minfz = GREAT;
00359         forAll (pf, pfi)
00360         {
00361             minfz = min(minfz, points[pf[pfi]].z());
00362         }
00363 
00364         if (minfz > minz)
00365         {
00366             pFaces[CYLINDERHEAD][0].append(pf);
00367         }
00368         else
00369         {
00370             newLinerFaces.append(pf);
00371         }
00372     }
00373 
00374     if (pFaces[LINER][0].size() != newLinerFaces.size())
00375     {
00376         Info<< "Transfered " << pFaces[LINER][0].size() - newLinerFaces.size()
00377             << " faces from liner region to cylinder head" << endl;
00378         pFaces[LINER][0] = newLinerFaces;
00379     }
00380 
00381     SLList<face> newCylinderHeadFaces;
00382 
00383     for
00384     (
00385         SLList<face>::iterator iter = pFaces[CYLINDERHEAD][0].begin();
00386         iter != pFaces[CYLINDERHEAD][0].end();
00387         ++iter
00388     )
00389     {
00390         const face& pf = iter();
00391 
00392         scalar minfz = GREAT;
00393         forAll (pf, pfi)
00394         {
00395             minfz = min(minfz, points[pf[pfi]].z());
00396         }
00397 
00398         if (minfz < zHeadMin)
00399         {
00400             pFaces[LINER][0].append(pf);
00401         }
00402         else
00403         {
00404             newCylinderHeadFaces.append(pf);
00405         }
00406     }
00407 
00408     if (pFaces[CYLINDERHEAD][0].size() != newCylinderHeadFaces.size())
00409     {
00410         Info<< "Transfered faces from cylinder-head region to linder" << endl;
00411         pFaces[CYLINDERHEAD][0] = newCylinderHeadFaces;
00412     }
00413 }
00414 
00415 
00416 // Count the number of non-zero sized patches
00417 label nPatches = 0;
00418 for (int bci=0; bci<nBCs; bci++)
00419 {
00420     forAll (pFaces[bci], rgi)
00421     {
00422         if (pFaces[bci][rgi].size())
00423         {
00424             nPatches++;
00425         }
00426     }
00427 }
00428 
00429 
00430 // Sort-out the 2D/3D wedge geometry
00431 if (pFaces[WEDGE].size() && pFaces[WEDGE][0].size())
00432 {
00433     if (pFaces[WEDGE][0].size() == cellShapes.size())
00434     {
00435         // Distribute the points to be +/- 2.5deg from the x-z plane
00436 
00437         scalar tanTheta = Foam::tan(2.5*mathematicalConstant::pi/180.0);
00438 
00439         SLList<face>::iterator iterf = pFaces[WEDGE][0].begin();
00440         SLList<face>::iterator iterb = pFaces[WEDGE][1].begin();
00441         for
00442         (
00443             ;
00444             iterf != pFaces[WEDGE][0].end(), iterb != pFaces[WEDGE][1].end();
00445             ++iterf, ++iterb
00446         )
00447         {
00448             for (direction d=0; d<4; d++)
00449             {
00450                 points[iterf()[d]].y() = -tanTheta*points[iterf()[d]].x();
00451                 points[iterb()[d]].y() =  tanTheta*points[iterb()[d]].x();
00452             }
00453         }
00454     }
00455     else
00456     {
00457         pFaces[CYCLIC].setSize(1);
00458         pFaces[CYCLIC][0] = pFaces[WEDGE][0];
00459         for
00460         (
00461             SLList<face>::iterator iterb = pFaces[WEDGE][1].begin();
00462             iterb != pFaces[WEDGE][1].end();
00463             ++iterb
00464         )
00465         {
00466             pFaces[CYCLIC][0].append(iterb());
00467         }
00468 
00469         pFaces[WEDGE].clear();
00470         nPatches--;
00471     }
00472 }
00473 
00474 
00475 // Build the patches
00476 
00477 faceListList boundary(nPatches);
00478 wordList patchNames(nPatches);
00479 wordList patchTypes(nPatches);
00480 word defaultFacesName = "defaultFaces";
00481 word defaultFacesType = emptyPolyPatch::typeName;
00482 wordList patchPhysicalTypes(nPatches);
00483 
00484 label nAddedPatches = 0;
00485 
00486 for (int bci=0; bci<nBCs; bci++)
00487 {
00488     forAll (pFaces[bci], rgi)
00489     {
00490         if (pFaces[bci][rgi].size())
00491         {
00492             patchTypes[nAddedPatches] = *kivaPatchTypes[bci];
00493 
00494             patchNames[nAddedPatches] = kivaPatchNames[bci];
00495 
00496             if (pFaces[bci].size() > 1)
00497             {
00498                 patchNames[nAddedPatches] += name(rgi+1);
00499             }
00500 
00501             boundary[nAddedPatches] = pFaces[bci][rgi];
00502             nAddedPatches++;
00503         }
00504     }
00505 }
00506 
00507 
00508 // Remove unused points and update cell and face labels accordingly
00509 
00510 labelList pointLabels(nPoints, -1);
00511 
00512 // Scan cells for used points
00513 forAll (cellShapes, celli)
00514 {
00515     forAll (cellShapes[celli], i)
00516     {
00517         pointLabels[cellShapes[celli][i]] = 1;
00518     }
00519 }
00520 
00521 // Create addressing for used points and pack points array
00522 label newPointi = 0;
00523 forAll (pointLabels, pointi)
00524 {
00525     if (pointLabels[pointi] != -1)
00526     {
00527         pointLabels[pointi] = newPointi;
00528         points[newPointi++] = points[pointi];
00529     }
00530 }
00531 points.setSize(newPointi);
00532 
00533 // Reset cell point labels
00534 forAll (cellShapes, celli)
00535 {
00536     cellShape& cs = cellShapes[celli];
00537 
00538     forAll (cs, i)
00539     {
00540         cs[i] = pointLabels[cs[i]];
00541     }
00542 }
00543 
00544 // Reset boundary-face point labels
00545 forAll (boundary, patchi)
00546 {
00547     forAll (boundary[patchi], facei)
00548     {
00549         face& f = boundary[patchi][facei];
00550 
00551         forAll (f, i)
00552         {
00553             f[i] = pointLabels[f[i]];
00554         }
00555     }
00556 }
00557 
00558 preservePatchTypes
00559 (
00560     runTime,
00561     runTime.constant(),
00562     polyMesh::defaultRegion,
00563     patchNames,
00564     patchTypes,
00565     defaultFacesName,
00566     defaultFacesType,
00567     patchPhysicalTypes
00568 );
00569 
00570 // Build the mesh and write it out
00571 polyMesh pShapeMesh
00572 (
00573     IOobject
00574     (
00575         polyMesh::defaultRegion,
00576         runTime.constant(),
00577         runTime
00578     ),
00579     xferMove(points),
00580     cellShapes,
00581     boundary,
00582     patchNames,
00583     patchTypes,
00584     defaultFacesName,
00585     defaultFacesType,
00586     patchPhysicalTypes
00587 );
00588 
00589 Info << "Writing polyMesh" << endl;
00590 pShapeMesh.write();
00591 
00592 fileName czPath
00593 (
00594     runTime.path()/runTime.constant()/polyMesh::defaultRegion/"cellZoning"
00595 );
00596 Info << "Writing cell zoning info to file: " << czPath << endl;
00597 OFstream cz(czPath);
00598 
00599 cz << cellZoning << endl;
00600 
00601 
00602 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines