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 #include <OpenFOAM/ListOps.H>
00058 #include <OpenFOAM/argList.H>
00059 #include <OpenFOAM/Time.H>
00060 #include <finiteVolume/fvMesh.H>
00061 #include <finiteVolume/volFields.H>
00062 #include <OpenFOAM/emptyPolyPatch.H>
00063 #include <OpenFOAM/symmetryPolyPatch.H>
00064 #include <OpenFOAM/wallPolyPatch.H>
00065 #include <OpenFOAM/SortableList.H>
00066 #include <meshTools/cellSet.H>
00067
00068 #include <libccmio/ccmio.h>
00069 #include <vector>
00070
00071 using namespace Foam;
00072
00073
00074
00075 static char const kDefaultState[] = "default";
00076 static int const kVertOffset = 2;
00077
00078
00079
00080 labelList getInternalFaceOrder
00081 (
00082 const cellList& cells,
00083 const labelList& owner,
00084 const labelList& neighbour
00085 )
00086 {
00087 labelList oldToNew(owner.size(), -1);
00088
00089
00090 label newFaceI = 0;
00091
00092 forAll(cells, cellI)
00093 {
00094 const labelList& cFaces = cells[cellI];
00095
00096 SortableList<label> nbr(cFaces.size());
00097
00098 forAll(cFaces, i)
00099 {
00100 label faceI = cFaces[i];
00101
00102 label nbrCellI = neighbour[faceI];
00103
00104 if (nbrCellI != -1)
00105 {
00106
00107 if (nbrCellI == cellI)
00108 {
00109 nbrCellI = owner[faceI];
00110 }
00111
00112 if (cellI < nbrCellI)
00113 {
00114
00115 nbr[i] = nbrCellI;
00116 }
00117 else
00118 {
00119
00120 nbr[i] = -1;
00121 }
00122 }
00123 else
00124 {
00125
00126 nbr[i] = -1;
00127 }
00128 }
00129
00130 nbr.sort();
00131
00132 forAll(nbr, i)
00133 {
00134 if (nbr[i] != -1)
00135 {
00136 oldToNew[cFaces[nbr.indices()[i]]] = newFaceI++;
00137 }
00138 }
00139 }
00140
00141
00142 for (label faceI = newFaceI; faceI < owner.size(); faceI++)
00143 {
00144 oldToNew[faceI] = faceI;
00145 }
00146
00147 return oldToNew;
00148 }
00149
00150
00151 void storeCellInZone
00152 (
00153 const label cellI,
00154 const label cellType,
00155 Map<label>& typeToZone,
00156 List<DynamicList<label> >& zoneCells
00157 )
00158 {
00159 if (cellType >= 0)
00160 {
00161 Map<label>::iterator zoneFnd = typeToZone.find(cellType);
00162
00163 if (zoneFnd == typeToZone.end())
00164 {
00165
00166 zoneCells.setSize(zoneCells.size() + 1);
00167
00168 label zoneI = zoneCells.size()-1;
00169
00170 Info<< "Mapping type " << cellType << " to Foam cellZone "
00171 << zoneI << endl;
00172 typeToZone.insert(cellType, zoneI);
00173
00174 zoneCells[zoneI].append(cellI);
00175 }
00176 else
00177 {
00178
00179 zoneCells[zoneFnd()].append(cellI);
00180 }
00181 }
00182 }
00183
00184
00185 void CheckError(CCMIOError const &err, const Foam::string& str)
00186 {
00187 if (err != kCCMIONoErr)
00188 {
00189 FatalErrorIn("CheckError")
00190 << str << exit(FatalError);
00191 }
00192 }
00193
00194
00195 void ReadVertices
00196 (
00197 CCMIOError &err,
00198 CCMIOID &vertices,
00199 labelList& foamPointMap,
00200 pointField& foamPoints
00201 )
00202 {
00203
00204
00205
00206
00207
00208
00209
00210 CCMIOSize nVertices;
00211 CCMIOEntitySize(&err, vertices, &nVertices, NULL);
00212
00213 List<int> mapData(nVertices, 0);
00214 List<float> verts(3*nVertices, 0);
00215
00216 int offset = 0;
00217 int offsetPlusSize = offset+nVertices;
00218
00219 int dims = 1;
00220 float scale;
00221 CCMIOID mapID;
00222 CCMIOReadVerticesf(&err, vertices, &dims, &scale, &mapID, verts.begin(),
00223 offset, offsetPlusSize);
00224 CCMIOReadMap(&err, mapID, mapData.begin(), offset, offsetPlusSize);
00225
00226
00227
00228
00229
00230
00231
00232
00233
00234 foamPoints.setSize(nVertices);
00235 foamPoints = vector::zero;
00236 foamPointMap.setSize(nVertices);
00237
00238 forAll(foamPointMap, i)
00239 {
00240 foamPointMap[i] = mapData[i];
00241 for (direction cmpt = 0; cmpt < dims; cmpt++)
00242 {
00243 foamPoints[i][cmpt] = verts[dims*i + cmpt]*scale;
00244 }
00245 }
00246 }
00247
00248
00249 void ReadProblem
00250 (
00251 CCMIOError &err,
00252 CCMIOID& problem,
00253
00254 const Map<label>& prostarToFoamPatch,
00255 Map<word>& foamCellTypeNames,
00256 wordList& foamPatchTypes,
00257 wordList& foamPatchNames
00258 )
00259 {
00260
00261 CCMIOID next;
00262 int i = 0;
00263 while
00264 (
00265 CCMIONextEntity(NULL, problem, kCCMIOCellType, &i, &next)
00266 == kCCMIONoErr
00267 )
00268 {
00269 char *name;
00270 int size;
00271 int cellType;
00272
00273
00274
00275
00276
00277
00278 if
00279 (
00280 CCMIOReadOptstr(NULL, next, "MaterialType", &size, NULL)
00281 == kCCMIONoErr
00282 )
00283 {
00284 name = new char[size + 1];
00285 CCMIOReadOptstr(&err, next, "MaterialType", &size, name);
00286 CCMIOGetEntityIndex(&err, next, &cellType);
00287
00288 foamCellTypeNames.insert(cellType, name);
00289 Pout<< "Celltype:" << cellType << " name:" << name << endl;
00290
00291 delete [] name;
00292 }
00293 }
00294
00295
00296
00297
00298 CCMIOID boundary;
00299 label regionI = 0;
00300 int k = 0;
00301 while
00302 (
00303 CCMIONextEntity(NULL, problem, kCCMIOBoundaryRegion, &k, &boundary)
00304 == kCCMIONoErr
00305 )
00306 {
00307
00308 label foamPatchI = -1;
00309
00310
00311
00312 int prostarI = -1;
00313 if
00314 (
00315 CCMIOReadOpti(NULL, boundary, "ProstarRegionNumber", &prostarI)
00316 == kCCMIONoErr
00317 )
00318 {
00319 Pout<< "For region:" << regionI
00320 << " found ProstarRegionNumber:" << prostarI << endl;
00321 }
00322 else
00323 {
00324 prostarI = regionI;
00325
00326 Pout<< "For region:" << regionI
00327 << "did not find ProstarRegionNumber entry. Assuming "
00328 << prostarI << endl;
00329 }
00330
00331
00332 if (prostarToFoamPatch.found(prostarI))
00333 {
00334 foamPatchI = prostarToFoamPatch[prostarI];
00335
00336
00337
00338 int size;
00339 if
00340 (
00341 CCMIOReadOptstr(NULL, boundary, "BoundaryType", &size, NULL)
00342 == kCCMIONoErr
00343 )
00344 {
00345 char* s = new char[size + 1];
00346 CCMIOReadOptstr(NULL, boundary, "BoundaryType", &size, s);
00347 s[size] = '\0';
00348 foamPatchTypes[foamPatchI] = string::validate<word>(string(s));
00349 delete [] s;
00350 }
00351
00352
00353
00354
00355
00356
00357
00358
00359
00360
00361 if
00362 (
00363 CCMIOReadOptstr(NULL, boundary, "BoundaryName", &size, NULL)
00364 == kCCMIONoErr
00365 )
00366 {
00367 char* name = new char[size + 1];
00368 CCMIOReadOptstr(NULL, boundary, "BoundaryName", &size, name);
00369 name[size] = '\0';
00370 foamPatchNames[foamPatchI] = string::validate<word>(string(name));
00371 delete [] name;
00372 }
00373 else if
00374 (
00375 CCMIOReadOptstr(NULL, boundary, "Label", &size, NULL)
00376 == kCCMIONoErr
00377 )
00378 {
00379 char* name = new char[size + 1];
00380 CCMIOReadOptstr(NULL, boundary, "Label", &size, name);
00381 name[size] = '\0';
00382 foamPatchNames[foamPatchI] = string::validate<word>(string(name));
00383 delete [] name;
00384 }
00385 else
00386 {
00387 foamPatchNames[foamPatchI] =
00388 foamPatchTypes[foamPatchI]
00389 + Foam::name(foamPatchI);
00390 Pout<< "Made up name:" << foamPatchNames[foamPatchI]
00391 << endl;
00392 }
00393
00394 Pout<< "Read patch:" << foamPatchI
00395 << " name:" << foamPatchNames[foamPatchI]
00396 << " foamPatchTypes:" << foamPatchTypes[foamPatchI]
00397 << endl;
00398 }
00399
00400 regionI++;
00401 }
00402 }
00403
00404
00405 void ReadCells
00406 (
00407 CCMIOError &err,
00408 CCMIOID& topology,
00409 labelList& foamCellMap,
00410 labelList& foamCellType,
00411 Map<label>& prostarToFoamPatch,
00412 DynamicList<label>& foamPatchSizes,
00413 DynamicList<label>& foamPatchStarts,
00414 labelList& foamFaceMap,
00415 labelList& foamOwner,
00416 labelList& foamNeighbour,
00417 faceList& foamFaces
00418 )
00419 {
00420
00421
00422
00423
00424
00425 CCMIOID id;
00426 CCMIOGetEntity(&err, topology, kCCMIOCells, 0, &id);
00427 CCMIOSize nCells;
00428 CCMIOEntitySize(&err, id, &nCells, NULL);
00429
00430 std::vector<int> mapData(nCells);
00431 std::vector<int> cellType(nCells);
00432
00433 CCMIOID mapID;
00434 CCMIOReadCells(&err, id, &mapID, &cellType[0], 0, nCells);
00435 CCMIOReadMap(&err, mapID, &mapData[0], 0, nCells);
00436 CheckError(err, "Error reading cells");
00437
00438 foamCellMap.setSize(nCells);
00439 foamCellType.setSize(nCells);
00440 forAll(foamCellMap, i)
00441 {
00442 foamCellMap[i] = mapData[i];
00443 foamCellType[i] = cellType[i];
00444 }
00445
00446
00447
00448
00449
00450 CCMIOGetEntity(&err, topology, kCCMIOInternalFaces, 0, &id);
00451 CCMIOSize nInternalFaces;
00452 CCMIOEntitySize(&err, id, &nInternalFaces, NULL);
00453 Pout<< "nInternalFaces:" << nInternalFaces << endl;
00454
00455
00456 label foamNFaces = nInternalFaces;
00457 int index = 0;
00458 while
00459 (
00460 CCMIONextEntity(NULL, topology, kCCMIOBoundaryFaces, &index, &id)
00461 == kCCMIONoErr
00462 )
00463 {
00464 CCMIOSize size;
00465 CCMIOEntitySize(&err, id, &size, NULL);
00466
00467 Pout<< "Read kCCMIOBoundaryFaces entry with " << size
00468 << " faces." << endl;
00469
00470 foamPatchStarts.append(foamNFaces);
00471 foamPatchSizes.append(size);
00472 foamNFaces += size;
00473 }
00474 foamPatchStarts.shrink();
00475 foamPatchSizes.shrink();
00476
00477 Pout<< "patchSizes:" << foamPatchSizes << endl;
00478 Pout<< "patchStarts:" << foamPatchStarts << endl;
00479 Pout<< "nFaces:" << foamNFaces << endl;
00480
00481 mapData.resize(nInternalFaces);
00482 CCMIOGetEntity(&err, topology, kCCMIOInternalFaces, 0, &id);
00483 CCMIOSize size;
00484 CCMIOReadFaces(&err, id, kCCMIOInternalFaces, NULL, &size, NULL,
00485 kCCMIOStart, kCCMIOEnd);
00486 std::vector<int> faces(size);
00487 CCMIOReadFaces(&err, id, kCCMIOInternalFaces, &mapID, NULL, &faces[0],
00488 kCCMIOStart, kCCMIOEnd);
00489 std::vector<int> faceCells(2*nInternalFaces);
00490 CCMIOReadFaceCells(&err, id, kCCMIOInternalFaces, &faceCells[0],
00491 kCCMIOStart, kCCMIOEnd);
00492 CCMIOReadMap(&err, mapID, &mapData[0], kCCMIOStart, kCCMIOEnd);
00493 CheckError(err, "Error reading internal faces");
00494
00495
00496 foamFaceMap.setSize(foamNFaces);
00497 foamFaces.setSize(foamNFaces);
00498 foamOwner.setSize(foamNFaces);
00499 foamNeighbour.setSize(foamNFaces);
00500
00501 unsigned int pos = 0;
00502
00503 for (unsigned int faceI = 0; faceI < nInternalFaces; faceI++)
00504 {
00505 foamFaceMap[faceI] = mapData[faceI];
00506 foamOwner[faceI] = faceCells[2*faceI];
00507 foamNeighbour[faceI] = faceCells[2*faceI+1];
00508 face& f = foamFaces[faceI];
00509
00510 f.setSize(faces[pos++]);
00511 forAll(f, fp)
00512 {
00513 f[fp] = faces[pos++];
00514 }
00515 }
00516
00517
00518
00519
00520 index = 0;
00521 label regionI = 0;
00522 while
00523 (
00524 CCMIONextEntity(NULL, topology, kCCMIOBoundaryFaces, &index, &id)
00525 == kCCMIONoErr
00526 )
00527 {
00528 CCMIOSize nFaces;
00529 CCMIOEntitySize(&err, id, &nFaces, NULL);
00530
00531 mapData.resize(nFaces);
00532 faceCells.resize(nFaces);
00533 CCMIOReadFaces(&err, id, kCCMIOBoundaryFaces, NULL, &size, NULL,
00534 kCCMIOStart, kCCMIOEnd);
00535 faces.resize(size);
00536 CCMIOReadFaces(&err, id, kCCMIOBoundaryFaces, &mapID, NULL, &faces[0],
00537 kCCMIOStart, kCCMIOEnd);
00538 CCMIOReadFaceCells(&err, id, kCCMIOBoundaryFaces, &faceCells[0],
00539 kCCMIOStart, kCCMIOEnd);
00540 CCMIOReadMap(&err, mapID, &mapData[0], kCCMIOStart, kCCMIOEnd);
00541 CheckError(err, "Error reading boundary faces");
00542
00543
00544 int prostarI;
00545 if
00546 (
00547 CCMIOReadOpti(NULL, id, "ProstarRegionNumber", &prostarI)
00548 == kCCMIONoErr
00549 )
00550 {
00551 Pout<< "For region:" << regionI
00552 << " found ProstarRegionNumber:" << prostarI << endl;
00553 }
00554 else
00555 {
00556 prostarI = regionI;
00557
00558 Pout<< "For region:" << regionI
00559 << " did not find ProstarRegionNumber entry. Assuming "
00560 << prostarI << endl;
00561 }
00562 prostarToFoamPatch.insert(prostarI, regionI);
00563
00564
00565 Pout<< "region:" << regionI
00566 << " ProstarRegionNumber:" << prostarI
00567 << " foamPatchStart:"
00568 << foamPatchStarts[regionI]
00569 << " size:"
00570 << foamPatchSizes[regionI]
00571 << endl;
00572
00573
00574 unsigned int pos = 0;
00575
00576 for (unsigned int i = 0; i < nFaces; i++)
00577 {
00578 label foamFaceI = foamPatchStarts[regionI] + i;
00579
00580 foamFaceMap[foamFaceI] = mapData[i];
00581 foamOwner[foamFaceI] = faceCells[i];
00582 foamNeighbour[foamFaceI] = -1;
00583 face& f = foamFaces[foamFaceI];
00584
00585 f.setSize(faces[pos++]);
00586 forAll(f, fp)
00587 {
00588 f[fp] = faces[pos++];
00589 }
00590 }
00591 Pout<< endl;
00592
00593 regionI++;
00594 }
00595 }
00596
00597
00598
00599
00600 int main(int argc, char *argv[])
00601 {
00602 argList::noParallel();
00603 argList::validArgs.append("ccm26 file");
00604
00605 # include <OpenFOAM/setRootCase.H>
00606 # include <OpenFOAM/createTime.H>
00607
00608
00609
00610
00611
00612
00613 pointField foamPoints;
00614 labelList foamPointMap;
00615
00616
00617 labelList foamCellType;
00618 labelList foamCellMap;
00619
00620
00621 Map<label> prostarToFoamPatch;
00622 DynamicList<label> foamPatchSizes;
00623 DynamicList<label> foamPatchStarts;
00624
00625 labelList foamFaceMap;
00626 labelList foamOwner;
00627 labelList foamNeighbour;
00628 faceList foamFaces;
00629
00630
00631
00632 Map<word> foamCellTypeNames;
00633 wordList foamPatchTypes;
00634 wordList foamPatchNames;
00635
00636 {
00637 fileName ccmFile(args.additionalArgs()[0]);
00638
00639 if (!isFile(ccmFile))
00640 {
00641 FatalErrorIn(args.executable())
00642 << "Cannot read file " << ccmFile
00643 << exit(FatalError);
00644 }
00645
00646 word ccmExt = ccmFile.ext();
00647
00648 if (ccmExt != "ccm" && ccmExt != "ccmg")
00649 {
00650 FatalErrorIn(args.executable())
00651 << "Illegal extension " << ccmExt << " for file " << ccmFile
00652 << nl << "Allowed extensions are '.ccm', '.ccmg'"
00653 << exit(FatalError);
00654 }
00655
00656
00657
00658
00659 CCMIOID root;
00660 CCMIOError err = CCMIOOpenFile(NULL, ccmFile.c_str(), kCCMIORead, &root);
00661
00662
00663
00664
00665 CCMIOID state;
00666 int stateI = 0;
00667 CCMIONextEntity(&err, root, kCCMIOState, &stateI, &state);
00668 CheckError(err, "Error opening state");
00669
00670 unsigned int size;
00671 CCMIOEntityDescription(&err, state, &size, NULL);
00672 char *desc = new char[size + 1];
00673 CCMIOEntityDescription(&err, state, NULL, desc);
00674 Pout<< "Reading state '" << kDefaultState << "' (" << desc << ")"
00675 << endl;
00676 delete [] desc;
00677
00678
00679
00680 int i = 0;
00681 CCMIOID processor;
00682 CCMIONextEntity(&err, state, kCCMIOProcessor, &i, &processor);
00683 CCMIOID solution, vertices, topology;
00684 CCMIOReadProcessor
00685 (
00686 &err,
00687 processor,
00688 &vertices,
00689 &topology,
00690 NULL,
00691 &solution
00692 );
00693
00694 if (err != kCCMIONoErr)
00695 {
00696
00697 err = kCCMIONoErr;
00698 CCMIOReadProcessor
00699 (
00700 &err,
00701 processor,
00702 &vertices,
00703 &topology,
00704 NULL,
00705 NULL
00706 );
00707 if (err != kCCMIONoErr)
00708 {
00709 FatalErrorIn(args.executable())
00710 << "Could not read the file."
00711 << exit(FatalError);
00712 }
00713 }
00714
00715 ReadVertices(err, vertices, foamPointMap, foamPoints);
00716
00717 Pout<< "nPoints:" << foamPoints.size() << endl
00718 << "bounding box:" << boundBox(foamPoints) << endl
00719 << endl;
00720
00721 ReadCells
00722 (
00723 err,
00724 topology,
00725 foamCellMap,
00726 foamCellType,
00727 prostarToFoamPatch,
00728 foamPatchSizes,
00729 foamPatchStarts,
00730 foamFaceMap,
00731 foamOwner,
00732 foamNeighbour,
00733 foamFaces
00734 );
00735
00736 Pout<< "nCells:" << foamCellMap.size() << endl
00737 << "nFaces:" << foamOwner.size() << endl
00738 << "nPatches:" << foamPatchStarts.size() << endl
00739 << "nInternalFaces:" << foamPatchStarts[0] << endl
00740 << endl;
00741
00742
00743
00744 foamPatchTypes.setSize(foamPatchStarts.size());
00745 foamPatchNames.setSize(foamPatchStarts.size());
00746
00747 forAll(foamPatchNames, i)
00748 {
00749 foamPatchNames[i] = word("patch") + Foam::name(i);
00750 foamPatchTypes[i] = "patch";
00751 }
00752
00753
00754
00755 CCMIOID problem;
00756 int problemI = 0;
00757 CCMIONextEntity
00758 (
00759 &err,
00760 root,
00761 kCCMIOProblemDescription,
00762 &problemI,
00763 &problem
00764 );
00765 CheckError(err, "Error stepping to first problem description");
00766
00767 if (CCMIOIsValidEntity(problem))
00768 {
00769 ReadProblem
00770 (
00771 err,
00772 problem,
00773 prostarToFoamPatch,
00774
00775 foamCellTypeNames,
00776 foamPatchTypes,
00777 foamPatchNames
00778 );
00779 }
00780
00781
00782 CCMIOCloseFile(&err, vertices);
00783 CCMIOCloseFile(&err, topology);
00784 CCMIOCloseFile(&err, solution);
00785 CCMIOCloseFile(&err, root);
00786 }
00787
00788
00789 Pout<< "foamPatchNames:" << foamPatchNames << endl;
00790
00791
00792 Pout<< "foamOwner : min:" << min(foamOwner)
00793 << " max:" << max(foamOwner)
00794 << nl
00795 << "foamNeighbour : min:" << min(foamNeighbour)
00796 << " max:" << max(foamNeighbour)
00797 << nl
00798 << "foamCellType : min:" << min(foamCellType)
00799 << " max:" << max(foamCellType)
00800 << nl << endl;
00801
00802
00803
00804
00805
00806
00807
00808
00809
00810
00811
00812 {
00813 label maxCCMPointI = max(foamPointMap);
00814 labelList toFoamPoints(invert(maxCCMPointI+1, foamPointMap));
00815
00816 forAll(foamFaces, faceI)
00817 {
00818 inplaceRenumber(toFoamPoints, foamFaces[faceI]);
00819 }
00820 }
00821
00822
00823 {
00824 label maxCCMCellI = max(foamCellMap);
00825 labelList toFoamCells(invert(maxCCMCellI+1, foamCellMap));
00826
00827 inplaceRenumber(toFoamCells, foamOwner);
00828 inplaceRenumber(toFoamCells, foamNeighbour);
00829 }
00830
00831
00832
00833
00834
00835
00836
00837
00838
00839
00840
00841
00842
00843 forAll(foamNeighbour, faceI)
00844 {
00845 label nbr = foamNeighbour[faceI];
00846 label own = foamOwner[faceI];
00847
00848 if (nbr >= foamCellType.size() || own >= foamCellType.size())
00849 {
00850 FatalErrorIn(args.executable())
00851 << "face:" << faceI
00852 << " nbr:" << nbr
00853 << " own:" << own
00854 << " nCells:" << foamCellType.size()
00855 << exit(FatalError);
00856 }
00857
00858 if (nbr >= 0)
00859 {
00860 if (nbr < own)
00861 {
00862 foamOwner[faceI] = foamNeighbour[faceI];
00863 foamNeighbour[faceI] = own;
00864 foamFaces[faceI] = foamFaces[faceI].reverseFace();
00865 }
00866 }
00867
00868
00869
00870 const face& f = foamFaces[faceI];
00871
00872 forAll(f, fp)
00873 {
00874 if (f[fp] < 0 || f[fp] >= foamPoints.size())
00875 {
00876 FatalErrorIn(args.executable()) << "face:" << faceI
00877 << " f:" << f << abort(FatalError);
00878 }
00879 }
00880 }
00881
00882
00883
00884
00885
00886 {
00887
00888 cellList foamCells;
00889 primitiveMesh::calcCells
00890 (
00891 foamCells,
00892 foamOwner,
00893 foamNeighbour,
00894 foamCellType.size()
00895 );
00896
00897
00898 labelList oldToNew
00899 (
00900 getInternalFaceOrder
00901 (
00902 foamCells,
00903 foamOwner,
00904 foamNeighbour
00905 )
00906 );
00907
00908
00909 forAll(foamCells, cellI)
00910 {
00911 inplaceRenumber(oldToNew, foamCells[cellI]);
00912 }
00913
00914
00915 inplaceReorder(oldToNew, foamFaces);
00916 inplaceReorder(oldToNew, foamOwner);
00917 inplaceReorder(oldToNew, foamNeighbour);
00918 }
00919
00920
00921
00922
00923
00924 fvMesh mesh
00925 (
00926 IOobject
00927 (
00928 fvMesh::defaultRegion,
00929 runTime.constant(),
00930 runTime
00931 ),
00932 xferMove<pointField>(foamPoints),
00933 xferMove<faceList>(foamFaces),
00934 xferCopy<labelList>(foamOwner),
00935 xferMove<labelList>(foamNeighbour)
00936 );
00937
00938
00939 List<polyPatch*> newPatches(foamPatchNames.size());
00940
00941 label meshFaceI = foamPatchStarts[0];
00942
00943 forAll(newPatches, patchI)
00944 {
00945 const word& patchName = foamPatchNames[patchI];
00946 const word& patchType = foamPatchTypes[patchI];
00947
00948 Pout<< "Patch:" << patchName << " start at:" << meshFaceI
00949 << " size:" << foamPatchSizes[patchI]
00950 << " end at:" << meshFaceI+foamPatchSizes[patchI]
00951 << endl;
00952
00953 if (patchType == "wall")
00954 {
00955 newPatches[patchI] =
00956 new wallPolyPatch
00957 (
00958 patchName,
00959 foamPatchSizes[patchI],
00960 meshFaceI,
00961 patchI,
00962 mesh.boundaryMesh()
00963 );
00964 }
00965 else if (patchType == "symmetryplane")
00966 {
00967 newPatches[patchI] =
00968 new symmetryPolyPatch
00969 (
00970 patchName,
00971 foamPatchSizes[patchI],
00972 meshFaceI,
00973 patchI,
00974 mesh.boundaryMesh()
00975 );
00976 }
00977 else if (patchType == "empty")
00978 {
00979
00980 newPatches[patchI] =
00981 new emptyPolyPatch
00982 (
00983 patchName,
00984 foamPatchSizes[patchI],
00985 meshFaceI,
00986 patchI,
00987 mesh.boundaryMesh()
00988 );
00989 }
00990 else
00991 {
00992
00993
00994 newPatches[patchI] =
00995 new polyPatch
00996 (
00997 patchName,
00998 foamPatchSizes[patchI],
00999 meshFaceI,
01000 patchI,
01001 mesh.boundaryMesh()
01002 );
01003 }
01004
01005 meshFaceI += foamPatchSizes[patchI];
01006 }
01007
01008 if (meshFaceI != foamOwner.size())
01009 {
01010 FatalErrorIn(args.executable())
01011 << "meshFaceI:" << meshFaceI
01012 << " nFaces:" << foamOwner.size()
01013 << abort(FatalError);
01014 }
01015 mesh.addFvPatches(newPatches);
01016
01017
01018
01019
01020
01021
01022 label maxType = max(foamCellType);
01023 label minType = min(foamCellType);
01024
01025 if (maxType > minType)
01026 {
01027
01028 Map<label> typeToZone;
01029
01030 List<DynamicList<label> > zoneCells(0);
01031
01032 forAll(foamCellType, cellI)
01033 {
01034 storeCellInZone
01035 (
01036 cellI,
01037 foamCellType[cellI],
01038 typeToZone,
01039 zoneCells
01040 );
01041 }
01042
01043 mesh.cellZones().clear();
01044 mesh.cellZones().setSize(typeToZone.size());
01045
01046 label nValidCellZones = 0;
01047
01048 forAllConstIter(Map<label>, typeToZone, iter)
01049 {
01050 label type = iter.key();
01051 label zoneI = iter();
01052
01053 zoneCells[zoneI].shrink();
01054
01055 word zoneName = "cellZone_" + name(type);
01056
01057 Info<< "Writing zone " << type
01058 << " size " << zoneCells[zoneI].size()
01059 << " to cellZone "
01060 << zoneName << " and cellSet " << zoneName
01061 << endl;
01062
01063 cellSet cset(mesh, zoneName, zoneCells[zoneI]);
01064 cset.write();
01065
01066 mesh.cellZones().set
01067 (
01068 nValidCellZones,
01069 new cellZone
01070 (
01071 zoneName,
01072 zoneCells[zoneI],
01073 nValidCellZones,
01074 mesh.cellZones()
01075 )
01076 );
01077 nValidCellZones++;
01078 }
01079 mesh.cellZones().writeOpt() = IOobject::AUTO_WRITE;
01080 }
01081
01082
01083 Info<< "Writing mesh to " << mesh.objectRegistry::objectPath()
01084 << "..." << nl << endl;
01085
01086
01087
01088 volScalarField cellIdField
01089 (
01090 IOobject
01091 (
01092 "cellId",
01093 runTime.timeName(),
01094 mesh,
01095 IOobject::NO_READ,
01096 IOobject::AUTO_WRITE
01097 ),
01098 mesh,
01099 dimensionedScalar("cellId", dimless, 0.0)
01100 );
01101
01102 forAll(foamCellMap, cellI)
01103 {
01104 cellIdField[cellI] = foamCellMap[cellI];
01105 }
01106
01107
01108 volScalarField cellTypeField
01109 (
01110 IOobject
01111 (
01112 "cellType",
01113 runTime.timeName(),
01114 mesh,
01115 IOobject::NO_READ,
01116 IOobject::AUTO_WRITE
01117 ),
01118 mesh,
01119 dimensionedScalar("cellType", dimless, 0.0)
01120 );
01121
01122 forAll(foamCellType, cellI)
01123 {
01124 cellTypeField[cellI] = foamCellType[cellI];
01125 }
01126
01127 Info<< "Writing cellIds as volScalarField to " << cellIdField.objectPath()
01128 << "..." << nl << endl;
01129 mesh.write();
01130
01131 Info<< "End\n" << endl;
01132
01133 return 0;
01134 }
01135
01136
01137