/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright held by original author \\/ M anipulation | ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version. OpenFOAM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OpenFOAM; if not, write to the Free Software Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA Application splitMeshRegions Description Splits mesh into multiple regions. Each region is defined as a domain whose cells can all be reached by cell-face-cell walking without crossing - boundary faces - additional faces from faceSet (-blockedFaces faceSet). - any face in between differing cellZones (-cellZones) Output is: - volScalarField with regions as different scalars (-detectOnly) or - mesh with multiple regions or - mesh with cells put into cellZones (-makeCellZones) Note: - cellZonesOnly does not do a walk and uses the cellZones only. Use this if you don't mind having disconnected domains in a single region. This option requires all cells to be in one (and one only) cellZone. - cellZonesFileOnly behaves like -cellZonesOnly but reads the cellZones from the specified file. This allows one to explicitly specify the region distribution and still have multiple cellZones per region. - useCellZonesOnly does not do a walk and uses the cellZones only. Use this if you don't mind having disconnected domains in a single region. This option requires all cells to be in one (and one only) cellZone. - Should work in parallel. cellZones can differ on either side of processor boundaries in which case the faces get moved from processor patch to directMapped patch. Not very well tested. - If a cell zone gets split into more than one region it can detect the largest matching region (-sloppyCellZones). This will accept any region that covers more than 50% of the zone. It has to be a subset so cannot have any cells in any other zone. - writes maps like decomposePar back to original mesh: - pointRegionAddressing : for every point in this region the point in the original mesh - cellRegionAddressing : ,, cell ,, cell ,, - faceRegionAddressing : ,, face ,, face in the original mesh + 'turning index'. For a face in the same orientation this is the original facelabel+1, for a turned face this is -facelabel - 1 \*---------------------------------------------------------------------------*/ #include "SortableList.H" #include "argList.H" #include "regionSplit.H" #include "fvMeshSubset.H" #include "IOobjectList.H" #include "volFields.H" #include "faceSet.H" #include "cellSet.H" #include "directTopoChange.H" #include "mapPolyMesh.H" #include "removeCells.H" #include "EdgeMap.H" #include "syncTools.H" #include "ReadFields.H" #include "directMappedWallPolyPatch.H" #include "zeroGradientFvPatchFields.H" using namespace Foam; // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // template void addPatchFields(fvMesh& mesh, const word& patchFieldType) { HashTable flds ( mesh.objectRegistry::lookupClass() ); for ( typename HashTable::const_iterator iter = flds.begin(); iter != flds.end(); ++iter ) { const GeoField& fld = *iter(); typename GeoField::GeometricBoundaryField& bfld = const_cast ( fld.boundaryField() ); label sz = bfld.size(); bfld.setSize(sz + 1); bfld.set ( sz, GeoField::PatchFieldType::New ( patchFieldType, mesh.boundary()[sz], fld.dimensionedInternalField() ) ); } } // Remove last patch field template void trimPatchFields(fvMesh& mesh, const label nPatches) { HashTable flds ( mesh.objectRegistry::lookupClass() ); for ( typename HashTable::const_iterator iter = flds.begin(); iter != flds.end(); ++iter ) { const GeoField& fld = *iter(); const_cast ( fld.boundaryField() ).setSize(nPatches); } } // Reorder patch field template void reorderPatchFields(fvMesh& mesh, const labelList& oldToNew) { HashTable flds ( mesh.objectRegistry::lookupClass() ); for ( typename HashTable::const_iterator iter = flds.begin(); iter != flds.end(); ++iter ) { const GeoField& fld = *iter(); typename GeoField::GeometricBoundaryField& bfld = const_cast ( fld.boundaryField() ); bfld.reorder(oldToNew); } } // Adds patch if not yet there. Returns patchID. label addPatch(fvMesh& mesh, const polyPatch& patch) { polyBoundaryMesh& polyPatches = const_cast(mesh.boundaryMesh()); label patchI = polyPatches.findPatchID(patch.name()); if (patchI != -1) { if (polyPatches[patchI].type() == patch.type()) { // Already there return patchI; } else { FatalErrorIn("addPatch(fvMesh&, const polyPatch*)") << "Already have patch " << patch.name() << " but of type " << patch.type() << exit(FatalError); } } label insertPatchI = polyPatches.size(); label startFaceI = mesh.nFaces(); forAll (polyPatches, patchI) { const polyPatch& pp = polyPatches[patchI]; if (isA(pp)) { insertPatchI = patchI; startFaceI = pp.start(); break; } } // Below is all quite a hack. Feel free to change once there is a better // mechanism to insert and reorder patches. // Clear local fields and e.g. polyMesh parallelInfo. mesh.clearOut(); label sz = polyPatches.size(); fvBoundaryMesh& fvPatches = const_cast(mesh.boundary()); Info<< "Inserting patch " << patch.name() << " to slot " << insertPatchI << " out of " << polyPatches.size() << endl; // Add polyPatch at the end polyPatches.setSize(sz + 1); polyPatches.set ( sz, patch.clone ( polyPatches, insertPatchI, // index 0, // size startFaceI // start ) ); fvPatches.setSize(sz + 1); fvPatches.set ( sz, fvPatch::New ( polyPatches[sz], // point to newly added polyPatch mesh.boundary() ) ); addPatchFields ( mesh, calculatedFvPatchField::typeName ); addPatchFields ( mesh, calculatedFvPatchField::typeName ); addPatchFields ( mesh, calculatedFvPatchField::typeName ); addPatchFields ( mesh, calculatedFvPatchField::typeName ); addPatchFields ( mesh, calculatedFvPatchField::typeName ); // Surface fields addPatchFields ( mesh, calculatedFvPatchField::typeName ); addPatchFields ( mesh, calculatedFvPatchField::typeName ); addPatchFields ( mesh, calculatedFvPatchField::typeName ); addPatchFields ( mesh, calculatedFvPatchField::typeName ); addPatchFields ( mesh, calculatedFvPatchField::typeName ); // Create reordering list // patches before insert position stay as is labelList oldToNew(sz + 1); for (label i = 0; i < insertPatchI; i++) { oldToNew[i] = i; } // patches after insert position move one up for (label i = insertPatchI; i < sz; i++) { oldToNew[i] = i + 1; } // appended patch gets moved to insert position oldToNew[sz] = insertPatchI; // Shuffle into place polyPatches.reorder(oldToNew); fvPatches.reorder(oldToNew); reorderPatchFields(mesh, oldToNew); reorderPatchFields(mesh, oldToNew); reorderPatchFields(mesh, oldToNew); reorderPatchFields(mesh, oldToNew); reorderPatchFields(mesh, oldToNew); reorderPatchFields(mesh, oldToNew); reorderPatchFields(mesh, oldToNew); reorderPatchFields(mesh, oldToNew); reorderPatchFields(mesh, oldToNew); reorderPatchFields(mesh, oldToNew); return insertPatchI; } // Reorder and delete patches. Deleted. HJ, 22/May/2011 template void subsetVolFields ( const fvMesh& mesh, const fvMesh& subMesh, const labelList& cellMap, const labelList& faceMap, const labelList& patchMap ) { // Real patch map passed by argument // HJ, 22/May/2011 HashTable fields ( mesh.objectRegistry::lookupClass() ); forAllConstIter(typename HashTable, fields, iter) { const GeoField& fld = *iter(); Info<< "Mapping field " << fld.name() << endl; tmp tSubFld ( fvMeshSubset::meshToMesh ( fld, subMesh, patchMap, cellMap, faceMap ) ); // Hack: set value to 0 for introduced patches (since don't // get initialised. forAll (tSubFld().boundaryField(), patchI) { const fvPatchField& pfld = tSubFld().boundaryField()[patchI]; if ( isA > (pfld) ) { tSubFld().boundaryField()[patchI] == pTraits::zero; } } // Store on subMesh GeoField* subFld = tSubFld.ptr(); subFld->rename(fld.name()); subFld->writeOpt() = IOobject::AUTO_WRITE; subFld->store(); } } template void subsetSurfaceFields ( const fvMesh& mesh, const fvMesh& subMesh, const labelList& faceMap, const labelList& patchMap ) { // Real patch map passed by argument // HJ, 22/May/2011 HashTable fields ( mesh.objectRegistry::lookupClass() ); forAllConstIter(typename HashTable, fields, iter) { const GeoField& fld = *iter(); Info<< "Mapping field " << fld.name() << endl; tmp tSubFld ( fvMeshSubset::meshToMesh ( fld, subMesh, patchMap, faceMap ) ); // Hack: set value to 0 for introduced patches (since don't // get initialised. forAll (tSubFld().boundaryField(), patchI) { const fvsPatchField& pfld = tSubFld().boundaryField()[patchI]; if ( isA > (pfld) ) { tSubFld().boundaryField()[patchI] == pTraits::zero; } } // Store on subMesh GeoField* subFld = tSubFld.ptr(); subFld->rename(fld.name()); subFld->writeOpt() = IOobject::AUTO_WRITE; subFld->store(); } } // Select all cells not in the region labelList getNonRegionCells(const labelList& cellRegion, const label regionI) { DynamicList