ENH: subsetMesh: have -cellSet and -dict argument

This commit is contained in:
mattijs
2012-11-23 11:27:46 +00:00
parent 4909ff35ba
commit b4bb215a8c
16 changed files with 1189 additions and 21 deletions

View File

@ -1,3 +1,6 @@
cellSelection/cellSelection.C
cellSelection/badQualityCellSelection.C
cellSelection/outsideCellSelection.C
subsetMesh.C subsetMesh.C
EXE = $(FOAM_APPBIN)/subsetMesh EXE = $(FOAM_APPBIN)/subsetMesh

View File

@ -1,8 +1,11 @@
EXE_INC = \ EXE_INC = \
-IcellSelection \
-I$(LIB_SRC)/finiteVolume/lnInclude \ -I$(LIB_SRC)/finiteVolume/lnInclude \
-I$(LIB_SRC)/dynamicMesh/lnInclude \
-I$(LIB_SRC)/meshTools/lnInclude -I$(LIB_SRC)/meshTools/lnInclude
EXE_LIBS = \ EXE_LIBS = \
-lfiniteVolume \ -lfiniteVolume \
-ldynamicMesh \
-lmeshTools \ -lmeshTools \
-lgenericPatchFields -lgenericPatchFields

View File

@ -0,0 +1,95 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012 OpenFOAM Foundation
\\/ 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 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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "badQualityCellSelection.H"
#include "addToRunTimeSelectionTable.H"
#include "faceSet.H"
#include "polyMesh.H"
#include "motionSmoother.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace cellSelections
{
defineTypeNameAndDebug(badQualityCellSelection, 0);
addToRunTimeSelectionTable
(
cellSelection,
badQualityCellSelection,
dictionary
);
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::cellSelections::badQualityCellSelection::badQualityCellSelection
(
const word& name,
const polyMesh& mesh,
const dictionary& dict
)
:
cellSelection(name, mesh, dict)
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::cellSelections::badQualityCellSelection::~badQualityCellSelection()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::cellSelections::badQualityCellSelection::select
(
boolList& selectedCell
) const
{
//- Delete cell of any face in error
faceSet faces(mesh_, "meshQualityFaces", mesh_.nFaces()/100+1);
motionSmoother::checkMesh(false, mesh_, dict_, faces);
label nFaces = returnReduce(faces.size(), sumOp<label>());
if (nFaces > 0)
{
faces.sync(mesh_);
forAllConstIter(faceSet, faces, iter)
{
label faceI = iter.key();
selectedCell[mesh_.faceOwner()[faceI]] = false;
if (mesh_.isInternalFace(faceI))
{
selectedCell[mesh_.faceNeighbour()[faceI]] = false;
}
}
}
}
// ************************************************************************* //

View File

@ -0,0 +1,101 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012 OpenFOAM Foundation
\\/ 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 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 <http://www.gnu.org/licenses/>.
Class
Foam::cellSelections::badQualityCellSelection
Description
Deselect bad quality cells
SourceFiles
badQualityCellSelection.C
\*---------------------------------------------------------------------------*/
#ifndef badQualityCellSelection_H
#define badQualityCellSelection_H
#include "cellSelection.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace cellSelections
{
/*---------------------------------------------------------------------------*\
Class badQualityCellSelection Declaration
\*---------------------------------------------------------------------------*/
class badQualityCellSelection
:
public cellSelection
{
public:
//- Runtime type information
TypeName("badQuality");
// Constructors
//- Construct from dictionary
badQualityCellSelection
(
const word& name,
const polyMesh& mesh,
const dictionary& dict
);
//- Clone
autoPtr<cellSelection> clone() const
{
notImplemented("autoPtr<cellSelection> clone() const");
return autoPtr<cellSelection>(NULL);
}
//- Destructor
virtual ~badQualityCellSelection();
// Member Functions
virtual void select(boolList&) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace cellSelections
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,104 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012 OpenFOAM Foundation
\\/ 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 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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "cellSelection.H"
#include "polyMesh.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
defineTypeNameAndDebug(cellSelection, 0);
defineRunTimeSelectionTable(cellSelection, dictionary);
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::cellSelection::cellSelection
(
const word& name,
const polyMesh& mesh,
const dictionary& dict
)
:
name_(name),
mesh_(mesh),
dict_(dict)
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::cellSelection::~cellSelection()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::autoPtr<Foam::cellSelection> Foam::cellSelection::New
(
const word& name,
const polyMesh& mesh,
const dictionary& dict
)
{
const word sampleType(dict.lookup("type"));
dictionaryConstructorTable::iterator cstrIter =
dictionaryConstructorTablePtr_->find(sampleType);
if (cstrIter == dictionaryConstructorTablePtr_->end())
{
FatalErrorIn
(
"cellSelection::New"
"(const word&, const polyMesh&, const dictionary&)"
) << "Unknown cellSelection type "
<< sampleType << nl << nl
<< "Valid cellSelection types : " << endl
<< dictionaryConstructorTablePtr_->sortedToc()
<< exit(FatalError);
}
return autoPtr<cellSelection>(cstrIter()(name, mesh, dict));
}
Foam::label Foam::cellSelection::count(const boolList& lst)
{
label n = 0;
forAll(lst, i)
{
if (lst[i])
{
n++;
}
}
return returnReduce(n, sumOp<label>());
}
// ************************************************************************* //

View File

@ -0,0 +1,149 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012 OpenFOAM Foundation
\\/ 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 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 <http://www.gnu.org/licenses/>.
Class
Foam::cellSelection
Description
Cell selection methods in subsetMesh
SourceFiles
cellSelection.C
\*---------------------------------------------------------------------------*/
#ifndef cellSelection_H
#define cellSelection_H
#include "dictionary.H"
#include "typeInfo.H"
#include "runTimeSelectionTables.H"
#include "autoPtr.H"
#include "boolList.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
// Forward declaration of classes
class polyMesh;
/*---------------------------------------------------------------------------*\
Class cellSelection Declaration
\*---------------------------------------------------------------------------*/
class cellSelection
{
protected:
// Protected data
//- Name
const word name_;
//- Reference to mesh
const polyMesh& mesh_;
//- Input dictionary
const dictionary dict_;
public:
//- Runtime type information
TypeName("cellSelection");
// Declare run-time constructor selection table
declareRunTimeSelectionTable
(
autoPtr,
cellSelection,
dictionary,
(
const word& name,
const polyMesh& mesh,
const dictionary& dict
),
(name, mesh, dict)
);
// Constructors
//- Construct from dictionary
cellSelection
(
const word& name,
const polyMesh& mesh,
const dictionary& dict
);
//- Clone
autoPtr<cellSelection> clone() const
{
notImplemented("autoPtr<cellSelection> clone() const");
return autoPtr<cellSelection>(NULL);
}
// Selectors
//- Return a reference to the selected cellSelection
static autoPtr<cellSelection> New
(
const word& name,
const polyMesh& mesh,
const dictionary& dict
);
//- Destructor
virtual ~cellSelection();
// Member Functions
//- Count global number of selected elements
static label count(const boolList&);
const word& name() const
{
return name_;
}
virtual void select(boolList&) const = 0;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,406 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012 OpenFOAM Foundation
\\/ 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 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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "outsideCellSelection.H"
#include "addToRunTimeSelectionTable.H"
#include "faceSet.H"
#include "polyMesh.H"
#include "motionSmoother.H"
#include "regionSplit.H"
#include "syncTools.H"
#include "zeroGradientFvPatchFields.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace cellSelections
{
defineTypeNameAndDebug(outsideCellSelection, 0);
addToRunTimeSelectionTable(cellSelection, outsideCellSelection, dictionary);
}
}
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
Foam::tmp<Foam::volScalarField>
Foam::cellSelections::outsideCellSelection::generateField
(
const word& name,
const boolList& lst
) const
{
const fvMesh& mesh = dynamic_cast<const fvMesh&>(mesh_);
tmp<volScalarField> tfld
(
new volScalarField
(
IOobject
(
name,
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
mesh,
dimensionedScalar(name, dimless, 0),
zeroGradientFvPatchScalarField::typeName
)
);
scalarField& fld = tfld().internalField();
forAll(fld, celli)
{
fld[celli] = 1.0*lst[celli];
}
tfld().correctBoundaryConditions();
return tfld;
}
void Foam::cellSelections::outsideCellSelection::markRegionFaces
(
const boolList& selectedCell,
boolList& regionFace
) const
{
// Internal faces
const labelList& faceOwner = mesh_.faceOwner();
const labelList& faceNeighbour = mesh_.faceNeighbour();
forAll(faceNeighbour, faceI)
{
if
(
selectedCell[faceOwner[faceI]]
!= selectedCell[faceNeighbour[faceI]]
)
{
regionFace[faceI] = true;
}
}
// Swap neighbour selectedCell state
boolList nbrSelected;
syncTools::swapBoundaryCellList(mesh_, selectedCell, nbrSelected);
// Boundary faces
const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
forAll(pbm, patchI)
{
const polyPatch& pp = pbm[patchI];
const labelUList& faceCells = pp.faceCells();
forAll(faceCells, i)
{
label faceI = pp.start()+i;
label bFaceI = faceI-mesh_.nInternalFaces();
if
(
selectedCell[faceCells[i]]
!= selectedCell[nbrSelected[bFaceI]]
)
{
regionFace[faceI] = true;
}
}
}
}
Foam::boolList Foam::cellSelections::outsideCellSelection::findRegions
(
const bool verbose,
const regionSplit& cellRegion
) const
{
boolList keepRegion(cellRegion.nRegions(), false);
forAll(locationsInMesh_, i)
{
// Find the region containing the insidePoint
label cellI = mesh_.findCell(locationsInMesh_[i]);
label keepRegionI = -1;
label keepProcI = -1;
if (cellI != -1)
{
keepRegionI = cellRegion[cellI];
keepProcI = Pstream::myProcNo();
}
reduce(keepRegionI, maxOp<label>());
keepRegion[keepRegionI] = true;
if (verbose)
{
reduce(keepProcI, maxOp<label>());
Info<< "Found location " << locationsInMesh_[i]
<< " in cell " << cellI << " on processor " << keepProcI
<< " in global region " << keepRegionI
<< " out of " << cellRegion.nRegions() << " regions." << endl;
}
}
return keepRegion;
}
void Foam::cellSelections::outsideCellSelection::unselectOutsideRegions
(
boolList& selectedCell
) const
{
// Determine faces on the edge of selectedCell
boolList blockedFace(mesh_.nFaces(), false);
markRegionFaces(selectedCell, blockedFace);
// Determine regions
regionSplit cellRegion(mesh_, blockedFace);
// Determine regions containing locationsInMesh_
boolList keepRegion(findRegions(true, cellRegion));
// Go back to bool per cell
forAll(cellRegion, cellI)
{
if (!keepRegion[cellRegion[cellI]])
{
selectedCell[cellI] = false;
}
}
}
void Foam::cellSelections::outsideCellSelection::shrinkRegions
(
boolList& selectedCell
) const
{
// Select points on unselected cells and boundary
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
boolList boundaryPoint(mesh_.nPoints(), false);
const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
forAll(pbm, patchI)
{
const polyPatch& pp = pbm[patchI];
if (!pp.coupled() && !isA<emptyPolyPatch>(pp))
{
forAll(pp, i)
{
const face& f = pp[i];
forAll(f, fp)
{
boundaryPoint[f[fp]] = true;
}
}
}
}
forAll(selectedCell, cellI)
{
if (!selectedCell[cellI])
{
const labelList& cPoints = mesh_.cellPoints(cellI);
forAll(cPoints, i)
{
boundaryPoint[cPoints[i]] = true;
}
}
}
syncTools::syncPointList(mesh_, boundaryPoint, orEqOp<bool>(), false);
// Select all cells using these points
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
label nChanged = 0;
forAll(boundaryPoint, pointI)
{
if (boundaryPoint[pointI])
{
const labelList& pCells = mesh_.pointCells(pointI);
forAll(pCells, i)
{
label cellI = pCells[i];
if (selectedCell[cellI])
{
selectedCell[cellI] = false;
nChanged++;
}
}
}
}
}
void Foam::cellSelections::outsideCellSelection::erode
(
boolList& selectedCell
) const
{
//Info<< "Entering shrinkRegions:" << count(selectedCell) << endl;
//generateField("selectedCell_before", selectedCell)().write();
// Now erode and see which regions get disconnected
boolList shrunkSelectedCell(selectedCell);
for (label iter = 0; iter < nErode_; iter++)
{
shrinkRegions(shrunkSelectedCell);
}
//Info<< "After shrinking:" << count(shrunkSelectedCell) << endl;
//generateField("shrunkSelectedCell", shrunkSelectedCell)().write();
// Determine faces on the edge of shrunkSelectedCell
boolList blockedFace(mesh_.nFaces(), false);
markRegionFaces(shrunkSelectedCell, blockedFace);
// Find disconnected regions
regionSplit cellRegion(mesh_, blockedFace);
// Determine regions containing insidePoints
boolList keepRegion(findRegions(true, cellRegion));
// Extract cells in regions that are not to be kept.
boolList removeCell(mesh_.nCells(), false);
forAll(cellRegion, cellI)
{
if (shrunkSelectedCell[cellI] && !keepRegion[cellRegion[cellI]])
{
removeCell[cellI] = true;
}
}
//Info<< "removeCell before:" << count(removeCell) << endl;
//generateField("removeCell_before", removeCell)().write();
// Grow removeCell
for (label iter = 0; iter < nErode_; iter++)
{
// Grow selected cell in regions that are not for keeping
boolList boundaryPoint(mesh_.nPoints(), false);
forAll(removeCell, cellI)
{
if (removeCell[cellI])
{
const labelList& cPoints = mesh_.cellPoints(cellI);
forAll(cPoints, i)
{
boundaryPoint[cPoints[i]] = true;
}
}
}
syncTools::syncPointList(mesh_, boundaryPoint, orEqOp<bool>(), false);
// Select all cells using these points
label nChanged = 0;
forAll(boundaryPoint, pointI)
{
if (boundaryPoint[pointI])
{
const labelList& pCells = mesh_.pointCells(pointI);
forAll(pCells, i)
{
label cellI = pCells[i];
if (!removeCell[cellI])
{
removeCell[cellI] = true;
nChanged++;
}
}
}
}
}
//Info<< "removeCell after:" << count(removeCell) << endl;
//generateField("removeCell_after", removeCell)().write();
// Unmark removeCell
forAll(removeCell, cellI)
{
if (removeCell[cellI])
{
selectedCell[cellI] = false;
}
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::cellSelections::outsideCellSelection::outsideCellSelection
(
const word& name,
const polyMesh& mesh,
const dictionary& dict
)
:
cellSelection(name, mesh, dict),
locationsInMesh_(dict.lookup("locationsInMesh")),
nErode_(readLabel(dict.lookup("nErodeLayers")))
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::cellSelections::outsideCellSelection::~outsideCellSelection()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::cellSelections::outsideCellSelection::select
(
boolList& selectedCell
) const
{
// Unselect all disconnected regions
unselectOutsideRegions(selectedCell);
if (nErode_ > 0)
{
erode(selectedCell);
}
}
// ************************************************************************* //

View File

@ -0,0 +1,145 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012 OpenFOAM Foundation
\\/ 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 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 <http://www.gnu.org/licenses/>.
Class
Foam::cellSelections::outsideCellSelection
Description
Deselect cells not reachable from 'inside' points
SourceFiles
outsideCellSelection.C
\*---------------------------------------------------------------------------*/
#ifndef outsideCellSelection_H
#define outsideCellSelection_H
#include "cellSelection.H"
#include "pointField.H"
#include "boolList.H"
#include "volFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
class regionSplit;
namespace cellSelections
{
/*---------------------------------------------------------------------------*\
Class outsideCellSelection Declaration
\*---------------------------------------------------------------------------*/
class outsideCellSelection
:
public cellSelection
{
// Private data
//- Locations to keep
const pointField locationsInMesh_;
//- Number of layers to erode
const label nErode_;
// Private Member Functions
//- For debugging: generate volScalarField with 1.0 for all true
tmp<volScalarField> generateField
(
const word& name,
const boolList& lst
) const;
//- Mark faces inbetween selected and unselected elements
void markRegionFaces
(
const boolList& selectedCell,
boolList& regionFace
) const;
//- Determine for every disconnected region in the mesh whether
// it contains a locationInMesh
boolList findRegions(const bool verbose, const regionSplit&) const;
//- Unselect regions not containing a locationInMesh
void unselectOutsideRegions(boolList& selectedCell) const;
//- Unselect one layer of cells from selectedCell
void shrinkRegions(boolList& selectedCell) const;
//- Erode a given number of layers from selectedCell. Remove any
// region that gets disconnected that way.
void erode(boolList& selectedCell) const;
public:
//- Runtime type information
TypeName("outside");
// Constructors
//- Construct from dictionary
outsideCellSelection
(
const word& name,
const polyMesh& mesh,
const dictionary& dict
);
//- Clone
autoPtr<cellSelection> clone() const
{
notImplemented("autoPtr<cellSelection> clone() const");
return autoPtr<cellSelection>(NULL);
}
//- Destructor
virtual ~outsideCellSelection();
// Member Functions
//- Apply this selector
virtual void select(boolList&) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace cellSelections
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -2,7 +2,7 @@
========= | ========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | \\ / O peration |
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation \\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
@ -37,12 +37,12 @@ Description
#include "cellSet.H" #include "cellSet.H"
#include "IOobjectList.H" #include "IOobjectList.H"
#include "volFields.H" #include "volFields.H"
#include "cellSelection.H"
using namespace Foam; using namespace Foam;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
template<class Type> template<class Type>
void subsetVolFields void subsetVolFields
( (
@ -152,12 +152,21 @@ int main(int argc, char *argv[])
{ {
argList::addNote argList::addNote
( (
"select a mesh subset based on a cellSet" "select a mesh subset based on a provided cellSet and/or"
" selection criteria"
);
argList::addBoolOption
(
"dict",
"read mesh subset selection criteria"
" from system/subsetMeshDict"
);
argList::addOption
(
"cellSet",
"name",
"operates on specified cellSet name"
); );
#include "addOverwriteOption.H"
#include "addRegionOption.H"
argList::validArgs.append("cellSet");
argList::addOption argList::addOption
( (
"patch", "patch",
@ -165,23 +174,81 @@ int main(int argc, char *argv[])
"add exposed internal faces to specified patch instead of to " "add exposed internal faces to specified patch instead of to "
"'oldInternalFaces'" "'oldInternalFaces'"
); );
#include "addOverwriteOption.H"
#include "addRegionOption.H"
#include "setRootCase.H" #include "setRootCase.H"
#include "createTime.H" #include "createTime.H"
runTime.functionObjects().off(); runTime.functionObjects().off();
Foam::word meshRegionName = polyMesh::defaultRegion; word setName;
args.optionReadIfPresent("region", meshRegionName); const bool useCellSet = args.optionReadIfPresent("cellSet", setName);
const bool useDict = args.optionFound("dict");
const bool overwrite = args.optionFound("overwrite");
if (!useCellSet && !useDict)
{
FatalErrorIn(args.executable())
<< "No cells to operate on selected. Please supply at least one of "
<< "'-cellSet', '-dict'"
<< exit(FatalError);
}
#include "createNamedMesh.H" #include "createNamedMesh.H"
const word oldInstance = mesh.pointsInstance(); const word oldInstance = mesh.pointsInstance();
const word setName = args[1]; autoPtr<cellSet> currentSet;
const bool overwrite = args.optionFound("overwrite"); if (useCellSet)
{
// Load the cellSet
Info<< "Operating on cell set " << setName << nl << endl;
currentSet.reset(new cellSet(mesh, setName));
}
PtrList<cellSelection> selectors;
if (useDict)
{
Info<< "Reading selection criteria from subsetMeshDict" << nl << endl;
IOdictionary dict
(
IOobject
(
"subsetMeshDict",
mesh.time().system(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE
)
);
const dictionary& selectionsDict = dict.subDict("selections");
label n = 0;
forAllConstIter(dictionary, selectionsDict, iter)
{
if (iter().isDict())
{
n++;
}
}
selectors.setSize(n);
n = 0;
forAllConstIter(dictionary, selectionsDict, iter)
{
if (iter().isDict())
{
selectors.set
(
n++,
cellSelection::New(iter().keyword(), mesh, iter().dict())
);
}
}
}
Info<< "Reading cell set from " << setName << endl << endl;
// Create mesh subsetting engine // Create mesh subsetting engine
fvMeshSubset subsetter(mesh); fvMeshSubset subsetter(mesh);
@ -212,9 +279,54 @@ int main(int argc, char *argv[])
} }
cellSet currentSet(mesh, setName); // Select cells to operate on
// ~~~~~~~~~~~~~~~~~~~~~~~~~~
boolList selectedCell(mesh.nCells(), false);
if (currentSet.valid())
{
const cellSet& set = currentSet();
forAllConstIter(cellSet, set, iter)
{
selectedCell[iter.key()] = true;
}
}
else
{
selectedCell = true;
}
// Manipulate selectedCell according to dictionary
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
forAll(selectors, i)
{
Info<< "Applying cellSelector " << selectors[i].name() << endl;
selectors[i].select(selectedCell);
Info<< "After applying cellSelector " << selectors[i].name()
<< " have " << cellSelection::count(selectedCell)
<< " cells" << nl << endl;
}
// Back from selectedCells to region
{
Info<< "Final selection : " << cellSelection::count(selectedCell)
<< " cells" << nl << endl;
labelList cellRegion(mesh.nCells(), -1);
forAll(selectedCell, cellI)
{
if (selectedCell[cellI])
{
cellRegion[cellI] = 0;
}
}
subsetter.setLargeCellSubset(cellRegion, 0, patchI, true);
}
subsetter.setLargeCellSubset(currentSet, patchI, true);
IOobjectList objects(mesh, runTime.timeName()); IOobjectList objects(mesh, runTime.timeName());

View File

@ -0,0 +1,50 @@
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: dev |
| \\ / A nd | Web: http://www.openfoam.org |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "system";
object subsetMeshDict;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
//- Load snappyHexMeshDict settings so we can use some entries below.
#include "./snappyHexMeshDict"
selections
{
badQuality
{
// Remove any cells using a 'bad quality' face. Uses the mesh checks
// used by snappyHexMesh, cvMesh.
type badQuality;
// Use the quality criteria from the snappyHexMeshDict
${.meshQualityControls}
}
outside
{
// Remove any cells not reachable from provided locations
type outside;
//- Number of cell layers to erode mesh to detect holes in the mesh
// Set to 0 if not used.
nErodeLayers 3;
//- Define mesh location for keeping.
// In this case use the one from snappyHexMeshDict. This can
// optionally be a list of locations.
locationsInMesh (${.castellatedMeshControls.locationInMesh});
}
}
// ************************************************************************* //

View File

@ -9,7 +9,7 @@ application=`getApplication`
runApplication blockMesh runApplication blockMesh
runApplication topoSet runApplication topoSet
runApplication subsetMesh -overwrite c0 -patch floatingObject runApplication subsetMesh -overwrite -cellSet c0 -patch floatingObject
cp -r 0.org 0 > /dev/null 2>&1 cp -r 0.org 0 > /dev/null 2>&1
runApplication $application runApplication $application

View File

@ -13,7 +13,7 @@ runApplication blockMesh
runApplication topoSet runApplication topoSet
# create the obstacles - add obstacle patches to wallFilm patch # create the obstacles - add obstacle patches to wallFilm patch
runApplication subsetMesh c0 -patch wallFilm -overwrite runApplication subsetMesh -cellSet c0 -patch wallFilm -overwrite
# split the obstacle patches into cube[1-6]_patch[1-6] # split the obstacle patches into cube[1-6]_patch[1-6]
echo "running patchifyObstacles" echo "running patchifyObstacles"

View File

@ -14,7 +14,7 @@ runApplication surfaceFeatureExtract
runApplication blockMesh runApplication blockMesh
runApplication topoSet -dict system/topoSetDict-background runApplication topoSet -dict system/topoSetDict-background
mv log.topoSet log.topoSet.background mv log.topoSet log.topoSet.background
runApplication subsetMesh background -patch walls -overwrite runApplication subsetMesh -cellSet background -patch walls -overwrite
runApplication decomposePar runApplication decomposePar

View File

@ -10,7 +10,7 @@ application=`getApplication`
runApplication blockMesh runApplication blockMesh
runApplication topoSet runApplication topoSet
runApplication subsetMesh -overwrite c0 -patch movingBlock runApplication subsetMesh -overwrite -cellSet c0 -patch movingBlock
cp -r 0.org 0 > /dev/null 2>&1 cp -r 0.org 0 > /dev/null 2>&1
runApplication $application runApplication $application
./extractData log.$application ./extractData log.$application

View File

@ -8,7 +8,7 @@ cp -r 0.org 0 > /dev/null 2>&1
runApplication blockMesh runApplication blockMesh
#runApplication setSet -batch createObstacle.setSet #runApplication setSet -batch createObstacle.setSet
runApplication topoSet runApplication topoSet
runApplication subsetMesh -overwrite c0 -patch walls runApplication subsetMesh -overwrite -cellSet c0 -patch walls
runApplication setFields runApplication setFields
runApplication `getApplication` runApplication `getApplication`

View File

@ -9,7 +9,7 @@ application=`getApplication`
runApplication blockMesh runApplication blockMesh
runApplication topoSet runApplication topoSet
runApplication subsetMesh -overwrite c0 -patch floatingObject runApplication subsetMesh -overwrite -cellSet c0 -patch floatingObject
cp -r 0.org 0 > /dev/null 2>&1 cp -r 0.org 0 > /dev/null 2>&1
runApplication setFields runApplication setFields
runApplication $application runApplication $application