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

splitMeshRegions.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     splitMeshRegions
00026 
00027 Description
00028     Splits mesh into multiple regions.
00029 
00030     Each region is defined as a domain whose cells can all be reached by
00031     cell-face-cell walking without crossing
00032     - boundary faces
00033     - additional faces from faceset (-blockedFaces faceSet).
00034     - any face inbetween differing cellZones (-cellZones)
00035 
00036     Output is:
00037     - volScalarField with regions as different scalars (-detectOnly) or
00038     - mesh with multiple regions or
00039     - mesh with cells put into cellZones (-makeCellZones)
00040 
00041 Note
00042     - cellZonesOnly does not do a walk and uses the cellZones only. Use
00043     this if you don't mind having disconnected domains in a single region.
00044     This option requires all cells to be in one (and one only) cellZone.
00045 
00046     - cellZonesFileOnly behaves like -cellZonesOnly but reads the cellZones
00047     from the specified file. This allows one to explicitly specify the region
00048     distribution and still have multiple cellZones per region.
00049 
00050     - useCellZonesOnly does not do a walk and uses the cellZones only. Use
00051     this if you don't mind having disconnected domains in a single region.
00052     This option requires all cells to be in one (and one only) cellZone.
00053 
00054 
00055     - Should work in parallel.
00056     cellZones can differ on either side of processor boundaries in which case
00057     the faces get moved from processor patch to directMapped patch. Not
00058     very well tested.
00059 
00060     - If a cell zone gets split into more than one region it can detect
00061     the largest matching region (-sloppyCellZones). This will accept any
00062     region that covers more than 50% of the zone. It has to be a subset
00063     so cannot have any cells in any other zone.
00064 
00065     - writes maps like decomposePar back to original mesh:
00066         - pointRegionAddressing : for every point in this region the point in
00067         the original mesh
00068         - cellRegionAddressing  :   ,,      cell                ,,  cell    ,,
00069         - faceRegionAddressing  :   ,,      face                ,,  face in
00070         the original mesh + 'turning index'. For a face in the same orientation
00071         this is the original facelabel+1, for a turned face this is -facelabel-1
00072 
00073 Usage
00074 
00075     - splitMeshRegions [OPTIONS]
00076 
00077     @param -cellZones \n
00078     Split different cell zones.
00079 
00080     @param -detectOnly \n
00081     Do no processing.
00082 
00083     @param -blockedFaces <faceSet name>\n
00084     Split at give face set.
00085 
00086     @param -sloppyCellZones \n
00087     Try to match regions to existing cell zones.
00088 
00089     @param -makeCellZones \n
00090     Create mesh with cells in different cell zones.
00091 
00092     @param -overwrite \n
00093     Overwrite existing data.
00094 
00095     @param -case <dir>\n
00096     Case directory.
00097 
00098     @param -parallel \n
00099     Run in parallel.
00100 
00101     @param -help \n
00102     Display help message.
00103 
00104     @param -doc \n
00105     Display Doxygen API documentation page for this application.
00106 
00107     @param -srcDoc \n
00108     Display Doxygen source documentation page for this application.
00109 
00110     @param -insidePoint <point>\n
00111     Only keep region containing specified point.
00112 
00113     @param -largestOnly \n
00114     Only keep largest region.
00115 
00116 \*---------------------------------------------------------------------------*/
00117 
00118 #include <OpenFOAM/SortableList.H>
00119 #include <OpenFOAM/argList.H>
00120 #include <meshTools/regionSplit.H>
00121 #include <finiteVolume/fvMeshSubset.H>
00122 #include <OpenFOAM/IOobjectList.H>
00123 #include <finiteVolume/volFields.H>
00124 #include <meshTools/faceSet.H>
00125 #include <meshTools/cellSet.H>
00126 #include <dynamicMesh/polyTopoChange.H>
00127 #include <dynamicMesh/removeCells.H>
00128 #include <OpenFOAM/EdgeMap.H>
00129 #include <OpenFOAM/syncTools.H>
00130 #include <OpenFOAM/ReadFields.H>
00131 #include <meshTools/directMappedWallPolyPatch.H>
00132 #include <finiteVolume/zeroGradientFvPatchFields.H>
00133 
00134 using namespace Foam;
00135 
00136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00137 
00138 template<class GeoField>
00139 void addPatchFields(fvMesh& mesh, const word& patchFieldType)
00140 {
00141     HashTable<const GeoField*> flds
00142     (
00143         mesh.objectRegistry::lookupClass<GeoField>()
00144     );
00145 
00146     for
00147     (
00148         typename HashTable<const GeoField*>::const_iterator iter = flds.begin();
00149         iter != flds.end();
00150         ++iter
00151     )
00152     {
00153         const GeoField& fld = *iter();
00154 
00155         typename GeoField::GeometricBoundaryField& bfld =
00156             const_cast<typename GeoField::GeometricBoundaryField&>
00157             (
00158                 fld.boundaryField()
00159             );
00160 
00161         label sz = bfld.size();
00162         bfld.setSize(sz+1);
00163         bfld.set
00164         (
00165             sz,
00166             GeoField::PatchFieldType::New
00167             (
00168                 patchFieldType,
00169                 mesh.boundary()[sz],
00170                 fld.dimensionedInternalField()
00171             )
00172         );
00173     }
00174 }
00175 
00176 
00177 // Remove last patch field
00178 template<class GeoField>
00179 void trimPatchFields(fvMesh& mesh, const label nPatches)
00180 {
00181     HashTable<const GeoField*> flds
00182     (
00183         mesh.objectRegistry::lookupClass<GeoField>()
00184     );
00185 
00186     for
00187     (
00188         typename HashTable<const GeoField*>::const_iterator iter = flds.begin();
00189         iter != flds.end();
00190         ++iter
00191     )
00192     {
00193         const GeoField& fld = *iter();
00194 
00195         const_cast<typename GeoField::GeometricBoundaryField&>
00196         (
00197             fld.boundaryField()
00198         ).setSize(nPatches);
00199     }
00200 }
00201 
00202 
00203 // Reorder patch field
00204 template<class GeoField>
00205 void reorderPatchFields(fvMesh& mesh, const labelList& oldToNew)
00206 {
00207     HashTable<const GeoField*> flds
00208     (
00209         mesh.objectRegistry::lookupClass<GeoField>()
00210     );
00211 
00212     for
00213     (
00214         typename HashTable<const GeoField*>::const_iterator iter = flds.begin();
00215         iter != flds.end();
00216         ++iter
00217     )
00218     {
00219         const GeoField& fld = *iter();
00220 
00221         typename GeoField::GeometricBoundaryField& bfld =
00222             const_cast<typename GeoField::GeometricBoundaryField&>
00223             (
00224                 fld.boundaryField()
00225             );
00226 
00227         bfld.reorder(oldToNew);
00228     }
00229 }
00230 
00231 
00232 // Adds patch if not yet there. Returns patchID.
00233 label addPatch(fvMesh& mesh, const polyPatch& patch)
00234 {
00235     polyBoundaryMesh& polyPatches =
00236         const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
00237 
00238     label patchI = polyPatches.findPatchID(patch.name());
00239     if (patchI != -1)
00240     {
00241         if (polyPatches[patchI].type() == patch.type())
00242         {
00243             // Already there
00244             return patchI;
00245         }
00246         else
00247         {
00248             FatalErrorIn("addPatch(fvMesh&, const polyPatch*)")
00249                 << "Already have patch " << patch.name()
00250                 << " but of type " << patch.type()
00251                 << exit(FatalError);
00252         }
00253     }
00254 
00255 
00256     label insertPatchI = polyPatches.size();
00257     label startFaceI = mesh.nFaces();
00258 
00259     forAll(polyPatches, patchI)
00260     {
00261         const polyPatch& pp = polyPatches[patchI];
00262 
00263         if (isA<processorPolyPatch>(pp))
00264         {
00265             insertPatchI = patchI;
00266             startFaceI = pp.start();
00267             break;
00268         }
00269     }
00270 
00271 
00272     // Below is all quite a hack. Feel free to change once there is a better
00273     // mechanism to insert and reorder patches.
00274 
00275     // Clear local fields and e.g. polyMesh parallelInfo.
00276     mesh.clearOut();
00277 
00278     label sz = polyPatches.size();
00279 
00280     fvBoundaryMesh& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
00281 
00282     // Add polyPatch at the end
00283     polyPatches.setSize(sz+1);
00284     polyPatches.set
00285     (
00286         sz,
00287         patch.clone
00288         (
00289             polyPatches,
00290             insertPatchI,   //index
00291             0,              //size
00292             startFaceI      //start
00293         )
00294     );
00295     fvPatches.setSize(sz+1);
00296     fvPatches.set
00297     (
00298         sz,
00299         fvPatch::New
00300         (
00301             polyPatches[sz],  // point to newly added polyPatch
00302             mesh.boundary()
00303         )
00304     );
00305 
00306     addPatchFields<volScalarField>
00307     (
00308         mesh,
00309         calculatedFvPatchField<scalar>::typeName
00310     );
00311     addPatchFields<volVectorField>
00312     (
00313         mesh,
00314         calculatedFvPatchField<vector>::typeName
00315     );
00316     addPatchFields<volSphericalTensorField>
00317     (
00318         mesh,
00319         calculatedFvPatchField<sphericalTensor>::typeName
00320     );
00321     addPatchFields<volSymmTensorField>
00322     (
00323         mesh,
00324         calculatedFvPatchField<symmTensor>::typeName
00325     );
00326     addPatchFields<volTensorField>
00327     (
00328         mesh,
00329         calculatedFvPatchField<tensor>::typeName
00330     );
00331 
00332     // Surface fields
00333 
00334     addPatchFields<surfaceScalarField>
00335     (
00336         mesh,
00337         calculatedFvPatchField<scalar>::typeName
00338     );
00339     addPatchFields<surfaceVectorField>
00340     (
00341         mesh,
00342         calculatedFvPatchField<vector>::typeName
00343     );
00344     addPatchFields<surfaceSphericalTensorField>
00345     (
00346         mesh,
00347         calculatedFvPatchField<sphericalTensor>::typeName
00348     );
00349     addPatchFields<surfaceSymmTensorField>
00350     (
00351         mesh,
00352         calculatedFvPatchField<symmTensor>::typeName
00353     );
00354     addPatchFields<surfaceTensorField>
00355     (
00356         mesh,
00357         calculatedFvPatchField<tensor>::typeName
00358     );
00359 
00360     // Create reordering list
00361     // patches before insert position stay as is
00362     labelList oldToNew(sz+1);
00363     for (label i = 0; i < insertPatchI; i++)
00364     {
00365         oldToNew[i] = i;
00366     }
00367     // patches after insert position move one up
00368     for (label i = insertPatchI; i < sz; i++)
00369     {
00370         oldToNew[i] = i+1;
00371     }
00372     // appended patch gets moved to insert position
00373     oldToNew[sz] = insertPatchI;
00374 
00375     // Shuffle into place
00376     polyPatches.reorder(oldToNew);
00377     fvPatches.reorder(oldToNew);
00378 
00379     reorderPatchFields<volScalarField>(mesh, oldToNew);
00380     reorderPatchFields<volVectorField>(mesh, oldToNew);
00381     reorderPatchFields<volSphericalTensorField>(mesh, oldToNew);
00382     reorderPatchFields<volSymmTensorField>(mesh, oldToNew);
00383     reorderPatchFields<volTensorField>(mesh, oldToNew);
00384     reorderPatchFields<surfaceScalarField>(mesh, oldToNew);
00385     reorderPatchFields<surfaceVectorField>(mesh, oldToNew);
00386     reorderPatchFields<surfaceSphericalTensorField>(mesh, oldToNew);
00387     reorderPatchFields<surfaceSymmTensorField>(mesh, oldToNew);
00388     reorderPatchFields<surfaceTensorField>(mesh, oldToNew);
00389 
00390     return insertPatchI;
00391 }
00392 
00393 
00394 // Reorder and delete patches.
00395 void reorderPatches
00396 (
00397     fvMesh& mesh,
00398     const labelList& oldToNew,
00399     const label nNewPatches
00400 )
00401 {
00402     polyBoundaryMesh& polyPatches =
00403         const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
00404     fvBoundaryMesh& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
00405 
00406     // Shuffle into place
00407     polyPatches.reorder(oldToNew);
00408     fvPatches.reorder(oldToNew);
00409 
00410     reorderPatchFields<volScalarField>(mesh, oldToNew);
00411     reorderPatchFields<volVectorField>(mesh, oldToNew);
00412     reorderPatchFields<volSphericalTensorField>(mesh, oldToNew);
00413     reorderPatchFields<volSymmTensorField>(mesh, oldToNew);
00414     reorderPatchFields<volTensorField>(mesh, oldToNew);
00415     reorderPatchFields<surfaceScalarField>(mesh, oldToNew);
00416     reorderPatchFields<surfaceVectorField>(mesh, oldToNew);
00417     reorderPatchFields<surfaceSphericalTensorField>(mesh, oldToNew);
00418     reorderPatchFields<surfaceSymmTensorField>(mesh, oldToNew);
00419     reorderPatchFields<surfaceTensorField>(mesh, oldToNew);
00420 
00421     // Remove last.
00422     polyPatches.setSize(nNewPatches);
00423     fvPatches.setSize(nNewPatches);
00424     trimPatchFields<volScalarField>(mesh, nNewPatches);
00425     trimPatchFields<volVectorField>(mesh, nNewPatches);
00426     trimPatchFields<volSphericalTensorField>(mesh, nNewPatches);
00427     trimPatchFields<volSymmTensorField>(mesh, nNewPatches);
00428     trimPatchFields<volTensorField>(mesh, nNewPatches);
00429     trimPatchFields<surfaceScalarField>(mesh, nNewPatches);
00430     trimPatchFields<surfaceVectorField>(mesh, nNewPatches);
00431     trimPatchFields<surfaceSphericalTensorField>(mesh, nNewPatches);
00432     trimPatchFields<surfaceSymmTensorField>(mesh, nNewPatches);
00433     trimPatchFields<surfaceTensorField>(mesh, nNewPatches);
00434 }
00435 
00436 
00437 template<class GeoField>
00438 void subsetVolFields
00439 (
00440     const fvMesh& mesh,
00441     const fvMesh& subMesh,
00442     const labelList& cellMap,
00443     const labelList& faceMap,
00444     const labelHashSet& addedPatches
00445 )
00446 {
00447     const labelList patchMap(identity(mesh.boundaryMesh().size()));
00448 
00449     HashTable<const GeoField*> fields
00450     (
00451         mesh.objectRegistry::lookupClass<GeoField>()
00452     );
00453     forAllConstIter(typename HashTable<const GeoField*>, fields, iter)
00454     {
00455         const GeoField& fld = *iter();
00456 
00457         Info<< "Mapping field " << fld.name() << endl;
00458 
00459         tmp<GeoField> tSubFld
00460         (
00461             fvMeshSubset::interpolate
00462             (
00463                 fld,
00464                 subMesh,
00465                 patchMap,
00466                 cellMap,
00467                 faceMap
00468             )
00469         );
00470 
00471         // Hack: set value to 0 for introduced patches (since don't
00472         //       get initialised.
00473         forAll(tSubFld().boundaryField(), patchI)
00474         {
00475             if (addedPatches.found(patchI))
00476             {
00477                 tSubFld().boundaryField()[patchI] ==
00478                     pTraits<typename GeoField::value_type>::zero;
00479             }
00480         }
00481 
00482         // Store on subMesh
00483         GeoField* subFld = tSubFld.ptr();
00484         subFld->rename(fld.name());
00485         subFld->writeOpt() = IOobject::AUTO_WRITE;
00486         subFld->store();
00487     }
00488 }
00489 
00490 
00491 template<class GeoField>
00492 void subsetSurfaceFields
00493 (
00494     const fvMesh& mesh,
00495     const fvMesh& subMesh,
00496     const labelList& faceMap,
00497     const labelHashSet& addedPatches
00498 )
00499 {
00500     const labelList patchMap(identity(mesh.boundaryMesh().size()));
00501 
00502     HashTable<const GeoField*> fields
00503     (
00504         mesh.objectRegistry::lookupClass<GeoField>()
00505     );
00506     forAllConstIter(typename HashTable<const GeoField*>, fields, iter)
00507     {
00508         const GeoField& fld = *iter();
00509 
00510         Info<< "Mapping field " << fld.name() << endl;
00511 
00512         tmp<GeoField> tSubFld
00513         (
00514             fvMeshSubset::interpolate
00515             (
00516                 fld,
00517                 subMesh,
00518                 patchMap,
00519                 faceMap
00520             )
00521         );
00522 
00523         // Hack: set value to 0 for introduced patches (since don't
00524         //       get initialised.
00525         forAll(tSubFld().boundaryField(), patchI)
00526         {
00527             if (addedPatches.found(patchI))
00528             {
00529                 tSubFld().boundaryField()[patchI] ==
00530                     pTraits<typename GeoField::value_type>::zero;
00531             }
00532         }
00533 
00534         // Store on subMesh
00535         GeoField* subFld = tSubFld.ptr();
00536         subFld->rename(fld.name());
00537         subFld->writeOpt() = IOobject::AUTO_WRITE;
00538         subFld->store();
00539     }
00540 }
00541 
00542 // Select all cells not in the region
00543 labelList getNonRegionCells(const labelList& cellRegion, const label regionI)
00544 {
00545     DynamicList<label> nonRegionCells(cellRegion.size());
00546     forAll(cellRegion, cellI)
00547     {
00548         if (cellRegion[cellI] != regionI)
00549         {
00550             nonRegionCells.append(cellI);
00551         }
00552     }
00553     return nonRegionCells.shrink();
00554 }
00555 
00556 
00557 // Get per region-region interface the sizes. If sumParallel sums sizes.
00558 // Returns interfaces as straight list for looping in identical order.
00559 void getInterfaceSizes
00560 (
00561     const polyMesh& mesh,
00562     const labelList& cellRegion,
00563     const bool sumParallel,
00564 
00565     edgeList& interfaces,
00566     EdgeMap<label>& interfaceSizes
00567 )
00568 {
00569 
00570     // Internal faces
00571     // ~~~~~~~~~~~~~~
00572 
00573     forAll(mesh.faceNeighbour(), faceI)
00574     {
00575         label ownRegion = cellRegion[mesh.faceOwner()[faceI]];
00576         label neiRegion = cellRegion[mesh.faceNeighbour()[faceI]];
00577 
00578         if (ownRegion != neiRegion)
00579         {
00580             edge interface
00581             (
00582                 min(ownRegion, neiRegion),
00583                 max(ownRegion, neiRegion)
00584             );
00585 
00586             EdgeMap<label>::iterator iter = interfaceSizes.find(interface);
00587 
00588             if (iter != interfaceSizes.end())
00589             {
00590                 iter()++;
00591             }
00592             else
00593             {
00594                 interfaceSizes.insert(interface, 1);
00595             }
00596         }
00597     }
00598 
00599     // Boundary faces
00600     // ~~~~~~~~~~~~~~
00601 
00602     // Neighbour cellRegion.
00603     labelList coupledRegion(mesh.nFaces()-mesh.nInternalFaces());
00604 
00605     forAll(coupledRegion, i)
00606     {
00607         label cellI = mesh.faceOwner()[i+mesh.nInternalFaces()];
00608         coupledRegion[i] = cellRegion[cellI];
00609     }
00610     syncTools::swapBoundaryFaceList(mesh, coupledRegion, false);
00611 
00612     forAll(coupledRegion, i)
00613     {
00614         label faceI = i+mesh.nInternalFaces();
00615         label ownRegion = cellRegion[mesh.faceOwner()[faceI]];
00616         label neiRegion = coupledRegion[i];
00617 
00618         if (ownRegion != neiRegion)
00619         {
00620             edge interface
00621             (
00622                 min(ownRegion, neiRegion),
00623                 max(ownRegion, neiRegion)
00624             );
00625 
00626             EdgeMap<label>::iterator iter = interfaceSizes.find(interface);
00627 
00628             if (iter != interfaceSizes.end())
00629             {
00630                 iter()++;
00631             }
00632             else
00633             {
00634                 interfaceSizes.insert(interface, 1);
00635             }
00636         }
00637     }
00638 
00639 
00640     if (sumParallel && Pstream::parRun())
00641     {
00642         if (Pstream::master())
00643         {
00644             // Receive and add to my sizes
00645             for
00646             (
00647                 int slave=Pstream::firstSlave();
00648                 slave<=Pstream::lastSlave();
00649                 slave++
00650             )
00651             {
00652                 IPstream fromSlave(Pstream::blocking, slave);
00653 
00654                 EdgeMap<label> slaveSizes(fromSlave);
00655 
00656                 forAllConstIter(EdgeMap<label>, slaveSizes, slaveIter)
00657                 {
00658                     EdgeMap<label>::iterator masterIter =
00659                         interfaceSizes.find(slaveIter.key());
00660 
00661                     if (masterIter != interfaceSizes.end())
00662                     {
00663                         masterIter() += slaveIter();
00664                     }
00665                     else
00666                     {
00667                         interfaceSizes.insert(slaveIter.key(), slaveIter());
00668                     }
00669                 }
00670             }
00671 
00672             // Distribute
00673             for
00674             (
00675                 int slave=Pstream::firstSlave();
00676                 slave<=Pstream::lastSlave();
00677                 slave++
00678             )
00679             {
00680                 // Receive the edges using shared points from the slave.
00681                 OPstream toSlave(Pstream::blocking, slave);
00682                 toSlave << interfaceSizes;
00683             }
00684         }
00685         else
00686         {
00687             // Send to master
00688             {
00689                 OPstream toMaster(Pstream::blocking, Pstream::masterNo());
00690                 toMaster << interfaceSizes;
00691             }
00692             // Receive from master
00693             {
00694                 IPstream fromMaster(Pstream::blocking, Pstream::masterNo());
00695                 fromMaster >> interfaceSizes;
00696             }
00697         }
00698     }
00699 
00700     // Make sure all processors have interfaces in same order
00701     interfaces = interfaceSizes.toc();
00702     if (sumParallel)
00703     {
00704         Pstream::scatter(interfaces);
00705     }
00706 }
00707 
00708 
00709 // Create mesh for region.
00710 autoPtr<mapPolyMesh> createRegionMesh
00711 (
00712     const labelList& cellRegion,
00713     const EdgeMap<label>& interfaceToPatch,
00714     const fvMesh& mesh,
00715     const label regionI,
00716     const word& regionName,
00717     autoPtr<fvMesh>& newMesh
00718 )
00719 {
00720     // Create dummy system/fv*
00721     {
00722         IOobject io
00723         (
00724             "fvSchemes",
00725             mesh.time().system(),
00726             regionName,
00727             mesh,
00728             IOobject::NO_READ,
00729             IOobject::NO_WRITE,
00730             false
00731         );
00732 
00733         Info<< "Testing:" << io.objectPath() << endl;
00734 
00735         if (!io.headerOk())
00736         // if (!exists(io.objectPath()))
00737         {
00738             Info<< "Writing dummy " << regionName/io.name() << endl;
00739             dictionary dummyDict;
00740             dictionary divDict;
00741             dummyDict.add("divSchemes", divDict);
00742             dictionary gradDict;
00743             dummyDict.add("gradSchemes", gradDict);
00744             dictionary laplDict;
00745             dummyDict.add("laplacianSchemes", laplDict);
00746 
00747             IOdictionary(io, dummyDict).regIOobject::write();
00748         }
00749     }
00750     {
00751         IOobject io
00752         (
00753             "fvSolution",
00754             mesh.time().system(),
00755             regionName,
00756             mesh,
00757             IOobject::NO_READ,
00758             IOobject::NO_WRITE,
00759             false
00760         );
00761 
00762         if (!io.headerOk())
00763         //if (!exists(io.objectPath()))
00764         {
00765             Info<< "Writing dummy " << regionName/io.name() << endl;
00766             dictionary dummyDict;
00767             IOdictionary(io, dummyDict).regIOobject::write();
00768         }
00769     }
00770 
00771 
00772     // Neighbour cellRegion.
00773     labelList coupledRegion(mesh.nFaces()-mesh.nInternalFaces());
00774 
00775     forAll(coupledRegion, i)
00776     {
00777         label cellI = mesh.faceOwner()[i+mesh.nInternalFaces()];
00778         coupledRegion[i] = cellRegion[cellI];
00779     }
00780     syncTools::swapBoundaryFaceList(mesh, coupledRegion, false);
00781 
00782 
00783     // Topology change container. Start off from existing mesh.
00784     polyTopoChange meshMod(mesh);
00785 
00786     // Cell remover engine
00787     removeCells cellRemover(mesh);
00788 
00789     // Select all but region cells
00790     labelList cellsToRemove(getNonRegionCells(cellRegion, regionI));
00791 
00792     // Find out which faces will get exposed. Note that this
00793     // gets faces in mesh face order. So both regions will get same
00794     // face in same order (important!)
00795     labelList exposedFaces = cellRemover.getExposedFaces(cellsToRemove);
00796 
00797     labelList exposedPatchIDs(exposedFaces.size());
00798     forAll(exposedFaces, i)
00799     {
00800         label faceI = exposedFaces[i];
00801 
00802         label ownRegion = cellRegion[mesh.faceOwner()[faceI]];
00803         label neiRegion = -1;
00804 
00805         if (mesh.isInternalFace(faceI))
00806         {
00807             neiRegion = cellRegion[mesh.faceNeighbour()[faceI]];
00808         }
00809         else
00810         {
00811             neiRegion = coupledRegion[faceI-mesh.nInternalFaces()];
00812         }
00813 
00814         label otherRegion = -1;
00815 
00816         if (ownRegion == regionI && neiRegion != regionI)
00817         {
00818             otherRegion = neiRegion;
00819         }
00820         else if (ownRegion != regionI && neiRegion == regionI)
00821         {
00822             otherRegion = ownRegion;
00823         }
00824         else
00825         {
00826             FatalErrorIn("createRegionMesh(..)")
00827                 << "Exposed face:" << faceI
00828                 << " fc:" << mesh.faceCentres()[faceI]
00829                 << " has owner region " << ownRegion
00830                 << " and neighbour region " << neiRegion
00831                 << " when handling region:" << regionI
00832                 << exit(FatalError);
00833         }
00834 
00835         if (otherRegion != -1)
00836         {
00837             edge interface(regionI, otherRegion);
00838 
00839             // Find the patch.
00840             if (regionI < otherRegion)
00841             {
00842                 exposedPatchIDs[i] = interfaceToPatch[interface];
00843             }
00844             else
00845             {
00846                 exposedPatchIDs[i] = interfaceToPatch[interface]+1;
00847             }
00848         }
00849     }
00850 
00851     // Remove faces
00852     cellRemover.setRefinement
00853     (
00854         cellsToRemove,
00855         exposedFaces,
00856         exposedPatchIDs,
00857         meshMod
00858     );
00859 
00860     autoPtr<mapPolyMesh> map = meshMod.makeMesh
00861     (
00862         newMesh,
00863         IOobject
00864         (
00865             regionName,
00866             mesh.time().timeName(),
00867             mesh.time(),
00868             IOobject::NO_READ,
00869             IOobject::AUTO_WRITE
00870         ),
00871         mesh
00872     );
00873 
00874     return map;
00875 }
00876 
00877 
00878 void createAndWriteRegion
00879 (
00880     const fvMesh& mesh,
00881     const labelList& cellRegion,
00882     const wordList& regionNames,
00883     const EdgeMap<label>& interfaceToPatch,
00884     const label regionI,
00885     const word& newMeshInstance
00886 )
00887 {
00888     Info<< "Creating mesh for region " << regionI
00889         << ' ' << regionNames[regionI] << endl;
00890 
00891     autoPtr<fvMesh> newMesh;
00892     autoPtr<mapPolyMesh> map = createRegionMesh
00893     (
00894         cellRegion,
00895         interfaceToPatch,
00896         mesh,
00897         regionI,
00898         regionNames[regionI],
00899         newMesh
00900     );
00901 
00902 
00903     // Make map of all added patches
00904     labelHashSet addedPatches(2*interfaceToPatch.size());
00905     forAllConstIter(EdgeMap<label>, interfaceToPatch, iter)
00906     {
00907         addedPatches.insert(iter());
00908         addedPatches.insert(iter()+1);
00909     }
00910 
00911     Info<< "Mapping fields" << endl;
00912 
00913     // Map existing fields
00914     newMesh().updateMesh(map());
00915 
00916     // Add subsetted fields
00917     subsetVolFields<volScalarField>
00918     (
00919         mesh,
00920         newMesh(),
00921         map().cellMap(),
00922         map().faceMap(),
00923         addedPatches
00924     );
00925     subsetVolFields<volVectorField>
00926     (
00927         mesh,
00928         newMesh(),
00929         map().cellMap(),
00930         map().faceMap(),
00931         addedPatches
00932     );
00933     subsetVolFields<volSphericalTensorField>
00934     (
00935         mesh,
00936         newMesh(),
00937         map().cellMap(),
00938         map().faceMap(),
00939         addedPatches
00940     );
00941     subsetVolFields<volSymmTensorField>
00942     (
00943         mesh,
00944         newMesh(),
00945         map().cellMap(),
00946         map().faceMap(),
00947         addedPatches
00948     );
00949     subsetVolFields<volTensorField>
00950     (
00951         mesh,
00952         newMesh(),
00953         map().cellMap(),
00954         map().faceMap(),
00955         addedPatches
00956     );
00957 
00958     subsetSurfaceFields<surfaceScalarField>
00959     (
00960         mesh,
00961         newMesh(),
00962         map().faceMap(),
00963         addedPatches
00964     );
00965     subsetSurfaceFields<surfaceVectorField>
00966     (
00967         mesh,
00968         newMesh(),
00969         map().faceMap(),
00970         addedPatches
00971     );
00972     subsetSurfaceFields<surfaceSphericalTensorField>
00973     (
00974         mesh,
00975         newMesh(),
00976         map().faceMap(),
00977         addedPatches
00978     );
00979     subsetSurfaceFields<surfaceSymmTensorField>
00980     (
00981         mesh,
00982         newMesh(),
00983         map().faceMap(),
00984         addedPatches
00985     );
00986     subsetSurfaceFields<surfaceTensorField>
00987     (
00988         mesh,
00989         newMesh(),
00990         map().faceMap(),
00991         addedPatches
00992     );
00993 
00994 
00995     const polyBoundaryMesh& newPatches = newMesh().boundaryMesh();
00996     newPatches.checkParallelSync(true);
00997 
00998     // Delete empty patches
00999     // ~~~~~~~~~~~~~~~~~~~~
01000 
01001     // Create reordering list to move patches-to-be-deleted to end
01002     labelList oldToNew(newPatches.size(), -1);
01003     label newI = 0;
01004 
01005     Info<< "Deleting empty patches" << endl;
01006 
01007     // Assumes all non-proc boundaries are on all processors!
01008     forAll(newPatches, patchI)
01009     {
01010         const polyPatch& pp = newPatches[patchI];
01011 
01012         if (!isA<processorPolyPatch>(pp))
01013         {
01014             if (returnReduce(pp.size(), sumOp<label>()) > 0)
01015             {
01016                 oldToNew[patchI] = newI++;
01017             }
01018         }
01019     }
01020 
01021     // Same for processor patches (but need no reduction)
01022     forAll(newPatches, patchI)
01023     {
01024         const polyPatch& pp = newPatches[patchI];
01025 
01026         if (isA<processorPolyPatch>(pp) && pp.size())
01027         {
01028             oldToNew[patchI] = newI++;
01029         }
01030     }
01031 
01032     const label nNewPatches = newI;
01033 
01034     // Move all deleteable patches to the end
01035     forAll(oldToNew, patchI)
01036     {
01037         if (oldToNew[patchI] == -1)
01038         {
01039             oldToNew[patchI] = newI++;
01040         }
01041     }
01042 
01043     reorderPatches(newMesh(), oldToNew, nNewPatches);
01044 
01045 
01046     Info<< "Writing new mesh" << endl;
01047 
01048     newMesh().setInstance(newMeshInstance);
01049     newMesh().write();
01050 
01051     // Write addressing files like decomposePar
01052     Info<< "Writing addressing to base mesh" << endl;
01053 
01054     labelIOList pointProcAddressing
01055     (
01056         IOobject
01057         (
01058             "pointRegionAddressing",
01059             newMesh().facesInstance(),
01060             newMesh().meshSubDir,
01061             newMesh(),
01062             IOobject::NO_READ,
01063             IOobject::NO_WRITE,
01064             false
01065         ),
01066         map().pointMap()
01067     );
01068     Info<< "Writing map " << pointProcAddressing.name()
01069         << " from region" << regionI
01070         << " points back to base mesh." << endl;
01071     pointProcAddressing.write();
01072 
01073     labelIOList faceProcAddressing
01074     (
01075         IOobject
01076         (
01077             "faceRegionAddressing",
01078             newMesh().facesInstance(),
01079             newMesh().meshSubDir,
01080             newMesh(),
01081             IOobject::NO_READ,
01082             IOobject::NO_WRITE,
01083             false
01084         ),
01085         newMesh().nFaces()
01086     );
01087     forAll(faceProcAddressing, faceI)
01088     {
01089         // face + turning index. (see decomposePar)
01090         // Is the face pointing in the same direction?
01091         label oldFaceI = map().faceMap()[faceI];
01092 
01093         if
01094         (
01095             map().cellMap()[newMesh().faceOwner()[faceI]]
01096          == mesh.faceOwner()[oldFaceI]
01097         )
01098         {
01099             faceProcAddressing[faceI] = oldFaceI+1;
01100         }
01101         else
01102         {
01103             faceProcAddressing[faceI] = -(oldFaceI+1);
01104         }
01105     }
01106     Info<< "Writing map " << faceProcAddressing.name()
01107         << " from region" << regionI
01108         << " faces back to base mesh." << endl;
01109     faceProcAddressing.write();
01110 
01111     labelIOList cellProcAddressing
01112     (
01113         IOobject
01114         (
01115             "cellRegionAddressing",
01116             newMesh().facesInstance(),
01117             newMesh().meshSubDir,
01118             newMesh(),
01119             IOobject::NO_READ,
01120             IOobject::NO_WRITE,
01121             false
01122         ),
01123         map().cellMap()
01124     );
01125     Info<< "Writing map " <<cellProcAddressing.name()
01126         << " from region" << regionI
01127         << " cells back to base mesh." << endl;
01128     cellProcAddressing.write();
01129 }
01130 
01131 
01132 // Create for every region-region interface with non-zero size two patches.
01133 // First one is for minimumregion to maximumregion.
01134 // Note that patches get created in same order on all processors (if parallel)
01135 // since looping over synchronised 'interfaces'.
01136 EdgeMap<label> addRegionPatches
01137 (
01138     fvMesh& mesh,
01139     const labelList& cellRegion,
01140     const label nCellRegions,
01141     const edgeList& interfaces,
01142     const EdgeMap<label>& interfaceSizes,
01143     const wordList& regionNames
01144 )
01145 {
01146     // Check that all patches are present in same order.
01147     mesh.boundaryMesh().checkParallelSync(true);
01148 
01149     Info<< nl << "Adding patches" << nl << endl;
01150 
01151     EdgeMap<label> interfaceToPatch(nCellRegions);
01152 
01153     forAll(interfaces, interI)
01154     {
01155         const edge& e = interfaces[interI];
01156 
01157         if (interfaceSizes[e] > 0)
01158         {
01159             const word inter1 = regionNames[e[0]] + "_to_" + regionNames[e[1]];
01160             const word inter2 = regionNames[e[1]] + "_to_" + regionNames[e[0]];
01161 
01162             directMappedWallPolyPatch patch1
01163             (
01164                 inter1,
01165                 0,                  // overridden
01166                 0,                  // overridden
01167                 0,                  // overridden
01168                 regionNames[e[1]],  // sampleRegion
01169                 directMappedPatchBase::NEARESTPATCHFACE,
01170                 inter2,             // samplePatch
01171                 point::zero,        // offset
01172                 mesh.boundaryMesh()
01173             );
01174 
01175             label patchI = addPatch(mesh, patch1);
01176 
01177             directMappedWallPolyPatch patch2
01178             (
01179                 inter2,
01180                 0,
01181                 0,
01182                 0,
01183                 regionNames[e[0]],  // sampleRegion
01184                 directMappedPatchBase::NEARESTPATCHFACE,
01185                 inter1,
01186                 point::zero,        // offset
01187                 mesh.boundaryMesh()
01188             );
01189             addPatch(mesh, patch2);
01190 
01191             Info<< "For interface between region " << e[0]
01192                 << " and " << e[1] << " added patch " << patchI
01193                 << " " << mesh.boundaryMesh()[patchI].name()
01194                 << endl;
01195 
01196             interfaceToPatch.insert(e, patchI);
01197         }
01198     }
01199     return interfaceToPatch;
01200 }
01201 
01202 
01203 // Find region that covers most of cell zone
01204 label findCorrespondingRegion
01205 (
01206     const labelList& existingZoneID,    // per cell the (unique) zoneID
01207     const labelList& cellRegion,
01208     const label nCellRegions,
01209     const label zoneI,
01210     const label minOverlapSize
01211 )
01212 {
01213     // Per region the number of cells in zoneI
01214     labelList cellsInZone(nCellRegions, 0);
01215 
01216     forAll(cellRegion, cellI)
01217     {
01218         if (existingZoneID[cellI] == zoneI)
01219         {
01220             cellsInZone[cellRegion[cellI]]++;
01221         }
01222     }
01223 
01224     Pstream::listCombineGather(cellsInZone, plusEqOp<label>());
01225     Pstream::listCombineScatter(cellsInZone);
01226 
01227     // Pick region with largest overlap of zoneI
01228     label regionI = findMax(cellsInZone);
01229 
01230 
01231     if (cellsInZone[regionI] < minOverlapSize)
01232     {
01233         // Region covers too little of zone. Not good enough.
01234         regionI = -1;
01235     }
01236     else
01237     {
01238         // Check that region contains no cells that aren't in cellZone.
01239         forAll(cellRegion, cellI)
01240         {
01241             if (cellRegion[cellI] == regionI && existingZoneID[cellI] != zoneI)
01242             {
01243                 // cellI in regionI but not in zoneI
01244                 regionI = -1;
01245                 break;
01246             }
01247         }
01248         // If one in error, all should be in error. Note that branch gets taken
01249         // on all procs.
01250         reduce(regionI, minOp<label>());
01251     }
01252 
01253     return regionI;
01254 }
01255 
01256 
01258 //label findCorrespondingRegion
01259 //(
01260 //    const cellZoneMesh& cellZones,
01261 //    const labelList& existingZoneID,    // per cell the (unique) zoneID
01262 //    const labelList& cellRegion,
01263 //    const label nCellRegions,
01264 //    const label zoneI
01265 //)
01266 //{
01267 //    // Region corresponding to zone. Start off with special value: no
01268 //    // corresponding region.
01269 //    label regionI = labelMax;
01270 //
01271 //    const cellZone& cz = cellZones[zoneI];
01272 //
01273 //    if (cz.empty())
01274 //    {
01275 //        // My local portion is empty. Maps to any empty cellZone. Mark with
01276 //        // special value which can get overwritten by other processors.
01277 //        regionI = -1;
01278 //    }
01279 //    else
01280 //    {
01281 //        regionI = cellRegion[cz[0]];
01282 //
01283 //        forAll(cz, i)
01284 //        {
01285 //            if (cellRegion[cz[i]] != regionI)
01286 //            {
01287 //                regionI = labelMax;
01288 //                break;
01289 //            }
01290 //        }
01291 //    }
01292 //
01293 //    // Determine same zone over all processors.
01294 //    reduce(regionI, maxOp<label>());
01295 //
01296 //
01297 //    // 2. All of region present?
01298 //
01299 //    if (regionI == labelMax)
01300 //    {
01301 //        regionI = -1;
01302 //    }
01303 //    else if (regionI != -1)
01304 //    {
01305 //        forAll(cellRegion, cellI)
01306 //        {
01307 //            if (cellRegion[cellI] == regionI && existingZoneID[cellI] != zoneI)
01308 //            {
01309 //                // cellI in regionI but not in zoneI
01310 //                regionI = -1;
01311 //                break;
01312 //            }
01313 //        }
01314 //        // If one in error, all should be in error. Note that branch gets taken
01315 //        // on all procs.
01316 //        reduce(regionI, minOp<label>());
01317 //    }
01318 //
01319 //    return regionI;
01320 //}
01321 
01322 
01323 // Get zone per cell
01324 // - non-unique zoning
01325 // - coupled zones
01326 void getZoneID
01327 (
01328     const polyMesh& mesh,
01329     const cellZoneMesh& cellZones,
01330     labelList& zoneID,
01331     labelList& neiZoneID
01332 )
01333 {
01334     // Existing zoneID
01335     zoneID.setSize(mesh.nCells());
01336     zoneID = -1;
01337 
01338     forAll(cellZones, zoneI)
01339     {
01340         const cellZone& cz = cellZones[zoneI];
01341 
01342         forAll(cz, i)
01343         {
01344             label cellI = cz[i];
01345             if (zoneID[cellI] == -1)
01346             {
01347                 zoneID[cellI] = zoneI;
01348             }
01349             else
01350             {
01351                 FatalErrorIn("getZoneID(..)")
01352                     << "Cell " << cellI << " with cell centre "
01353                     << mesh.cellCentres()[cellI]
01354                     << " is multiple zones. This is not allowed." << endl
01355                     << "It is in zone " << cellZones[zoneID[cellI]].name()
01356                     << " and in zone " << cellZones[zoneI].name()
01357                     << exit(FatalError);
01358             }
01359         }
01360     }
01361 
01362     // Neighbour zoneID.
01363     neiZoneID.setSize(mesh.nFaces()-mesh.nInternalFaces());
01364 
01365     forAll(neiZoneID, i)
01366     {
01367         neiZoneID[i] = zoneID[mesh.faceOwner()[i+mesh.nInternalFaces()]];
01368     }
01369     syncTools::swapBoundaryFaceList(mesh, neiZoneID, false);
01370 }
01371 
01372 
01373 void matchRegions
01374 (
01375     const bool sloppyCellZones,
01376     const polyMesh& mesh,
01377 
01378     const label nCellRegions,
01379     const labelList& cellRegion,
01380 
01381     labelList& regionToZone,
01382     wordList& regionNames,
01383     labelList& zoneToRegion
01384 )
01385 {
01386     const cellZoneMesh& cellZones = mesh.cellZones();
01387 
01388     regionToZone.setSize(nCellRegions, -1);
01389     regionNames.setSize(nCellRegions);
01390     zoneToRegion.setSize(cellZones.size(), -1);
01391 
01392     // Get current per cell zoneID
01393     labelList zoneID(mesh.nCells(), -1);
01394     labelList neiZoneID(mesh.nFaces()-mesh.nInternalFaces());
01395     getZoneID(mesh, cellZones, zoneID, neiZoneID);
01396 
01397     // Sizes per cellzone
01398     labelList zoneSizes(cellZones.size(), 0);
01399     {
01400         List<wordList> zoneNames(Pstream::nProcs());
01401         zoneNames[Pstream::myProcNo()] = cellZones.names();
01402         Pstream::gatherList(zoneNames);
01403         Pstream::scatterList(zoneNames);
01404 
01405         forAll(zoneNames, procI)
01406         {
01407             if (zoneNames[procI] != zoneNames[0])
01408             {
01409                 FatalErrorIn("matchRegions(..)")
01410                     << "cellZones not synchronised across processors." << endl
01411                     << "Master has cellZones " << zoneNames[0] << endl
01412                     << "Processor " << procI
01413                     << " has cellZones " << zoneNames[procI]
01414                     << exit(FatalError);
01415             }
01416         }
01417 
01418         forAll(cellZones, zoneI)
01419         {
01420             zoneSizes[zoneI] = returnReduce
01421             (
01422                 cellZones[zoneI].size(),
01423                 sumOp<label>()
01424             );
01425         }
01426     }
01427 
01428 
01429     if (sloppyCellZones)
01430     {
01431         Info<< "Trying to match regions to existing cell zones;"
01432             << " region can be subset of cell zone." << nl << endl;
01433 
01434         forAll(cellZones, zoneI)
01435         {
01436             label regionI = findCorrespondingRegion
01437             (
01438                 zoneID,
01439                 cellRegion,
01440                 nCellRegions,
01441                 zoneI,
01442                 label(0.5*zoneSizes[zoneI]) // minimum overlap
01443             );
01444 
01445             if (regionI != -1)
01446             {
01447                 Info<< "Sloppily matched region " << regionI
01448                     //<< " size " << regionSizes[regionI]
01449                     << " to zone " << zoneI << " size " << zoneSizes[zoneI]
01450                     << endl;
01451                 zoneToRegion[zoneI] = regionI;
01452                 regionToZone[regionI] = zoneI;
01453                 regionNames[regionI] = cellZones[zoneI].name();
01454             }
01455         }
01456     }
01457     else
01458     {
01459         Info<< "Trying to match regions to existing cell zones." << nl << endl;
01460 
01461         forAll(cellZones, zoneI)
01462         {
01463             label regionI = findCorrespondingRegion
01464             (
01465                 zoneID,
01466                 cellRegion,
01467                 nCellRegions,
01468                 zoneI,
01469                 1               // minimum overlap
01470             );
01471 
01472             if (regionI != -1)
01473             {
01474                 zoneToRegion[zoneI] = regionI;
01475                 regionToZone[regionI] = zoneI;
01476                 regionNames[regionI] = cellZones[zoneI].name();
01477             }
01478         }
01479     }
01480     // Allocate region names for unmatched regions.
01481     forAll(regionToZone, regionI)
01482     {
01483         if (regionToZone[regionI] == -1)
01484         {
01485             regionNames[regionI] = "domain" + Foam::name(regionI);
01486         }
01487     }
01488 }
01489 
01490 
01491 void writeCellToRegion(const fvMesh& mesh, const labelList& cellRegion)
01492 {
01493     // Write to manual decomposition option
01494     {
01495         labelIOList cellToRegion
01496         (
01497             IOobject
01498             (
01499                 "cellToRegion",
01500                 mesh.facesInstance(),
01501                 mesh,
01502                 IOobject::NO_READ,
01503                 IOobject::NO_WRITE,
01504                 false
01505             ),
01506             cellRegion
01507         );
01508         cellToRegion.write();
01509 
01510         Info<< "Writing region per cell file (for manual decomposition) to "
01511             << cellToRegion.objectPath() << nl << endl;
01512     }
01513     // Write for postprocessing
01514     {
01515         volScalarField cellToRegion
01516         (
01517             IOobject
01518             (
01519                 "cellToRegion",
01520                 mesh.time().timeName(),
01521                 mesh,
01522                 IOobject::NO_READ,
01523                 IOobject::NO_WRITE,
01524                 false
01525             ),
01526             mesh,
01527             dimensionedScalar("zero", dimless, 0),
01528             zeroGradientFvPatchScalarField::typeName
01529         );
01530         forAll(cellRegion, cellI)
01531         {
01532             cellToRegion[cellI] = cellRegion[cellI];
01533         }
01534         cellToRegion.write();
01535 
01536         Info<< "Writing region per cell as volScalarField to "
01537             << cellToRegion.objectPath() << nl << endl;
01538     }
01539 }
01540 
01541 
01542 
01543 // Main program:
01544 
01545 int main(int argc, char *argv[])
01546 {
01547     argList::validOptions.insert("cellZones", "");
01548     argList::validOptions.insert("cellZonesOnly", "");
01549     argList::validOptions.insert("cellZonesFileOnly", "cellZonesName");
01550     argList::validOptions.insert("blockedFaces", "faceSet");
01551     argList::validOptions.insert("makeCellZones", "");
01552     argList::validOptions.insert("largestOnly", "");
01553     argList::validOptions.insert("insidePoint", "point");
01554     argList::validOptions.insert("overwrite", "");
01555     argList::validOptions.insert("detectOnly", "");
01556     argList::validOptions.insert("sloppyCellZones", "");
01557 
01558 #   include <OpenFOAM/setRootCase.H>
01559 #   include <OpenFOAM/createTime.H>
01560     runTime.functionObjects().off();
01561 #   include <OpenFOAM/createMesh.H>
01562     const word oldInstance = mesh.pointsInstance();
01563 
01564     word blockedFacesName;
01565     if (args.optionFound("blockedFaces"))
01566     {
01567         blockedFacesName = args.option("blockedFaces");
01568         Info<< "Reading blocked internal faces from faceSet "
01569             << blockedFacesName << nl << endl;
01570     }
01571 
01572     bool makeCellZones    = args.optionFound("makeCellZones");
01573     bool largestOnly      = args.optionFound("largestOnly");
01574     bool insidePoint      = args.optionFound("insidePoint");
01575     bool useCellZones     = args.optionFound("cellZones");
01576     bool useCellZonesOnly = args.optionFound("cellZonesOnly");
01577     bool useCellZonesFile = args.optionFound("cellZonesFileOnly");
01578     bool overwrite        = args.optionFound("overwrite");
01579     bool detectOnly       = args.optionFound("detectOnly");
01580     bool sloppyCellZones  = args.optionFound("sloppyCellZones");
01581 
01582     if
01583     (
01584         (useCellZonesOnly || useCellZonesFile)
01585      && (
01586             blockedFacesName != word::null
01587          || useCellZones
01588         )
01589     )
01590     {
01591         FatalErrorIn(args.executable())
01592             << "You cannot specify both -cellZonesOnly or -cellZonesFileOnly"
01593             << " (which specify complete split)"
01594             << " in combination with -blockedFaces or -cellZones"
01595             << " (which imply a split based on topology)"
01596             << exit(FatalError);
01597     }
01598 
01599 
01600 
01601     if (insidePoint && largestOnly)
01602     {
01603         FatalErrorIn(args.executable())
01604             << "You cannot specify both -largestOnly"
01605             << " (keep region with most cells)"
01606             << " and -insidePoint (keep region containing point)"
01607             << exit(FatalError);
01608     }
01609 
01610 
01611     const cellZoneMesh& cellZones = mesh.cellZones();
01612 
01613     // Existing zoneID
01614     labelList zoneID(mesh.nCells(), -1);
01615     // Neighbour zoneID.
01616     labelList neiZoneID(mesh.nFaces()-mesh.nInternalFaces());
01617     getZoneID(mesh, cellZones, zoneID, neiZoneID);
01618 
01619 
01620 
01621     // Determine per cell the region it belongs to
01622     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
01623 
01624     // cellRegion is the labelList with the region per cell.
01625     labelList cellRegion;
01626     // Region per zone
01627     labelList regionToZone;
01628     // Name of region
01629     wordList regionNames;
01630     // Zone to region
01631     labelList zoneToRegion;
01632 
01633     label nCellRegions = 0;
01634     if (useCellZonesOnly)
01635     {
01636         Info<< "Using current cellZones to split mesh into regions."
01637             << " This requires all"
01638             << " cells to be in one and only one cellZone." << nl << endl;
01639 
01640         label unzonedCellI = findIndex(zoneID, -1);
01641         if (unzonedCellI != -1)
01642         {
01643             FatalErrorIn(args.executable())
01644                 << "For the cellZonesOnly option all cells "
01645                 << "have to be in a cellZone." << endl
01646                 << "Cell " << unzonedCellI
01647                 << " at" << mesh.cellCentres()[unzonedCellI]
01648                 << " is not in a cellZone. There might be more unzoned cells."
01649                 << exit(FatalError);
01650         }
01651         cellRegion = zoneID;
01652         nCellRegions = gMax(cellRegion)+1;
01653         regionToZone.setSize(nCellRegions);
01654         regionNames.setSize(nCellRegions);
01655         zoneToRegion.setSize(cellZones.size(), -1);
01656         for (label regionI = 0; regionI < nCellRegions; regionI++)
01657         {
01658             regionToZone[regionI] = regionI;
01659             zoneToRegion[regionI] = regionI;
01660             regionNames[regionI] = cellZones[regionI].name();
01661         }
01662     }
01663     else if (useCellZonesFile)
01664     {
01665         const word zoneFile = args.option("cellZonesFileOnly");
01666         Info<< "Reading split from cellZones file " << zoneFile << endl
01667             << "This requires all"
01668             << " cells to be in one and only one cellZone." << nl << endl;
01669 
01670         cellZoneMesh newCellZones
01671         (
01672             IOobject
01673             (
01674                 zoneFile,
01675                 mesh.facesInstance(),
01676                 polyMesh::meshSubDir,
01677                 mesh,
01678                 IOobject::MUST_READ,
01679                 IOobject::NO_WRITE,
01680                 false
01681             ),
01682             mesh
01683         );
01684 
01685         labelList newZoneID(mesh.nCells(), -1);
01686         labelList newNeiZoneID(mesh.nFaces()-mesh.nInternalFaces());
01687         getZoneID(mesh, newCellZones, newZoneID, newNeiZoneID);
01688 
01689         label unzonedCellI = findIndex(newZoneID, -1);
01690         if (unzonedCellI != -1)
01691         {
01692             FatalErrorIn(args.executable())
01693                 << "For the cellZonesFileOnly option all cells "
01694                 << "have to be in a cellZone." << endl
01695                 << "Cell " << unzonedCellI
01696                 << " at" << mesh.cellCentres()[unzonedCellI]
01697                 << " is not in a cellZone. There might be more unzoned cells."
01698                 << exit(FatalError);
01699         }
01700         cellRegion = newZoneID;
01701         nCellRegions = gMax(cellRegion)+1;
01702         zoneToRegion.setSize(newCellZones.size(), -1);
01703         regionToZone.setSize(nCellRegions);
01704         regionNames.setSize(nCellRegions);
01705         for (label regionI = 0; regionI < nCellRegions; regionI++)
01706         {
01707             regionToZone[regionI] = regionI;
01708             zoneToRegion[regionI] = regionI;
01709             regionNames[regionI] = newCellZones[regionI].name();
01710         }
01711     }
01712     else
01713     {
01714         // Determine connected regions
01715         // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
01716 
01717         // Mark additional faces that are blocked
01718         boolList blockedFace;
01719 
01720         // Read from faceSet
01721         if (blockedFacesName.size())
01722         {
01723             faceSet blockedFaceSet(mesh, blockedFacesName);
01724             Info<< "Read "
01725                 << returnReduce(blockedFaceSet.size(), sumOp<label>())
01726                 << " blocked faces from set " << blockedFacesName << nl << endl;
01727 
01728             blockedFace.setSize(mesh.nFaces(), false);
01729 
01730             forAllConstIter(faceSet, blockedFaceSet, iter)
01731             {
01732                 blockedFace[iter.key()] = true;
01733             }
01734         }
01735 
01736         // Imply from differing cellZones
01737         if (useCellZones)
01738         {
01739             blockedFace.setSize(mesh.nFaces(), false);
01740 
01741             for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
01742             {
01743                 label own = mesh.faceOwner()[faceI];
01744                 label nei = mesh.faceNeighbour()[faceI];
01745 
01746                 if (zoneID[own] != zoneID[nei])
01747                 {
01748                     blockedFace[faceI] = true;
01749                 }
01750             }
01751 
01752             // Different cellZones on either side of processor patch.
01753             forAll(neiZoneID, i)
01754             {
01755                 label faceI = i+mesh.nInternalFaces();
01756 
01757                 if (zoneID[mesh.faceOwner()[faceI]] != neiZoneID[i])
01758                 {
01759                     blockedFace[faceI] = true;
01760                 }
01761             }
01762         }
01763 
01764         // Do a topological walk to determine regions
01765         regionSplit regions(mesh, blockedFace);
01766         nCellRegions = regions.nRegions();
01767         cellRegion.transfer(regions);
01768 
01769         // Make up region names. If possible match them to existing zones.
01770         matchRegions
01771         (
01772             sloppyCellZones,
01773             mesh,
01774             nCellRegions,
01775             cellRegion,
01776 
01777             regionToZone,
01778             regionNames,
01779             zoneToRegion
01780         );
01781     }
01782 
01783     Info<< endl << "Number of regions:" << nCellRegions << nl << endl;
01784 
01785 
01786     // Write decomposition to file
01787     writeCellToRegion(mesh, cellRegion);
01788 
01789 
01790 
01791     // Sizes per region
01792     // ~~~~~~~~~~~~~~~~
01793 
01794     labelList regionSizes(nCellRegions, 0);
01795 
01796     forAll(cellRegion, cellI)
01797     {
01798         regionSizes[cellRegion[cellI]]++;
01799     }
01800     forAll(regionSizes, regionI)
01801     {
01802         reduce(regionSizes[regionI], sumOp<label>());
01803     }
01804 
01805     Info<< "Region\tCells" << nl
01806         << "------\t-----" << endl;
01807 
01808     forAll(regionSizes, regionI)
01809     {
01810         Info<< regionI << '\t' << regionSizes[regionI] << nl;
01811     }
01812     Info<< endl;
01813 
01814 
01815 
01816     // Print region to zone
01817     Info<< "Region\tZone\tName" << nl
01818         << "------\t----\t----" << endl;
01819     forAll(regionToZone, regionI)
01820     {
01821         Info<< regionI << '\t' << regionToZone[regionI] << '\t'
01822             << regionNames[regionI] << nl;
01823     }
01824     Info<< endl;
01825 
01826 
01827 
01828     // Since we're going to mess with patches make sure all non-processor ones
01829     // are on all processors.
01830     mesh.boundaryMesh().checkParallelSync(true);
01831 
01832 
01833     // Sizes of interface between regions. From pair of regions to number of
01834     // faces.
01835     edgeList interfaces;
01836     EdgeMap<label> interfaceSizes;
01837     getInterfaceSizes
01838     (
01839         mesh,
01840         cellRegion,
01841         true,      // sum in parallel?
01842 
01843         interfaces,
01844         interfaceSizes
01845     );
01846 
01847     Info<< "Sizes inbetween regions:" << nl << nl
01848         << "Region\tRegion\tFaces" << nl
01849         << "------\t------\t-----" << endl;
01850 
01851     forAll(interfaces, interI)
01852     {
01853         const edge& e = interfaces[interI];
01854 
01855         Info<< e[0] << '\t' << e[1] << '\t' << interfaceSizes[e] << nl;
01856     }
01857     Info<< endl;
01858 
01859 
01860     if (detectOnly)
01861     {
01862         return 0;
01863     }
01864 
01865 
01866     // Read objects in time directory
01867     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
01868 
01869     IOobjectList objects(mesh, runTime.timeName());
01870 
01871     // Read vol fields.
01872 
01873     PtrList<volScalarField> vsFlds;
01874     ReadFields(mesh, objects, vsFlds);
01875 
01876     PtrList<volVectorField> vvFlds;
01877     ReadFields(mesh, objects, vvFlds);
01878 
01879     PtrList<volSphericalTensorField> vstFlds;
01880     ReadFields(mesh, objects, vstFlds);
01881 
01882     PtrList<volSymmTensorField> vsymtFlds;
01883     ReadFields(mesh, objects, vsymtFlds);
01884 
01885     PtrList<volTensorField> vtFlds;
01886     ReadFields(mesh, objects, vtFlds);
01887 
01888     // Read surface fields.
01889 
01890     PtrList<surfaceScalarField> ssFlds;
01891     ReadFields(mesh, objects, ssFlds);
01892 
01893     PtrList<surfaceVectorField> svFlds;
01894     ReadFields(mesh, objects, svFlds);
01895 
01896     PtrList<surfaceSphericalTensorField> sstFlds;
01897     ReadFields(mesh, objects, sstFlds);
01898 
01899     PtrList<surfaceSymmTensorField> ssymtFlds;
01900     ReadFields(mesh, objects, ssymtFlds);
01901 
01902     PtrList<surfaceTensorField> stFlds;
01903     ReadFields(mesh, objects, stFlds);
01904 
01905     Info<< endl;
01906 
01907 
01908     // Remove any demand-driven fields ('S', 'V' etc)
01909     mesh.clearOut();
01910 
01911 
01912     if (nCellRegions == 1)
01913     {
01914         Info<< "Only one region. Doing nothing." << endl;
01915     }
01916     else if (makeCellZones)
01917     {
01918         Info<< "Putting cells into cellZones instead of splitting mesh."
01919             << endl;
01920 
01921         // Check if region overlaps with existing zone. If so keep.
01922 
01923         for (label regionI = 0; regionI < nCellRegions; regionI++)
01924         {
01925             label zoneI = regionToZone[regionI];
01926 
01927             if (zoneI != -1)
01928             {
01929                 Info<< "    Region " << regionI << " : corresponds to existing"
01930                     << " cellZone "
01931                     << zoneI << ' ' << cellZones[zoneI].name() << endl;
01932             }
01933             else
01934             {
01935                 // Create new cellZone.
01936                 labelList regionCells = findIndices(cellRegion, regionI);
01937 
01938                 word zoneName = "region" + Foam::name(regionI);
01939 
01940                 zoneI = cellZones.findZoneID(zoneName);
01941 
01942                 if (zoneI == -1)
01943                 {
01944                     zoneI = cellZones.size();
01945                     mesh.cellZones().setSize(zoneI+1);
01946                     mesh.cellZones().set
01947                     (
01948                         zoneI,
01949                         new cellZone
01950                         (
01951                             zoneName,           //name
01952                             regionCells,        //addressing
01953                             zoneI,              //index
01954                             cellZones           //cellZoneMesh
01955                         )
01956                     );
01957                 }
01958                 else
01959                 {
01960                     mesh.cellZones()[zoneI].clearAddressing();
01961                     mesh.cellZones()[zoneI] = regionCells;
01962                 }
01963                 Info<< "    Region " << regionI << " : created new cellZone "
01964                     << zoneI << ' ' << cellZones[zoneI].name() << endl;
01965             }
01966         }
01967         mesh.cellZones().writeOpt() = IOobject::AUTO_WRITE;
01968 
01969         if (!overwrite)
01970         {
01971             runTime++;
01972             mesh.setInstance(runTime.timeName());
01973         }
01974         else
01975         {
01976             mesh.setInstance(oldInstance);
01977         }
01978 
01979         Info<< "Writing cellZones as new mesh to time " << runTime.timeName()
01980             << nl << endl;
01981 
01982         mesh.write();
01983 
01984 
01985         // Write cellSets for convenience
01986         // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
01987 
01988         Info<< "Writing cellSets corresponding to cellZones." << nl << endl;
01989 
01990         forAll(cellZones, zoneI)
01991         {
01992             const cellZone& cz = cellZones[zoneI];
01993 
01994             cellSet(mesh, cz.name(), cz).write();
01995         }
01996     }
01997     else
01998     {
01999         // Add patches for interfaces
02000         // ~~~~~~~~~~~~~~~~~~~~~~~~~~
02001 
02002         // Add all possible patches. Empty ones get filtered later on.
02003         Info<< nl << "Adding patches" << nl << endl;
02004 
02005         EdgeMap<label> interfaceToPatch
02006         (
02007             addRegionPatches
02008             (
02009                 mesh,
02010                 cellRegion,
02011                 nCellRegions,
02012                 interfaces,
02013                 interfaceSizes,
02014                 regionNames
02015             )
02016         );
02017 
02018 
02019         if (!overwrite)
02020         {
02021             runTime++;
02022         }
02023 
02024 
02025         // Create regions
02026         // ~~~~~~~~~~~~~~
02027 
02028         if (insidePoint)
02029         {
02030             point insidePoint(args.optionLookup("insidePoint")());
02031 
02032             label regionI = -1;
02033 
02034             label cellI = mesh.findCell(insidePoint);
02035 
02036             Info<< nl << "Found point " << insidePoint << " in cell " << cellI
02037                 << endl;
02038 
02039             if (cellI != -1)
02040             {
02041                 regionI = cellRegion[cellI];
02042             }
02043 
02044             reduce(regionI, maxOp<label>());
02045 
02046             Info<< nl
02047                 << "Subsetting region " << regionI
02048                 << " containing point " << insidePoint << endl;
02049 
02050             if (regionI == -1)
02051             {
02052                 FatalErrorIn(args.executable())
02053                     << "Point " << insidePoint
02054                     << " is not inside the mesh." << nl
02055                     << "Bounding box of the mesh:" << mesh.bounds()
02056                     << exit(FatalError);
02057             }
02058 
02059             createAndWriteRegion
02060             (
02061                 mesh,
02062                 cellRegion,
02063                 regionNames,
02064                 interfaceToPatch,
02065                 regionI,
02066                 (overwrite ? oldInstance : runTime.timeName())
02067             );
02068         }
02069         else if (largestOnly)
02070         {
02071             label regionI = findMax(regionSizes);
02072 
02073             Info<< nl
02074                 << "Subsetting region " << regionI
02075                 << " of size " << regionSizes[regionI] << endl;
02076 
02077             createAndWriteRegion
02078             (
02079                 mesh,
02080                 cellRegion,
02081                 regionNames,
02082                 interfaceToPatch,
02083                 regionI,
02084                 (overwrite ? oldInstance : runTime.timeName())
02085             );
02086         }
02087         else
02088         {
02089             // Split all
02090             for (label regionI = 0; regionI < nCellRegions; regionI++)
02091             {
02092                 Info<< nl
02093                     << "Region " << regionI << nl
02094                     << "-------- " << endl;
02095 
02096                 createAndWriteRegion
02097                 (
02098                     mesh,
02099                     cellRegion,
02100                     regionNames,
02101                     interfaceToPatch,
02102                     regionI,
02103                     (overwrite ? oldInstance : runTime.timeName())
02104                 );
02105             }
02106         }
02107     }
02108 
02109     Info<< "End\n" << endl;
02110 
02111     return 0;
02112 }
02113 
02114 
02115 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines