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

redistributeMeshPar.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 Application
00025     redistributeMeshPar
00026 
00027 Description
00028     Redistributes existing decomposed mesh and fields according to the current
00029     settings in the decomposeParDict file.
00030 
00031     Must be run on maximum number of source and destination processors.
00032     Balances mesh and writes new mesh to new time directory.
00033 
00034     Can also work like decomposePar:
00035 
00036         # Create empty processors
00037         mkdir processor0
00038                 ..
00039         mkdir processorN
00040 
00041         # Copy undecomposed polyMesh
00042         cp -r constant processor0
00043 
00044         # Distribute
00045         mpirun -np ddd redistributeMeshPar -parallel
00046 
00047 Note
00048     You might want to unset FOAM_SIGFPE and FOAM_SETNAN since
00049     patchfields that hold additional data might not be initialised
00050     (since mapped from 0 faces)
00051 
00052 Usage
00053     - redistributeMeshPar [OPTIONS]
00054 
00055     @param -mergeTol <number>\n
00056     Relative merge distance.
00057 
00058     @param -case <dir>\n
00059     Case directory.
00060 
00061     @param -parallel \n
00062     Run in parallel.
00063 
00064     @param -help \n
00065     Display help message.
00066 
00067     @param -doc \n
00068     Display Doxygen API documentation page for this application.
00069 
00070     @param -srcDoc \n
00071     Display Doxygen source documentation page for this application.
00072 
00073 \*---------------------------------------------------------------------------*/
00074 
00075 #include <finiteVolume/fvMesh.H>
00076 #include <decompositionMethods/decompositionMethod.H>
00077 #include <OpenFOAM/PstreamReduceOps.H>
00078 #include <finiteVolume/fvCFD.H>
00079 #include <dynamicMesh/fvMeshDistribute.H>
00080 #include <OpenFOAM/mapDistributePolyMesh.H>
00081 #include <OpenFOAM/IOobjectList.H>
00082 
00083 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00084 
00085 // Tolerance (as fraction of the bounding box). Needs to be fairly lax since
00086 // usually meshes get written with limited precision (6 digits)
00087 static const scalar defaultMergeTol = 1E-6;
00088 
00089 
00090 // Read mesh if available. Otherwise create empty mesh with same non-proc
00091 // patches as proc0 mesh. Requires all processors to have all patches
00092 // (and in same order).
00093 autoPtr<fvMesh> createMesh
00094 (
00095     const Time& runTime,
00096     const word& regionName,
00097     const fileName& instDir,
00098     const bool haveMesh
00099 )
00100 {
00101     Pout<< "Create mesh for time = "
00102         << runTime.timeName() << nl << endl;
00103 
00104     IOobject io
00105     (
00106         regionName,
00107         instDir,
00108         runTime,
00109         IOobject::MUST_READ
00110     );
00111 
00112     if (!haveMesh)
00113     {
00114         // Create dummy mesh. Only used on procs that don't have mesh.
00115         fvMesh dummyMesh
00116         (
00117             io,
00118             xferCopy(pointField()),
00119             xferCopy(faceList()),
00120             xferCopy(labelList()),
00121             xferCopy(labelList()),
00122             false
00123         );
00124         Pout<< "Writing dummy mesh to " << dummyMesh.polyMesh::objectPath()
00125             << endl;
00126         dummyMesh.write();
00127     }
00128 
00129     Pout<< "Reading mesh from " << io.objectPath() << endl;
00130     autoPtr<fvMesh> meshPtr(new fvMesh(io));
00131     fvMesh& mesh = meshPtr();
00132 
00133 
00134     // Determine patches.
00135     if (Pstream::master())
00136     {
00137         // Send patches
00138         for
00139         (
00140             int slave=Pstream::firstSlave();
00141             slave<=Pstream::lastSlave();
00142             slave++
00143         )
00144         {
00145             OPstream toSlave(Pstream::blocking, slave);
00146             toSlave << mesh.boundaryMesh();
00147         }
00148     }
00149     else
00150     {
00151         // Receive patches
00152         IPstream fromMaster(Pstream::blocking, Pstream::masterNo());
00153         PtrList<entry> patchEntries(fromMaster);
00154 
00155         if (haveMesh)
00156         {
00157             // Check master names against mine
00158 
00159             const polyBoundaryMesh& patches = mesh.boundaryMesh();
00160 
00161             forAll(patchEntries, patchI)
00162             {
00163                 const entry& e = patchEntries[patchI];
00164                 const word type(e.dict().lookup("type"));
00165                 const word& name = e.keyword();
00166 
00167                 if (type == processorPolyPatch::typeName)
00168                 {
00169                     break;
00170                 }
00171 
00172                 if (patchI >= patches.size())
00173                 {
00174                     FatalErrorIn
00175                     (
00176                         "createMesh(const Time&, const fileName&, const bool)"
00177                     )   << "Non-processor patches not synchronised."
00178                         << endl
00179                         << "Processor " << Pstream::myProcNo()
00180                         << " has only " << patches.size()
00181                         << " patches, master has "
00182                         << patchI
00183                         << exit(FatalError);
00184                 }
00185 
00186                 if
00187                 (
00188                     type != patches[patchI].type()
00189                  || name != patches[patchI].name()
00190                 )
00191                 {
00192                     FatalErrorIn
00193                     (
00194                         "createMesh(const Time&, const fileName&, const bool)"
00195                     )   << "Non-processor patches not synchronised."
00196                         << endl
00197                         << "Master patch " << patchI
00198                         << " name:" << type
00199                         << " type:" << type << endl
00200                         << "Processor " << Pstream::myProcNo()
00201                         << " patch " << patchI
00202                         << " has name:" << patches[patchI].name()
00203                         << " type:" << patches[patchI].type()
00204                         << exit(FatalError);
00205                 }
00206             }
00207         }
00208         else
00209         {
00210             // Add patch
00211             List<polyPatch*> patches(patchEntries.size());
00212             label nPatches = 0;
00213 
00214             forAll(patchEntries, patchI)
00215             {
00216                 const entry& e = patchEntries[patchI];
00217                 const word type(e.dict().lookup("type"));
00218                 const word& name = e.keyword();
00219 
00220                 if (type == processorPolyPatch::typeName)
00221                 {
00222                     break;
00223                 }
00224 
00225                 Pout<< "Adding patch:" << nPatches
00226                     << " name:" << name
00227                     << " type:" << type << endl;
00228 
00229                 dictionary patchDict(e.dict());
00230                 patchDict.remove("nFaces");
00231                 patchDict.add("nFaces", 0);
00232                 patchDict.remove("startFace");
00233                 patchDict.add("startFace", 0);
00234 
00235                 patches[patchI] = polyPatch::New
00236                 (
00237                     name,
00238                     patchDict,
00239                     nPatches++,
00240                     mesh.boundaryMesh()
00241                 ).ptr();
00242             }
00243             patches.setSize(nPatches);
00244             mesh.addFvPatches(patches, false);  // no parallel comms
00245 
00247             //meshPtr().write();
00248         }
00249     }
00250 
00251     if (!haveMesh)
00252     {
00253         // We created a dummy mesh file above. Delete it.
00254         Pout<< "Removing dummy mesh " << io.objectPath()
00255             << endl;
00256         rmDir(io.objectPath());
00257     }
00258 
00259     // Force recreation of globalMeshData.
00260     mesh.clearOut();
00261     mesh.globalData();
00262 
00263     return meshPtr;
00264 }
00265 
00266 
00267 // Get merging distance when matching face centres
00268 scalar getMergeDistance
00269 (
00270     const argList& args,
00271     const Time& runTime,
00272     const boundBox& bb
00273 )
00274 {
00275     scalar mergeTol = defaultMergeTol;
00276     args.optionReadIfPresent("mergeTol", mergeTol);
00277 
00278     scalar writeTol =
00279         Foam::pow(scalar(10.0), -scalar(IOstream::defaultPrecision()));
00280 
00281     Info<< "Merge tolerance : " << mergeTol << nl
00282         << "Write tolerance : " << writeTol << endl;
00283 
00284     if (runTime.writeFormat() == IOstream::ASCII && mergeTol < writeTol)
00285     {
00286         FatalErrorIn("getMergeDistance")
00287             << "Your current settings specify ASCII writing with "
00288             << IOstream::defaultPrecision() << " digits precision." << endl
00289             << "Your merging tolerance (" << mergeTol << ") is finer than this."
00290             << endl
00291             << "Please change your writeFormat to binary"
00292             << " or increase the writePrecision" << endl
00293             << "or adjust the merge tolerance (-mergeTol)."
00294             << exit(FatalError);
00295     }
00296 
00297     scalar mergeDist = mergeTol * bb.mag();
00298 
00299     Info<< "Overall meshes bounding box : " << bb << nl
00300         << "Relative tolerance          : " << mergeTol << nl
00301         << "Absolute matching distance  : " << mergeDist << nl
00302         << endl;
00303 
00304     return mergeDist;
00305 }
00306 
00307 
00308 void printMeshData(Ostream& os, const polyMesh& mesh)
00309 {
00310     os  << "Number of points:           " << mesh.points().size() << nl
00311         << "          faces:            " << mesh.faces().size() << nl
00312         << "          internal faces:   " << mesh.faceNeighbour().size() << nl
00313         << "          cells:            " << mesh.cells().size() << nl
00314         << "          boundary patches: " << mesh.boundaryMesh().size() << nl
00315         << "          point zones:      " << mesh.pointZones().size() << nl
00316         << "          face zones:       " << mesh.faceZones().size() << nl
00317         << "          cell zones:       " << mesh.cellZones().size() << nl;
00318 }
00319 
00320 
00321 // Debugging: write volScalarField with decomposition for post processing.
00322 void writeDecomposition
00323 (
00324     const word& name,
00325     const fvMesh& mesh,
00326     const labelList& decomp
00327 )
00328 {
00329     Info<< "Writing wanted cell distribution to volScalarField " << name
00330         << " for postprocessing purposes." << nl << endl;
00331 
00332     volScalarField procCells
00333     (
00334         IOobject
00335         (
00336             name,
00337             mesh.time().timeName(),
00338             mesh,
00339             IOobject::NO_READ,
00340             IOobject::AUTO_WRITE,
00341             false                   // do not register
00342         ),
00343         mesh,
00344         dimensionedScalar(name, dimless, -1),
00345         zeroGradientFvPatchScalarField::typeName
00346     );
00347 
00348     forAll(procCells, cI)
00349     {
00350         procCells[cI] = decomp[cI];
00351     }
00352     procCells.write();
00353 }
00354 
00355 
00356 // Read vol or surface fields
00357 //template<class T, class Mesh>
00358 template<class GeoField>
00359 void readFields
00360 (
00361     const boolList& haveMesh,
00362     const fvMesh& mesh,
00363     const autoPtr<fvMeshSubset>& subsetterPtr,
00364     IOobjectList& allObjects,
00365     PtrList<GeoField>& fields
00366 )
00367 {
00368     //typedef GeometricField<T, fvPatchField, Mesh> fldType;
00369 
00370     // Get my objects of type
00371     IOobjectList objects(allObjects.lookupClass(GeoField::typeName));
00372     // Check that we all have all objects
00373     wordList objectNames = objects.toc();
00374     // Get master names
00375     wordList masterNames(objectNames);
00376     Pstream::scatter(masterNames);
00377 
00378     if (haveMesh[Pstream::myProcNo()] && objectNames != masterNames)
00379     {
00380         FatalErrorIn("readFields()")
00381             << "differing fields of type " << GeoField::typeName
00382             << " on processors." << endl
00383             << "Master has:" << masterNames << endl
00384             << Pstream::myProcNo() << " has:" << objectNames
00385             << abort(FatalError);
00386     }
00387 
00388     fields.setSize(masterNames.size());
00389 
00390     // Have master send all fields to processors that don't have a mesh
00391     if (Pstream::master())
00392     {
00393         forAll(masterNames, i)
00394         {
00395             const word& name = masterNames[i];
00396             IOobject& io = *objects[name];
00397             io.writeOpt() = IOobject::AUTO_WRITE;
00398 
00399             // Load field
00400             fields.set(i, new GeoField(io, mesh));
00401 
00402             // Create zero sized field and send
00403             if (subsetterPtr.valid())
00404             {
00405                 tmp<GeoField> tsubfld = subsetterPtr().interpolate(fields[i]);
00406 
00407                 // Send to all processors that don't have a mesh
00408                 for (label procI = 1; procI < Pstream::nProcs(); procI++)
00409                 {
00410                     if (!haveMesh[procI])
00411                     {
00412                         OPstream toProc(Pstream::blocking, procI);
00413                         toProc<< tsubfld();
00414                     }
00415                 }
00416             }
00417         }
00418     }
00419     else if (!haveMesh[Pstream::myProcNo()])
00420     {
00421         // Don't have mesh (nor fields). Receive empty field from master.
00422 
00423         forAll(masterNames, i)
00424         {
00425             const word& name = masterNames[i];
00426 
00427             // Receive field
00428             IPstream fromMaster(Pstream::blocking, Pstream::masterNo());
00429 
00430             fields.set
00431             (
00432                 i,
00433                 new GeoField
00434                 (
00435                     IOobject
00436                     (
00437                         name,
00438                         mesh.time().timeName(),
00439                         mesh,
00440                         IOobject::NO_READ,
00441                         IOobject::AUTO_WRITE
00442                     ),
00443                     mesh,
00444                     fromMaster
00445                 )
00446             );
00447 
00449             //fields[i].write();
00450         }
00451     }
00452     else
00453     {
00454         // Have mesh so just try to load
00455         forAll(masterNames, i)
00456         {
00457             const word& name = masterNames[i];
00458             IOobject& io = *objects[name];
00459             io.writeOpt() = IOobject::AUTO_WRITE;
00460 
00461             // Load field
00462             fields.set(i, new GeoField(io, mesh));
00463         }
00464     }
00465 }
00466 
00467 
00468 // Debugging: compare two fields.
00469 void compareFields
00470 (
00471     const scalar tolDim,
00472     const volVectorField& a,
00473     const volVectorField& b
00474 )
00475 {
00476     forAll(a, cellI)
00477     {
00478         if (mag(b[cellI] - a[cellI]) > tolDim)
00479         {
00480             FatalErrorIn
00481             (
00482                 "compareFields"
00483                 "(const scalar, const volVectorField&, const volVectorField&)"
00484             )   << "Did not map volVectorField correctly:" << nl
00485                 << "cell:" << cellI
00486                 << " transfer b:" << b[cellI]
00487                 << " real cc:" << a[cellI]
00488                 << abort(FatalError);
00489         }
00490     }
00491     forAll(a.boundaryField(), patchI)
00492     {
00493         // We have real mesh cellcentre and
00494         // mapped original cell centre.
00495 
00496         const fvPatchVectorField& aBoundary =
00497             a.boundaryField()[patchI];
00498 
00499         const fvPatchVectorField& bBoundary =
00500             b.boundaryField()[patchI];
00501 
00502         if (!bBoundary.coupled())
00503         {
00504             forAll(aBoundary, i)
00505             {
00506                 if (mag(aBoundary[i] - bBoundary[i]) > tolDim)
00507                 {
00508                     FatalErrorIn
00509                     (
00510                         "compareFields"
00511                         "(const scalar, const volVectorField&"
00512                         ", const volVectorField&)"
00513                     )   << "Did not map volVectorField correctly:"
00514                         << endl
00515                         << "patch:" << patchI << " patchFace:" << i
00516                         << " cc:" << endl
00517                         << "    real    :" << aBoundary[i] << endl
00518                         << "    mapped  :" << bBoundary[i] << endl
00519                         << abort(FatalError);
00520                 }
00521             }
00522         }
00523     }
00524 }
00525 
00526 
00527 // Main program:
00528 
00529 int main(int argc, char *argv[])
00530 {
00531 #   include <OpenFOAM/addRegionOption.H>
00532     argList::validOptions.insert("mergeTol", "relative merge distance");
00533 
00534 #   include <OpenFOAM/setRootCase.H>
00535 
00536     // Create processor directory if non-existing
00537     if (!Pstream::master() && !isDir(args.path()))
00538     {
00539         Pout<< "Creating case directory " << args.path() << endl;
00540         mkDir(args.path());
00541     }
00542 
00543 #   include <OpenFOAM/createTime.H>
00544 
00545     word regionName = polyMesh::defaultRegion;
00546     fileName meshSubDir;
00547 
00548     if (args.optionReadIfPresent("region", regionName))
00549     {
00550         meshSubDir = regionName/polyMesh::meshSubDir;
00551     }
00552     else
00553     {
00554         meshSubDir = polyMesh::meshSubDir;
00555     }
00556     Info<< "Using mesh subdirectory " << meshSubDir << nl << endl;
00557 
00558 
00559     // Get time instance directory. Since not all processors have meshes
00560     // just use the master one everywhere.
00561 
00562     fileName masterInstDir;
00563     if (Pstream::master())
00564     {
00565         masterInstDir = runTime.findInstance(meshSubDir, "points");
00566     }
00567     Pstream::scatter(masterInstDir);
00568 
00569     // Check who has a mesh
00570     const fileName meshPath = runTime.path()/masterInstDir/meshSubDir;
00571 
00572     Info<< "Found points in " << meshPath << nl << endl;
00573 
00574 
00575     boolList haveMesh(Pstream::nProcs(), false);
00576     haveMesh[Pstream::myProcNo()] = isDir(meshPath);
00577     Pstream::gatherList(haveMesh);
00578     Pstream::scatterList(haveMesh);
00579     Info<< "Per processor mesh availability : " << haveMesh << endl;
00580     const bool allHaveMesh = (findIndex(haveMesh, false) == -1);
00581 
00582     // Create mesh
00583     autoPtr<fvMesh> meshPtr = createMesh
00584     (
00585         runTime,
00586         regionName,
00587         masterInstDir,
00588         haveMesh[Pstream::myProcNo()]
00589     );
00590     fvMesh& mesh = meshPtr();
00591 
00592     Pout<< "Read mesh:" << endl;
00593     printMeshData(Pout, mesh);
00594     Pout<< endl;
00595 
00596     IOdictionary decompositionDict
00597     (
00598         IOobject
00599         (
00600             "decomposeParDict",
00601             runTime.system(),
00602             mesh,
00603             IOobject::MUST_READ,
00604             IOobject::NO_WRITE
00605         )
00606     );
00607 
00608     labelList finalDecomp;
00609 
00610     // Create decompositionMethod and new decomposition
00611     {
00612         autoPtr<decompositionMethod> decomposer
00613         (
00614             decompositionMethod::New
00615             (
00616                 decompositionDict,
00617                 mesh
00618             )
00619         );
00620 
00621         if (!decomposer().parallelAware())
00622         {
00623             WarningIn(args.executable())
00624                 << "You have selected decomposition method "
00625                 << decomposer().typeName
00626                 << " which does" << endl
00627                 << "not synchronise the decomposition across"
00628                 << " processor patches." << endl
00629                 << "    You might want to select a decomposition method which"
00630                 << " is aware of this. Continuing."
00631                 << endl;
00632         }
00633 
00634         finalDecomp = decomposer().decompose(mesh.cellCentres());
00635     }
00636 
00637     // Dump decomposition to volScalarField
00638     writeDecomposition("decomposition", mesh, finalDecomp);
00639 
00640 
00641     // Create 0 sized mesh to do all the generation of zero sized
00642     // fields on processors that have zero sized meshes. Note that this is
00643     // only nessecary on master but since polyMesh construction with
00644     // Pstream::parRun does parallel comms we have to do it on all
00645     // processors
00646     autoPtr<fvMeshSubset> subsetterPtr;
00647 
00648     if (!allHaveMesh)
00649     {
00650         // Find last non-processor patch.
00651         const polyBoundaryMesh& patches = mesh.boundaryMesh();
00652 
00653         label nonProcI = -1;
00654 
00655         forAll(patches, patchI)
00656         {
00657             if (isA<processorPolyPatch>(patches[patchI]))
00658             {
00659                 break;
00660             }
00661             nonProcI++;
00662         }
00663 
00664         if (nonProcI == -1)
00665         {
00666             FatalErrorIn(args.executable())
00667                 << "Cannot find non-processor patch on processor "
00668                 << Pstream::myProcNo() << endl
00669                 << " Current patches:" << patches.names() << abort(FatalError);
00670         }
00671 
00672         // Subset 0 cells, no parallel comms. This is used to create zero-sized
00673         // fields.
00674         subsetterPtr.reset(new fvMeshSubset(mesh));
00675         subsetterPtr().setLargeCellSubset(labelHashSet(0), nonProcI, false);
00676     }
00677 
00678 
00679     // Get original objects (before incrementing time!)
00680     IOobjectList objects(mesh, runTime.timeName());
00681     // We don't want to map the decomposition (mapping already tested when
00682     // mapping the cell centre field)
00683     IOobjectList::iterator iter = objects.find("decomposition");
00684     if (iter != objects.end())
00685     {
00686         objects.erase(iter);
00687     }
00688 
00689     PtrList<volScalarField> volScalarFields;
00690     readFields
00691     (
00692         haveMesh,
00693         mesh,
00694         subsetterPtr,
00695         objects,
00696         volScalarFields
00697     );
00698 
00699     PtrList<volVectorField> volVectorFields;
00700     readFields
00701     (
00702         haveMesh,
00703         mesh,
00704         subsetterPtr,
00705         objects,
00706         volVectorFields
00707     );
00708 
00709     PtrList<volSphericalTensorField> volSphereTensorFields;
00710     readFields
00711     (
00712         haveMesh,
00713         mesh,
00714         subsetterPtr,
00715         objects,
00716         volSphereTensorFields
00717     );
00718 
00719     PtrList<volSymmTensorField> volSymmTensorFields;
00720     readFields
00721     (
00722         haveMesh,
00723         mesh,
00724         subsetterPtr,
00725         objects,
00726         volSymmTensorFields
00727     );
00728 
00729     PtrList<volTensorField> volTensorFields;
00730     readFields
00731     (
00732         haveMesh,
00733         mesh,
00734         subsetterPtr,
00735         objects,
00736         volTensorFields
00737     );
00738 
00739     PtrList<surfaceScalarField> surfScalarFields;
00740     readFields
00741     (
00742         haveMesh,
00743         mesh,
00744         subsetterPtr,
00745         objects,
00746         surfScalarFields
00747     );
00748 
00749     PtrList<surfaceVectorField> surfVectorFields;
00750     readFields
00751     (
00752         haveMesh,
00753         mesh,
00754         subsetterPtr,
00755         objects,
00756         surfVectorFields
00757     );
00758 
00759     PtrList<surfaceSphericalTensorField> surfSphereTensorFields;
00760     readFields
00761     (
00762         haveMesh,
00763         mesh,
00764         subsetterPtr,
00765         objects,
00766         surfSphereTensorFields
00767     );
00768 
00769     PtrList<surfaceSymmTensorField> surfSymmTensorFields;
00770     readFields
00771     (
00772         haveMesh,
00773         mesh,
00774         subsetterPtr,
00775         objects,
00776         surfSymmTensorFields
00777     );
00778 
00779     PtrList<surfaceTensorField> surfTensorFields;
00780     readFields
00781     (
00782         haveMesh,
00783         mesh,
00784         subsetterPtr,
00785         objects,
00786         surfTensorFields
00787     );
00788 
00789 
00790     // Debugging: Create additional volField that will be mapped.
00791     // Used to test correctness of mapping
00792     volVectorField mapCc("mapCc", 1*mesh.C());
00793 
00794     // Global matching tolerance
00795     const scalar tolDim = getMergeDistance
00796     (
00797         args,
00798         runTime,
00799         mesh.globalData().bb()
00800     );
00801 
00802     // Mesh distribution engine
00803     fvMeshDistribute distributor(mesh, tolDim);
00804 
00805     Pout<< "Wanted distribution:"
00806         << distributor.countCells(finalDecomp) << nl << endl;
00807 
00808     // Do actual sending/receiving of mesh
00809     autoPtr<mapDistributePolyMesh> map = distributor.distribute(finalDecomp);
00810 
00812     //map().distributeFaceData(faceCc);
00813 
00814 
00815     // Print a bit
00816     Pout<< "After distribution mesh:" << endl;
00817     printMeshData(Pout, mesh);
00818     Pout<< endl;
00819 
00820     runTime++;
00821     Pout<< "Writing redistributed mesh to " << runTime.timeName() << nl << endl;
00822     mesh.write();
00823 
00824 
00825     // Debugging: test mapped cellcentre field.
00826     compareFields(tolDim, mesh.C(), mapCc);
00827 
00828     // Print nice message
00829     // ~~~~~~~~~~~~~~~~~~
00830 
00831     // Determine which processors remain so we can print nice final message.
00832     labelList nFaces(Pstream::nProcs());
00833     nFaces[Pstream::myProcNo()] = mesh.nFaces();
00834     Pstream::gatherList(nFaces);
00835     Pstream::scatterList(nFaces);
00836 
00837     Info<< nl
00838         << "You can pick up the redecomposed mesh from the polyMesh directory"
00839         << " in " << runTime.timeName() << "." << nl
00840         << "If you redecomposed the mesh to less processors you can delete"
00841         << nl
00842         << "the processor directories with 0 sized meshes in them." << nl
00843         << "Below is a sample set of commands to do this."
00844         << " Take care when issuing these" << nl
00845         << "commands." << nl << endl;
00846 
00847     forAll(nFaces, procI)
00848     {
00849         fileName procDir = "processor" + name(procI);
00850 
00851         if (nFaces[procI] == 0)
00852         {
00853             Info<< "    rm -r " << procDir.c_str() << nl;
00854         }
00855         else
00856         {
00857             fileName timeDir = procDir/runTime.timeName()/meshSubDir;
00858             fileName constDir = procDir/runTime.constant()/meshSubDir;
00859 
00860             Info<< "    rm -r " << constDir.c_str() << nl
00861                 << "    mv " << timeDir.c_str()
00862                 << ' ' << constDir.c_str() << nl;
00863         }
00864     }
00865     Info<< endl;
00866 
00867 
00868     Pout<< "End\n" << endl;
00869 
00870     return 0;
00871 }
00872 
00873 
00874 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines