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/IOobjectList.H>
00078 #include <finiteVolume/fvMesh.H>
00079 #include <dynamicMesh/polyTopoChange.H>
00080 #include <OpenFOAM/ReadFields.H>
00081 #include <finiteVolume/volFields.H>
00082 #include <finiteVolume/surfaceFields.H>
00083 #include <OpenFOAM/bandCompression.H>
00084 #include <meshTools/faceSet.H>
00085 #include <OpenFOAM/SortableList.H>
00086 #include <decompositionMethods/decompositionMethod.H>
00087 #include <finiteVolume/fvMeshSubset.H>
00088 #include <finiteVolume/zeroGradientFvPatchFields.H>
00089
00090 using namespace Foam;
00091
00092
00093
00094 label getBand(const labelList& owner, const labelList& neighbour)
00095 {
00096 label band = 0;
00097
00098 forAll(neighbour, faceI)
00099 {
00100 label diff = neighbour[faceI] - owner[faceI];
00101
00102 if (diff > band)
00103 {
00104 band = diff;
00105 }
00106 }
00107 return band;
00108 }
00109
00110
00111
00112 labelList regionBandCompression
00113 (
00114 const fvMesh& mesh,
00115 const labelList& cellToRegion
00116 )
00117 {
00118 Pout<< "Determining cell order:" << endl;
00119
00120 labelList cellOrder(cellToRegion.size());
00121
00122 label nRegions = max(cellToRegion)+1;
00123
00124 labelListList regionToCells(invertOneToMany(nRegions, cellToRegion));
00125
00126 label cellI = 0;
00127
00128 forAll(regionToCells, regionI)
00129 {
00130 Pout<< " region " << regionI << " starts at " << cellI << endl;
00131
00132
00133 fvMeshSubset subsetter(mesh);
00134 subsetter.setLargeCellSubset(cellToRegion, regionI);
00135 const fvMesh& subMesh = subsetter.subMesh();
00136 labelList subCellOrder(bandCompression(subMesh.cellCells()));
00137
00138 const labelList& cellMap = subsetter.cellMap();
00139
00140 forAll(subCellOrder, i)
00141 {
00142 cellOrder[cellI++] = cellMap[subCellOrder[i]];
00143 }
00144 }
00145 Pout<< endl;
00146
00147 return cellOrder;
00148 }
00149
00150
00151
00152
00153 labelList regionFaceOrder
00154 (
00155 const primitiveMesh& mesh,
00156 const labelList& cellOrder,
00157 const labelList& cellToRegion
00158 )
00159 {
00160 Pout<< "Determining face order:" << endl;
00161
00162 labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
00163
00164 labelList oldToNewFace(mesh.nFaces(), -1);
00165
00166 label newFaceI = 0;
00167
00168 label prevRegion = -1;
00169
00170 forAll (cellOrder, newCellI)
00171 {
00172 label oldCellI = cellOrder[newCellI];
00173
00174 if (cellToRegion[oldCellI] != prevRegion)
00175 {
00176 prevRegion = cellToRegion[oldCellI];
00177 Pout<< " region " << prevRegion << " internal faces start at "
00178 << newFaceI << endl;
00179 }
00180
00181 const cell& cFaces = mesh.cells()[oldCellI];
00182
00183 SortableList<label> nbr(cFaces.size());
00184
00185 forAll(cFaces, i)
00186 {
00187 label faceI = cFaces[i];
00188
00189 if (mesh.isInternalFace(faceI))
00190 {
00191
00192 label nbrCellI = reverseCellOrder[mesh.faceNeighbour()[faceI]];
00193 if (nbrCellI == newCellI)
00194 {
00195 nbrCellI = reverseCellOrder[mesh.faceOwner()[faceI]];
00196 }
00197
00198 if (cellToRegion[oldCellI] != cellToRegion[cellOrder[nbrCellI]])
00199 {
00200
00201 nbr[i] = -1;
00202 }
00203 else if (newCellI < nbrCellI)
00204 {
00205
00206 nbr[i] = nbrCellI;
00207 }
00208 else
00209 {
00210
00211 nbr[i] = -1;
00212 }
00213 }
00214 else
00215 {
00216
00217 nbr[i] = -1;
00218 }
00219 }
00220
00221 nbr.sort();
00222
00223 forAll(nbr, i)
00224 {
00225 if (nbr[i] != -1)
00226 {
00227 oldToNewFace[cFaces[nbr.indices()[i]]] = newFaceI++;
00228 }
00229 }
00230 }
00231
00232
00233 label nRegions = max(cellToRegion)+1;
00234 {
00235
00236 SortableList<label> sortKey(mesh.nFaces(), labelMax);
00237
00238 for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
00239 {
00240 label ownRegion = cellToRegion[mesh.faceOwner()[faceI]];
00241 label neiRegion = cellToRegion[mesh.faceNeighbour()[faceI]];
00242
00243 if (ownRegion != neiRegion)
00244 {
00245 sortKey[faceI] =
00246 min(ownRegion, neiRegion)*nRegions
00247 +max(ownRegion, neiRegion);
00248 }
00249 }
00250 sortKey.sort();
00251
00252
00253 label prevKey = -1;
00254 forAll(sortKey, i)
00255 {
00256 label key = sortKey[i];
00257
00258 if (key == labelMax)
00259 {
00260 break;
00261 }
00262
00263 if (prevKey != key)
00264 {
00265 Pout<< " faces inbetween region " << key/nRegions
00266 << " and " << key%nRegions
00267 << " start at " << newFaceI << endl;
00268 prevKey = key;
00269 }
00270
00271 oldToNewFace[sortKey.indices()[i]] = newFaceI++;
00272 }
00273 }
00274
00275
00276 for (label faceI = newFaceI; faceI < mesh.nFaces(); faceI++)
00277 {
00278 oldToNewFace[faceI] = faceI;
00279 }
00280
00281
00282
00283 forAll(oldToNewFace, faceI)
00284 {
00285 if (oldToNewFace[faceI] == -1)
00286 {
00287 FatalErrorIn
00288 (
00289 "polyDualMesh::getFaceOrder"
00290 "(const labelList&, const labelList&, const label) const"
00291 ) << "Did not determine new position"
00292 << " for face " << faceI
00293 << abort(FatalError);
00294 }
00295 }
00296 Pout<< endl;
00297
00298 return invert(mesh.nFaces(), oldToNewFace);
00299 }
00300
00301
00302
00303
00304
00305 autoPtr<mapPolyMesh> reorderMesh
00306 (
00307 polyMesh& mesh,
00308 const labelList& cellOrder,
00309 const labelList& faceOrder
00310 )
00311 {
00312 labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
00313 labelList reverseFaceOrder(invert(faceOrder.size(), faceOrder));
00314
00315 faceList newFaces(reorder(reverseFaceOrder, mesh.faces()));
00316 labelList newOwner
00317 (
00318 renumber
00319 (
00320 reverseCellOrder,
00321 reorder(reverseFaceOrder, mesh.faceOwner())
00322 )
00323 );
00324 labelList newNeighbour
00325 (
00326 renumber
00327 (
00328 reverseCellOrder,
00329 reorder(reverseFaceOrder, mesh.faceNeighbour())
00330 )
00331 );
00332
00333
00334 forAll(newNeighbour, faceI)
00335 {
00336 label own = newOwner[faceI];
00337 label nei = newNeighbour[faceI];
00338
00339 if (nei < own)
00340 {
00341 newFaces[faceI] = newFaces[faceI].reverseFace();
00342 Swap(newOwner[faceI], newNeighbour[faceI]);
00343 }
00344 }
00345
00346 const polyBoundaryMesh& patches = mesh.boundaryMesh();
00347 labelList patchSizes(patches.size());
00348 labelList patchStarts(patches.size());
00349 labelList oldPatchNMeshPoints(patches.size());
00350 labelListList patchPointMap(patches.size());
00351 forAll(patches, patchI)
00352 {
00353 patchSizes[patchI] = patches[patchI].size();
00354 patchStarts[patchI] = patches[patchI].start();
00355 oldPatchNMeshPoints[patchI] = patches[patchI].nPoints();
00356 patchPointMap[patchI] = identity(patches[patchI].nPoints());
00357 }
00358
00359 mesh.resetPrimitives
00360 (
00361 Xfer<pointField>::null(),
00362 xferMove(newFaces),
00363 xferMove(newOwner),
00364 xferMove(newNeighbour),
00365 patchSizes,
00366 patchStarts
00367 );
00368
00369 return autoPtr<mapPolyMesh>
00370 (
00371 new mapPolyMesh
00372 (
00373 mesh,
00374 mesh.nPoints(),
00375 mesh.nFaces(),
00376 mesh.nCells(),
00377 identity(mesh.nPoints()),
00378 List<objectMap>(0),
00379 faceOrder,
00380 List<objectMap>(0),
00381 List<objectMap>(0),
00382 List<objectMap>(0),
00383 cellOrder,
00384 List<objectMap>(0),
00385 List<objectMap>(0),
00386 List<objectMap>(0),
00387 List<objectMap>(0),
00388 identity(mesh.nPoints()),
00389 reverseFaceOrder,
00390 reverseCellOrder,
00391 labelHashSet(0),
00392 patchPointMap,
00393 labelListList(0),
00394 labelListList(0),
00395 labelListList(0),
00396 labelListList(0),
00397 pointField(0),
00398 patchStarts,
00399 oldPatchNMeshPoints
00400 )
00401 );
00402 }
00403
00404
00405
00406
00407 int main(int argc, char *argv[])
00408 {
00409 argList::validOptions.insert("blockOrder", "");
00410 argList::validOptions.insert("orderPoints", "");
00411 argList::validOptions.insert("writeMaps", "");
00412 argList::validOptions.insert("overwrite", "");
00413
00414 # include <OpenFOAM/addRegionOption.H>
00415 # include <OpenFOAM/addTimeOptions.H>
00416
00417 # include <OpenFOAM/setRootCase.H>
00418 # include <OpenFOAM/createTime.H>
00419 runTime.functionObjects().off();
00420
00421
00422 instantList Times = runTime.times();
00423
00424
00425 # include <OpenFOAM/checkTimeOptions.H>
00426
00427 runTime.setTime(Times[startTime], startTime);
00428
00429 # include <OpenFOAM/createNamedMesh.H>
00430 const word oldInstance = mesh.pointsInstance();
00431
00432 const bool blockOrder = args.optionFound("blockOrder");
00433 if (blockOrder)
00434 {
00435 Info<< "Ordering cells into regions (using decomposition);"
00436 << " ordering faces into region-internal and region-external." << nl
00437 << endl;
00438 }
00439
00440 const bool orderPoints = args.optionFound("orderPoints");
00441 if (orderPoints)
00442 {
00443 Info<< "Ordering points into internal and boundary points." << nl
00444 << endl;
00445 }
00446
00447 const bool writeMaps = args.optionFound("writeMaps");
00448
00449 if (writeMaps)
00450 {
00451 Info<< "Writing renumber maps (new to old) to polyMesh." << nl
00452 << endl;
00453 }
00454
00455 bool overwrite = args.optionFound("overwrite");
00456
00457 label band = getBand(mesh.faceOwner(), mesh.faceNeighbour());
00458
00459 Info<< "Mesh size: " << returnReduce(mesh.nCells(), sumOp<label>()) << nl
00460 << "Band before renumbering: "
00461 << returnReduce(band, maxOp<label>()) << nl << endl;
00462
00463
00464
00465 labelIOList cellProcAddressing
00466 (
00467 IOobject
00468 (
00469 "cellProcAddressing",
00470 mesh.facesInstance(),
00471 polyMesh::meshSubDir,
00472 mesh,
00473 IOobject::READ_IF_PRESENT
00474 ),
00475 labelList(0)
00476 );
00477
00478 labelIOList faceProcAddressing
00479 (
00480 IOobject
00481 (
00482 "faceProcAddressing",
00483 mesh.facesInstance(),
00484 polyMesh::meshSubDir,
00485 mesh,
00486 IOobject::READ_IF_PRESENT
00487 ),
00488 labelList(0)
00489 );
00490 labelIOList pointProcAddressing
00491 (
00492 IOobject
00493 (
00494 "pointProcAddressing",
00495 mesh.pointsInstance(),
00496 polyMesh::meshSubDir,
00497 mesh,
00498 IOobject::READ_IF_PRESENT
00499 ),
00500 labelList(0)
00501 );
00502 labelIOList boundaryProcAddressing
00503 (
00504 IOobject
00505 (
00506 "boundaryProcAddressing",
00507 mesh.pointsInstance(),
00508 polyMesh::meshSubDir,
00509 mesh,
00510 IOobject::READ_IF_PRESENT
00511 ),
00512 labelList(0)
00513 );
00514
00515
00516
00517 IOobjectList objects(mesh, runTime.timeName());
00518
00519
00520
00521 PtrList<volScalarField> vsFlds;
00522 ReadFields(mesh, objects, vsFlds);
00523
00524 PtrList<volVectorField> vvFlds;
00525 ReadFields(mesh, objects, vvFlds);
00526
00527 PtrList<volSphericalTensorField> vstFlds;
00528 ReadFields(mesh, objects, vstFlds);
00529
00530 PtrList<volSymmTensorField> vsymtFlds;
00531 ReadFields(mesh, objects, vsymtFlds);
00532
00533 PtrList<volTensorField> vtFlds;
00534 ReadFields(mesh, objects, vtFlds);
00535
00536
00537
00538 PtrList<surfaceScalarField> ssFlds;
00539 ReadFields(mesh, objects, ssFlds);
00540
00541 PtrList<surfaceVectorField> svFlds;
00542 ReadFields(mesh, objects, svFlds);
00543
00544 PtrList<surfaceSphericalTensorField> sstFlds;
00545 ReadFields(mesh, objects, sstFlds);
00546
00547 PtrList<surfaceSymmTensorField> ssymtFlds;
00548 ReadFields(mesh, objects, ssymtFlds);
00549
00550 PtrList<surfaceTensorField> stFlds;
00551 ReadFields(mesh, objects, stFlds);
00552
00553
00554 autoPtr<mapPolyMesh> map;
00555
00556 if (blockOrder)
00557 {
00558
00559
00560
00561
00562 IOdictionary decomposeDict
00563 (
00564 IOobject
00565 (
00566 "decomposeParDict",
00567 runTime.system(),
00568 mesh,
00569 IOobject::MUST_READ,
00570 IOobject::NO_WRITE
00571 )
00572 );
00573 autoPtr<decompositionMethod> decomposePtr = decompositionMethod::New
00574 (
00575 decomposeDict,
00576 mesh
00577 );
00578
00579 labelList cellToRegion(decomposePtr().decompose(mesh.cellCentres()));
00580
00581
00582 {
00583 volScalarField cellDist
00584 (
00585 IOobject
00586 (
00587 "cellDist",
00588 runTime.timeName(),
00589 mesh,
00590 IOobject::NO_READ,
00591 IOobject::NO_WRITE,
00592 false
00593 ),
00594 mesh,
00595 dimensionedScalar("cellDist", dimless, 0),
00596 zeroGradientFvPatchScalarField::typeName
00597 );
00598
00599 forAll(cellToRegion, cellI)
00600 {
00601 cellDist[cellI] = cellToRegion[cellI];
00602 }
00603
00604 cellDist.write();
00605
00606 Info<< nl << "Written decomposition as volScalarField to "
00607 << cellDist.name() << " for use in postprocessing."
00608 << nl << endl;
00609 }
00610
00611
00612
00613 labelList cellOrder(regionBandCompression(mesh, cellToRegion));
00614
00615
00616 labelList faceOrder
00617 (
00618 regionFaceOrder
00619 (
00620 mesh,
00621 cellOrder,
00622 cellToRegion
00623 )
00624 );
00625
00626 if (!overwrite)
00627 {
00628 runTime++;
00629 }
00630
00631
00632 map = reorderMesh(mesh, cellOrder, faceOrder);
00633 }
00634 else
00635 {
00636
00637
00638 polyTopoChange meshMod(mesh);
00639
00640 if (!overwrite)
00641 {
00642 runTime++;
00643 }
00644
00645
00646 map = meshMod.changeMesh
00647 (
00648 mesh,
00649 false,
00650 true,
00651 true,
00652 orderPoints
00653 );
00654 }
00655
00656
00657 mesh.updateMesh(map);
00658
00659
00660 if (cellProcAddressing.headerOk())
00661 {
00662 Info<< "Renumbering processor cell decomposition map "
00663 << cellProcAddressing.name() << endl;
00664
00665 cellProcAddressing = labelList
00666 (
00667 UIndirectList<label>(cellProcAddressing, map().cellMap())
00668 );
00669 }
00670 if (faceProcAddressing.headerOk())
00671 {
00672 Info<< "Renumbering processor face decomposition map "
00673 << faceProcAddressing.name() << endl;
00674
00675 faceProcAddressing = labelList
00676 (
00677 UIndirectList<label>(faceProcAddressing, map().faceMap())
00678 );
00679 }
00680 if (pointProcAddressing.headerOk())
00681 {
00682 Info<< "Renumbering processor point decomposition map "
00683 << pointProcAddressing.name() << endl;
00684
00685 pointProcAddressing = labelList
00686 (
00687 UIndirectList<label>(pointProcAddressing, map().pointMap())
00688 );
00689 }
00690
00691
00692
00693 if (map().hasMotionPoints())
00694 {
00695 mesh.movePoints(map().preMotionPoints());
00696 }
00697
00698
00699 band = getBand(mesh.faceOwner(), mesh.faceNeighbour());
00700
00701 Info<< "Band after renumbering: "
00702 << returnReduce(band, maxOp<label>()) << nl << endl;
00703
00704
00705 if (orderPoints)
00706 {
00707
00708
00709 (void)mesh.edges();
00710
00711 label nTotPoints = returnReduce
00712 (
00713 mesh.nPoints(),
00714 sumOp<label>()
00715 );
00716 label nTotIntPoints = returnReduce
00717 (
00718 mesh.nInternalPoints(),
00719 sumOp<label>()
00720 );
00721
00722 label nTotEdges = returnReduce
00723 (
00724 mesh.nEdges(),
00725 sumOp<label>()
00726 );
00727 label nTotIntEdges = returnReduce
00728 (
00729 mesh.nInternalEdges(),
00730 sumOp<label>()
00731 );
00732 label nTotInt0Edges = returnReduce
00733 (
00734 mesh.nInternal0Edges(),
00735 sumOp<label>()
00736 );
00737 label nTotInt1Edges = returnReduce
00738 (
00739 mesh.nInternal1Edges(),
00740 sumOp<label>()
00741 );
00742
00743 Info<< "Points:" << nl
00744 << " total : " << nTotPoints << nl
00745 << " internal: " << nTotIntPoints << nl
00746 << " boundary: " << nTotPoints-nTotIntPoints << nl
00747 << "Edges:" << nl
00748 << " total : " << nTotEdges << nl
00749 << " internal: " << nTotIntEdges << nl
00750 << " internal using 0 boundary points: "
00751 << nTotInt0Edges << nl
00752 << " internal using 1 boundary points: "
00753 << nTotInt1Edges-nTotInt0Edges << nl
00754 << " internal using 2 boundary points: "
00755 << nTotIntEdges-nTotInt1Edges << nl
00756 << " boundary: " << nTotEdges-nTotIntEdges << nl
00757 << endl;
00758 }
00759
00760
00761 if (overwrite)
00762 {
00763 mesh.setInstance(oldInstance);
00764 }
00765
00766 Info<< "Writing mesh to " << mesh.facesInstance() << endl;
00767
00768 mesh.write();
00769 if (cellProcAddressing.headerOk())
00770 {
00771 cellProcAddressing.instance() = mesh.facesInstance();
00772 cellProcAddressing.write();
00773 }
00774 if (faceProcAddressing.headerOk())
00775 {
00776 faceProcAddressing.instance() = mesh.facesInstance();
00777 faceProcAddressing.write();
00778 }
00779 if (pointProcAddressing.headerOk())
00780 {
00781 pointProcAddressing.instance() = mesh.facesInstance();
00782 pointProcAddressing.write();
00783 }
00784 if (boundaryProcAddressing.headerOk())
00785 {
00786 boundaryProcAddressing.instance() = mesh.facesInstance();
00787 boundaryProcAddressing.write();
00788 }
00789
00790
00791 if (writeMaps)
00792 {
00793 labelIOList
00794 (
00795 IOobject
00796 (
00797 "cellMap",
00798 mesh.facesInstance(),
00799 polyMesh::meshSubDir,
00800 mesh,
00801 IOobject::NO_READ,
00802 IOobject::NO_WRITE,
00803 false
00804 ),
00805 map().cellMap()
00806 ).write();
00807 labelIOList
00808 (
00809 IOobject
00810 (
00811 "faceMap",
00812 mesh.facesInstance(),
00813 polyMesh::meshSubDir,
00814 mesh,
00815 IOobject::NO_READ,
00816 IOobject::NO_WRITE,
00817 false
00818 ),
00819 map().faceMap()
00820 ).write();
00821 labelIOList
00822 (
00823 IOobject
00824 (
00825 "pointMap",
00826 mesh.facesInstance(),
00827 polyMesh::meshSubDir,
00828 mesh,
00829 IOobject::NO_READ,
00830 IOobject::NO_WRITE,
00831 false
00832 ),
00833 map().pointMap()
00834 ).write();
00835 }
00836
00837 Info<< "\nEnd.\n" << endl;
00838
00839 return 0;
00840 }
00841
00842
00843