/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation \\/ M anipulation | Copyright (C) 2016-2018 OpenCFD Ltd. ------------------------------------------------------------------------------- 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 3 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, see . Application subsetMesh Group grpMeshManipulationUtilities Description Create a mesh subset for a particular region of interest based on a cellSet or cellZone. See setSet/topoSet utilities on how to define select cells based on various shapes. Will subset all points, faces and cells needed to make a sub-mesh, but not preserve attached boundary types. \*---------------------------------------------------------------------------*/ #include "fvMeshSubset.H" #include "argList.H" #include "IOobjectList.H" #include "volFields.H" #include "topoDistanceData.H" #include "FaceCellWave.H" #include "cellSet.H" #include "faceSet.H" #include "pointSet.H" #include "ReadFields.H" #include "processorMeshes.H" using namespace Foam; // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // Get the exposed patchId or define the exposedPatchName in fvMeshSubset label getExposedPatchId(const polyMesh& mesh, const word& patchName) { const label patchId = mesh.boundaryMesh().findPatchID(patchName); if (patchId == -1) { fvMeshSubset::exposedPatchName = patchName; } Info<< "Adding exposed internal faces to " << (patchId == -1 ? "new" : "existing") << " patch \"" << patchName << "\"" << nl << endl; return patchId; } labelList nearestPatch(const polyMesh& mesh, const labelList& patchIDs) { const polyBoundaryMesh& pbm = mesh.boundaryMesh(); // Count number of faces in exposedPatchIDs label nFaces = 0; for (const label patchi : patchIDs) { nFaces += pbm[patchi].size(); } // Field on cells and faces. List cellData(mesh.nCells()); List faceData(mesh.nFaces()); // Start of changes labelList patchFaces(nFaces); List patchData(nFaces); nFaces = 0; for (const label patchi : patchIDs) { const polyPatch& pp = pbm[patchi]; forAll(pp, i) { patchFaces[nFaces] = pp.start()+i; patchData[nFaces] = topoDistanceData(patchi, 0); ++nFaces; } } // Propagate information inwards FaceCellWave deltaCalc ( mesh, patchFaces, patchData, faceData, cellData, mesh.globalData().nTotalCells()+1 ); // And extract labelList nearest(mesh.nFaces()); bool haveWarned = false; forAll(faceData, faceI) { if (!faceData[faceI].valid(deltaCalc.data())) { if (!haveWarned) { WarningInFunction << "Did not visit some faces, e.g. face " << faceI << " at " << mesh.faceCentres()[faceI] << nl << "Using patch " << patchIDs[0] << " as nearest" << endl; haveWarned = true; } nearest[faceI] = patchIDs[0]; } else { nearest[faceI] = faceData[faceI].data(); } } return nearest; } // // Subset field-type, availability information cached // in the availableFields hashtable. // template class PatchField, class GeoMesh> void subsetFields ( const fvMeshSubset& subsetter, HashTable& availableFields, PtrList>& subFields ) { typedef GeometricField FieldType; const word fieldType = FieldType::typeName; const wordList fieldNames = availableFields(fieldType).sortedToc(); subFields.setSize(fieldNames.size()); const fvMesh& baseMesh = subsetter.baseMesh(); label nFields = 0; for (const word& fieldName : fieldNames) { if (!nFields) { Info<< "Subsetting " << fieldType << " ("; } else { Info<< ' '; } Info<< fieldName; FieldType fld ( IOobject ( fieldName, baseMesh.time().timeName(), baseMesh, IOobject::MUST_READ, IOobject::NO_WRITE ), baseMesh ); subFields.set(nFields, subsetter.interpolate(fld)); // Subsetting adds 'subset' prefix - rename to match original. subFields[nFields].rename(fieldName); ++nFields; } if (nFields) { Info<< ')' << nl; } } template void subsetPointFields ( const fvMeshSubset& subsetter, const pointMesh& pMesh, HashTable& availableFields, PtrList>& subFields ) { typedef GeometricField FieldType; const word fieldType = FieldType::typeName; const wordList fieldNames = availableFields(fieldType).sortedToc(); subFields.setSize(fieldNames.size()); const fvMesh& baseMesh = subsetter.baseMesh(); label nFields = 0; for (const word& fieldName : fieldNames) { if (!nFields) { Info<< "Subsetting " << fieldType << " ("; } else { Info<< ' '; } Info<< fieldName; FieldType fld ( IOobject ( fieldName, baseMesh.time().timeName(), baseMesh, IOobject::MUST_READ, IOobject::NO_WRITE ), pMesh ); subFields.set(nFields, subsetter.interpolate(fld)); // Subsetting adds 'subset' prefix - rename to match original. subFields[nFields].rename(fieldName); ++nFields; } if (nFields) { Info<< ')' << nl; } } template void subsetDimensionedFields ( const fvMeshSubset& subsetter, HashTable& availableFields, PtrList>& subFields ) { typedef DimensionedField FieldType; const word fieldType = FieldType::typeName; const wordList fieldNames = availableFields(fieldType).sortedToc(); subFields.setSize(fieldNames.size()); const fvMesh& baseMesh = subsetter.baseMesh(); label nFields = 0; for (const word& fieldName : fieldNames) { if (!nFields) { Info<< "Subsetting " << fieldType << " ("; } else { Info<< ' '; } Info<< fieldName; FieldType fld ( IOobject ( fieldName, baseMesh.time().timeName(), baseMesh, IOobject::MUST_READ, IOobject::NO_WRITE ), baseMesh ); subFields.set(nFields, subsetter.interpolate(fld)); // Subsetting adds 'subset' prefix - rename to match original. subFields[nFields].rename(fieldName); ++nFields; } if (nFields) { Info<< ')' << nl; } } template void subsetTopoSets ( const fvMesh& mesh, const IOobjectList& objects, const labelList& map, const fvMesh& subMesh, PtrList& subSets ) { // Read original sets PtrList sets; ReadFields(objects, sets); subSets.setSize(sets.size()); forAll(sets, seti) { const TopoSet& set = sets[seti]; Info<< "Subsetting " << set.type() << " " << set.name() << endl; labelHashSet subset(2*min(set.size(), map.size())); forAll(map, i) { if (set.found(map[i])) { subset.insert(i); } } subSets.set ( seti, new TopoSet ( subMesh, set.name(), std::move(subset), IOobject::AUTO_WRITE ) ); } } int main(int argc, char *argv[]) { argList::addNote ( "Create a mesh subset for a particular region of interest based on a" " cellSet or cellZone(s) specified as the first command argument.\n" "See setSet/topoSet utilities on how to select cells based on" " various shapes." ); #include "addOverwriteOption.H" #include "addRegionOption.H" argList::addArgument("setOrZoneName"); argList::addOption ( "patch", "name", "Add exposed internal faces to specified patch " "instead of \"oldInternalFaces\"" ); argList::addOption ( "patches", "wordRes", "Add exposed internal faces to closest of specified patches " "instead of \"oldInternalFaces\"" ); argList::addBoolOption ( "zone", "Subset with cellZone(s) instead of cellSet. " "The command argument may be a list of words or regexs" ); argList::addOption ( "resultTime", "time", "Specify a time for the resulting mesh" ); argList::noFunctionObjects(); // Never use function objects #include "setRootCase.H" #include "createTime.H" #include "createNamedMesh.H" // arg[1] = word (cellSet) or wordRes (cellZone) // const word selectionName = args[1]; word meshInstance = mesh.pointsInstance(); word fieldsInstance = runTime.timeName(); const bool useCellZone = args.found("zone"); const bool overwrite = args.found("overwrite"); const bool specifiedInstance = args.readIfPresent ( "resultTime", fieldsInstance ); if (specifiedInstance) { // Set both mesh and field to this time meshInstance = fieldsInstance; } // Default exposed patch id labelList exposedPatchIDs(one(), -1); if (args.found("patches")) { const wordRes patchNames(args.getList("patches")); if (patchNames.size() == 1 && patchNames.first().isLiteral()) { exposedPatchIDs.first() = getExposedPatchId(mesh, patchNames.first()); } else { exposedPatchIDs = mesh.boundaryMesh().patchSet(patchNames).sortedToc(); Info<< "Adding exposed internal faces to nearest of patches " << flatOutput(patchNames) << nl << endl; if (exposedPatchIDs.empty()) { FatalErrorInFunction << nl << "No patches matched. Patches: " << mesh.boundaryMesh().names() << nl << exit(FatalError); } } } else if (args.found("patch")) { exposedPatchIDs.first() = getExposedPatchId(mesh, args.opt("patch")); } else { Info<< "Adding exposed internal faces to patch \"" << fvMeshSubset::exposedPatchName << "\" (created if necessary)" << nl << nl; } autoPtr cellSetPtr; // arg[1] can be a word (for cellSet) or wordRes (for cellZone) wordRes zoneNames; if (useCellZone) { List selectionNames = args.getList(1); zoneNames.transfer(selectionNames); Info<< "Using cellZone " << flatOutput(zoneNames) << nl << endl; if (mesh.cellZones().findIndex(zoneNames) == -1) { FatalErrorInFunction << "No cellZones found: " << flatOutput(zoneNames) << nl << nl << exit(FatalError); } } else { const word selectionName = args[1]; Info<< "Using cellSet " << selectionName << nl << endl; cellSetPtr = autoPtr::New(mesh, selectionName); } // Mesh subsetting engine fvMeshSubset subsetter(mesh); { bitSet selectedCells = ( cellSetPtr ? BitSetOps::create(mesh.nCells(), *cellSetPtr) : mesh.cellZones().selection(zoneNames) ); if (exposedPatchIDs.size() == 1) { // Single patch for exposed faces subsetter.setCellSubset ( selectedCells, exposedPatchIDs.first(), true ); } else { // The nearest patch per face labelList nearestExposedPatch(nearestPatch(mesh, exposedPatchIDs)); labelList exposedFaces ( subsetter.getExposedFaces(selectedCells, true) ); subsetter.setCellSubset ( selectedCells, exposedFaces, labelUIndList(nearestExposedPatch, exposedFaces)(), true ); } Info<< "Subset " << returnReduce(subsetter.subMesh().nCells(), sumOp