mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: allow new patch names in subsetMesh (issue #1019)
Previously had 3 possibilities for handling exposed internal faces
1. use default "oldInternalFaces"
2. specify -patch, to use the specified (existing) patch
3. specify -patches, to use the geometrically closest patches
Now relaxed the restriction on -patch to allow specification of a new
(not yet existing) patch name. This improves flexibility, but won't
catch typing mistakes.
Harmonize behaviour of -patches and -patch. When -patches is used to
specify a single, non-regex patch name, it now behaves identically to
-patch. Since the getList handling for options already allows special
treatment for single parameter lists, the following will work
identically:
subsetMesh -patch patch0
subsetMesh -patches patch0
subsetMesh -patches '( patch0 )'
In the future it might be reasonable to fully combine the behaviour of
'-patch' and '-patches' and treat them as aliases for each other.
ENH: support subsetMesh on a cellZone.
- when the '-zone' option is specified, the command argument is treated
as the name (or names) of cellZones to be selected instead of as the
name of the cellSet.
The command argument can be a single word, regex, or list of
word/regex.
Eg,
subsetMesh -zone -patch mypatch mixer
subsetMesh -zone -patch mypatch '(mixer "moving.*" )'
STYLE: simplify set handling and other code cleanup in subsetMesh
This commit is contained in:
@ -45,7 +45,6 @@ Description
|
|||||||
#include "volFields.H"
|
#include "volFields.H"
|
||||||
#include "topoDistanceData.H"
|
#include "topoDistanceData.H"
|
||||||
#include "FaceCellWave.H"
|
#include "FaceCellWave.H"
|
||||||
#include "BitOps.H"
|
|
||||||
#include "cellSet.H"
|
#include "cellSet.H"
|
||||||
#include "faceSet.H"
|
#include "faceSet.H"
|
||||||
#include "pointSet.H"
|
#include "pointSet.H"
|
||||||
@ -56,16 +55,33 @@ 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)
|
labelList nearestPatch(const polyMesh& mesh, const labelList& patchIDs)
|
||||||
{
|
{
|
||||||
const polyBoundaryMesh& pbm = mesh.boundaryMesh();
|
const polyBoundaryMesh& pbm = mesh.boundaryMesh();
|
||||||
|
|
||||||
// Count number of faces in exposedPatchIDs
|
// Count number of faces in exposedPatchIDs
|
||||||
label nFaces = 0;
|
label nFaces = 0;
|
||||||
forAll(patchIDs, i)
|
for (const label patchi : patchIDs)
|
||||||
{
|
{
|
||||||
const polyPatch& pp = pbm[patchIDs[i]];
|
nFaces += pbm[patchi].size();
|
||||||
nFaces += pp.size();
|
|
||||||
}
|
}
|
||||||
|
|
||||||
// Field on cells and faces.
|
// Field on cells and faces.
|
||||||
@ -76,16 +92,15 @@ labelList nearestPatch(const polyMesh& mesh, const labelList& patchIDs)
|
|||||||
labelList patchFaces(nFaces);
|
labelList patchFaces(nFaces);
|
||||||
List<topoDistanceData> patchData(nFaces);
|
List<topoDistanceData> patchData(nFaces);
|
||||||
nFaces = 0;
|
nFaces = 0;
|
||||||
forAll(patchIDs, i)
|
for (const label patchi : patchIDs)
|
||||||
{
|
{
|
||||||
label patchI = patchIDs[i];
|
const polyPatch& pp = pbm[patchi];
|
||||||
const polyPatch& pp = pbm[patchI];
|
|
||||||
|
|
||||||
forAll(pp, i)
|
forAll(pp, i)
|
||||||
{
|
{
|
||||||
patchFaces[nFaces] = pp.start()+i;
|
patchFaces[nFaces] = pp.start()+i;
|
||||||
patchData[nFaces] = topoDistanceData(patchI, 0);
|
patchData[nFaces] = topoDistanceData(patchi, 0);
|
||||||
nFaces++;
|
++nFaces;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
@ -150,15 +165,17 @@ void subsetFields
|
|||||||
|
|
||||||
const fvMesh& baseMesh = subsetter.baseMesh();
|
const fvMesh& baseMesh = subsetter.baseMesh();
|
||||||
|
|
||||||
if (fieldNames.empty())
|
label nFields = 0;
|
||||||
|
for (const word& fieldName : fieldNames)
|
||||||
{
|
{
|
||||||
return;
|
if (!nFields)
|
||||||
}
|
{
|
||||||
Info<< "Subsetting " << fieldType << " (";
|
Info<< "Subsetting " << fieldType << " (";
|
||||||
forAll(fieldNames, i)
|
}
|
||||||
{
|
else
|
||||||
const word& fieldName = fieldNames[i];
|
{
|
||||||
if (i) Info<< ' ';
|
Info<< ' ';
|
||||||
|
}
|
||||||
Info<< fieldName;
|
Info<< fieldName;
|
||||||
|
|
||||||
FieldType fld
|
FieldType fld
|
||||||
@ -174,12 +191,18 @@ void subsetFields
|
|||||||
baseMesh
|
baseMesh
|
||||||
);
|
);
|
||||||
|
|
||||||
subFields.set(i, subsetter.interpolate(fld));
|
subFields.set(nFields, subsetter.interpolate(fld));
|
||||||
|
|
||||||
// Subsetting adds 'subset' prefix - rename to match original.
|
// Subsetting adds 'subset' prefix - rename to match original.
|
||||||
subFields[i].rename(fieldName);
|
subFields[nFields].rename(fieldName);
|
||||||
|
|
||||||
|
++nFields;
|
||||||
|
}
|
||||||
|
|
||||||
|
if (nFields)
|
||||||
|
{
|
||||||
|
Info<< ')' << nl;
|
||||||
}
|
}
|
||||||
Info<< ")" << nl;
|
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
@ -200,15 +223,17 @@ void subsetPointFields
|
|||||||
|
|
||||||
const fvMesh& baseMesh = subsetter.baseMesh();
|
const fvMesh& baseMesh = subsetter.baseMesh();
|
||||||
|
|
||||||
if (fieldNames.empty())
|
label nFields = 0;
|
||||||
|
for (const word& fieldName : fieldNames)
|
||||||
{
|
{
|
||||||
return;
|
if (!nFields)
|
||||||
}
|
{
|
||||||
Info<< "Subsetting " << fieldType << " (";
|
Info<< "Subsetting " << fieldType << " (";
|
||||||
forAll(fieldNames, i)
|
}
|
||||||
{
|
else
|
||||||
const word& fieldName = fieldNames[i];
|
{
|
||||||
if (i) Info<< ' ';
|
Info<< ' ';
|
||||||
|
}
|
||||||
Info<< fieldName;
|
Info<< fieldName;
|
||||||
|
|
||||||
FieldType fld
|
FieldType fld
|
||||||
@ -224,12 +249,18 @@ void subsetPointFields
|
|||||||
pMesh
|
pMesh
|
||||||
);
|
);
|
||||||
|
|
||||||
subFields.set(i, subsetter.interpolate(fld));
|
subFields.set(nFields, subsetter.interpolate(fld));
|
||||||
|
|
||||||
// Subsetting adds 'subset' prefix - rename to match original.
|
// Subsetting adds 'subset' prefix - rename to match original.
|
||||||
subFields[i].rename(fieldName);
|
subFields[nFields].rename(fieldName);
|
||||||
|
|
||||||
|
++nFields;
|
||||||
|
}
|
||||||
|
|
||||||
|
if (nFields)
|
||||||
|
{
|
||||||
|
Info<< ')' << nl;
|
||||||
}
|
}
|
||||||
Info<< ")" << nl;
|
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
@ -249,15 +280,17 @@ void subsetDimensionedFields
|
|||||||
|
|
||||||
const fvMesh& baseMesh = subsetter.baseMesh();
|
const fvMesh& baseMesh = subsetter.baseMesh();
|
||||||
|
|
||||||
if (fieldNames.empty())
|
label nFields = 0;
|
||||||
|
for (const word& fieldName : fieldNames)
|
||||||
{
|
{
|
||||||
return;
|
if (!nFields)
|
||||||
}
|
{
|
||||||
Info<< "Subsetting " << fieldType << " (";
|
Info<< "Subsetting " << fieldType << " (";
|
||||||
forAll(fieldNames, i)
|
}
|
||||||
{
|
else
|
||||||
const word& fieldName = fieldNames[i];
|
{
|
||||||
if (i) Info<< ' ';
|
Info<< ' ';
|
||||||
|
}
|
||||||
Info<< fieldName;
|
Info<< fieldName;
|
||||||
|
|
||||||
FieldType fld
|
FieldType fld
|
||||||
@ -273,12 +306,18 @@ void subsetDimensionedFields
|
|||||||
baseMesh
|
baseMesh
|
||||||
);
|
);
|
||||||
|
|
||||||
subFields.set(i, subsetter.interpolate(fld));
|
subFields.set(nFields, subsetter.interpolate(fld));
|
||||||
|
|
||||||
// Subsetting adds 'subset' prefix - rename to match original.
|
// Subsetting adds 'subset' prefix - rename to match original.
|
||||||
subFields[i].rename(fieldName);
|
subFields[nFields].rename(fieldName);
|
||||||
|
|
||||||
|
++nFields;
|
||||||
|
}
|
||||||
|
|
||||||
|
if (nFields)
|
||||||
|
{
|
||||||
|
Info<< ')' << nl;
|
||||||
}
|
}
|
||||||
Info<< ")" << nl;
|
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
@ -297,40 +336,33 @@ void subsetTopoSets
|
|||||||
ReadFields<TopoSet>(objects, sets);
|
ReadFields<TopoSet>(objects, sets);
|
||||||
|
|
||||||
subSets.setSize(sets.size());
|
subSets.setSize(sets.size());
|
||||||
forAll(sets, i)
|
forAll(sets, seti)
|
||||||
{
|
{
|
||||||
TopoSet& set = sets[i];
|
const TopoSet& set = sets[seti];
|
||||||
|
|
||||||
Info<< "Subsetting " << set.type() << " " << set.name() << endl;
|
Info<< "Subsetting " << set.type() << " " << set.name() << endl;
|
||||||
|
|
||||||
// Map the data
|
labelHashSet subset(2*min(set.size(), map.size()));
|
||||||
bitSet isSet(set.maxSize(mesh));
|
|
||||||
forAllConstIters(set, iter)
|
for (const label id : map)
|
||||||
{
|
{
|
||||||
isSet.set(iter.key());
|
if (set.found(id))
|
||||||
}
|
|
||||||
label nSet = 0;
|
|
||||||
forAll(map, i)
|
|
||||||
{
|
|
||||||
if (isSet.test(map[i]))
|
|
||||||
{
|
{
|
||||||
++nSet;
|
subset.insert(id);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
subSets.set
|
subSets.set
|
||||||
(
|
(
|
||||||
i,
|
seti,
|
||||||
new TopoSet(subMesh, set.name(), nSet, IOobject::AUTO_WRITE)
|
new TopoSet
|
||||||
|
(
|
||||||
|
subMesh,
|
||||||
|
set.name(),
|
||||||
|
std::move(subset),
|
||||||
|
IOobject::AUTO_WRITE
|
||||||
|
)
|
||||||
);
|
);
|
||||||
TopoSet& subSet = subSets[i];
|
|
||||||
forAll(map, i)
|
|
||||||
{
|
|
||||||
if (isSet.test(map[i]))
|
|
||||||
{
|
|
||||||
subSet.insert(i);
|
|
||||||
}
|
|
||||||
}
|
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
@ -340,25 +372,32 @@ 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 cellSet or on cellZone(s) specified "
|
||||||
|
"as the first command argument."
|
||||||
);
|
);
|
||||||
|
|
||||||
#include "addOverwriteOption.H"
|
#include "addOverwriteOption.H"
|
||||||
#include "addRegionOption.H"
|
#include "addRegionOption.H"
|
||||||
argList::addArgument("cellSet");
|
argList::addArgument("setOrZoneName");
|
||||||
argList::addOption
|
argList::addOption
|
||||||
(
|
(
|
||||||
"patch",
|
"patch",
|
||||||
"name",
|
"name",
|
||||||
"Add exposed internal faces to specified patch instead of to "
|
"Add exposed internal faces to specified patch "
|
||||||
"'oldInternalFaces'"
|
"instead of \"oldInternalFaces\""
|
||||||
);
|
);
|
||||||
argList::addOption
|
argList::addOption
|
||||||
(
|
(
|
||||||
"patches",
|
"patches",
|
||||||
"names",
|
"wordRes",
|
||||||
"Add exposed internal faces to nearest of specified patches"
|
"Add exposed internal faces to closest of specified patches "
|
||||||
" instead of to 'oldInternalFaces'"
|
"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
|
argList::addOption
|
||||||
(
|
(
|
||||||
@ -374,11 +413,13 @@ int main(int argc, char *argv[])
|
|||||||
|
|
||||||
#include "createNamedMesh.H"
|
#include "createNamedMesh.H"
|
||||||
|
|
||||||
const word selectionName = args[1];
|
// arg[1] = word (cellSet) or wordRes (cellZone)
|
||||||
|
// const word selectionName = args[1];
|
||||||
|
|
||||||
word meshInstance = mesh.pointsInstance();
|
word meshInstance = mesh.pointsInstance();
|
||||||
word fieldsInstance = runTime.timeName();
|
word fieldsInstance = runTime.timeName();
|
||||||
|
|
||||||
|
const bool useCellZone = args.found("zone");
|
||||||
const bool overwrite = args.found("overwrite");
|
const bool overwrite = args.found("overwrite");
|
||||||
const bool specifiedInstance = args.readIfPresent
|
const bool specifiedInstance = args.readIfPresent
|
||||||
(
|
(
|
||||||
@ -391,8 +432,6 @@ int main(int argc, char *argv[])
|
|||||||
meshInstance = fieldsInstance;
|
meshInstance = fieldsInstance;
|
||||||
}
|
}
|
||||||
|
|
||||||
Info<< "Reading cell set from " << selectionName << nl << endl;
|
|
||||||
|
|
||||||
|
|
||||||
// Default exposed patch id
|
// Default exposed patch id
|
||||||
labelList exposedPatchIDs(one(), -1);
|
labelList exposedPatchIDs(one(), -1);
|
||||||
@ -401,51 +440,81 @@ int main(int argc, char *argv[])
|
|||||||
{
|
{
|
||||||
const wordRes patchNames(args.getList<wordRe>("patches"));
|
const wordRes patchNames(args.getList<wordRe>("patches"));
|
||||||
|
|
||||||
exposedPatchIDs = mesh.boundaryMesh().patchSet(patchNames).sortedToc();
|
if (patchNames.size() == 1 && !patchNames.first().isPattern())
|
||||||
|
|
||||||
Info<< "Adding exposed internal faces to nearest of patches "
|
|
||||||
<< flatOutput(patchNames) << nl << endl;
|
|
||||||
|
|
||||||
if (exposedPatchIDs.empty())
|
|
||||||
{
|
{
|
||||||
FatalErrorInFunction
|
exposedPatchIDs.first() =
|
||||||
<< nl << "No patches matched. Patches: "
|
getExposedPatchId(mesh, patchNames.first());
|
||||||
<< mesh.boundaryMesh().names()
|
}
|
||||||
<< exit(FatalError);
|
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"))
|
else if (args.found("patch"))
|
||||||
{
|
{
|
||||||
const word patchName = args["patch"];
|
exposedPatchIDs.first() =
|
||||||
|
getExposedPatchId(mesh, args["patch"]);
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
Info<< "Adding exposed internal faces to patch \""
|
||||||
|
<< fvMeshSubset::exposedPatchName
|
||||||
|
<< "\" (created if necessary)" << nl
|
||||||
|
<< nl;
|
||||||
|
}
|
||||||
|
|
||||||
exposedPatchIDs.first() = mesh.boundaryMesh().findPatchID(patchName);
|
|
||||||
|
|
||||||
Info<< "Adding exposed internal faces to patch " << patchName
|
autoPtr<cellSet> cellSetPtr;
|
||||||
<< nl << endl;
|
|
||||||
|
|
||||||
if (exposedPatchIDs.first() == -1)
|
// arg[1] can be a word (for cellSet) or wordRes (for cellZone)
|
||||||
|
|
||||||
|
wordRes zoneNames;
|
||||||
|
if (useCellZone)
|
||||||
|
{
|
||||||
|
List<wordRe> selectionNames = args.getList<wordRe>(1);
|
||||||
|
zoneNames.transfer(selectionNames);
|
||||||
|
|
||||||
|
Info<< "Using cellZone " << flatOutput(zoneNames) << nl << endl;
|
||||||
|
|
||||||
|
if (mesh.cellZones().findIndex(zoneNames) == -1)
|
||||||
{
|
{
|
||||||
FatalErrorInFunction
|
FatalErrorInFunction
|
||||||
<< nl << "No such patch. Patches: "
|
<< "No cellZones found: " << flatOutput(zoneNames) << nl << nl
|
||||||
<< mesh.boundaryMesh().names()
|
|
||||||
<< exit(FatalError);
|
<< exit(FatalError);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
else
|
else
|
||||||
{
|
{
|
||||||
Info<< "Adding exposed internal faces to a patch called"
|
const word selectionName = args[1];
|
||||||
<< " \"oldInternalFaces\" (created if necessary)" << endl
|
|
||||||
<< endl;
|
Info<< "Using cellSet " << selectionName << nl << endl;
|
||||||
|
|
||||||
|
cellSetPtr = autoPtr<cellSet>::New(mesh, selectionName);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
// Mesh subsetting engine
|
// Mesh subsetting engine
|
||||||
fvMeshSubset subsetter(mesh);
|
fvMeshSubset subsetter(mesh);
|
||||||
|
|
||||||
cellSet currentSet(mesh, selectionName);
|
|
||||||
|
|
||||||
{
|
{
|
||||||
bitSet selectedCells = BitSetOps::create(mesh.nCells(), currentSet);
|
bitSet selectedCells =
|
||||||
|
(
|
||||||
|
cellSetPtr
|
||||||
|
? BitSetOps::create(mesh.nCells(), *cellSetPtr)
|
||||||
|
: mesh.cellZones().selection(zoneNames)
|
||||||
|
);
|
||||||
|
|
||||||
if (exposedPatchIDs.size() == 1)
|
if (exposedPatchIDs.size() == 1)
|
||||||
{
|
{
|
||||||
@ -475,6 +544,12 @@ int main(int argc, char *argv[])
|
|||||||
true
|
true
|
||||||
);
|
);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
Info<< "Subset "
|
||||||
|
<< returnReduce(subsetter.subMesh().nCells(), sumOp<label>())
|
||||||
|
<< " of "
|
||||||
|
<< returnReduce(mesh.nCells(), sumOp<label>())
|
||||||
|
<< " cells" << nl << nl;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
@ -485,39 +560,39 @@ int main(int argc, char *argv[])
|
|||||||
// Read vol fields and subset
|
// Read vol fields and subset
|
||||||
// ~~~~~~~~~~~~~~~~~~~~~~~~~~
|
// ~~~~~~~~~~~~~~~~~~~~~~~~~~
|
||||||
|
|
||||||
PtrList<volScalarField> scalarFlds;
|
PtrList<volScalarField> vScalarFlds;
|
||||||
subsetFields(subsetter, availableFields, scalarFlds);
|
subsetFields(subsetter, availableFields, vScalarFlds);
|
||||||
|
|
||||||
PtrList<volVectorField> vectorFlds;
|
PtrList<volVectorField> vVectorFlds;
|
||||||
subsetFields(subsetter, availableFields, vectorFlds);
|
subsetFields(subsetter, availableFields, vVectorFlds);
|
||||||
|
|
||||||
PtrList<volSphericalTensorField> sphTensorFlds;
|
PtrList<volSphericalTensorField> vSphTensorFlds;
|
||||||
subsetFields(subsetter, availableFields, sphTensorFlds);
|
subsetFields(subsetter, availableFields, vSphTensorFlds);
|
||||||
|
|
||||||
PtrList<volSymmTensorField> symmTensorFlds;
|
PtrList<volSymmTensorField> vSymmTensorFlds;
|
||||||
subsetFields(subsetter, availableFields, symmTensorFlds);
|
subsetFields(subsetter, availableFields, vSymmTensorFlds);
|
||||||
|
|
||||||
PtrList<volTensorField> tensorFlds;
|
PtrList<volTensorField> vTensorFlds;
|
||||||
subsetFields(subsetter, availableFields, tensorFlds);
|
subsetFields(subsetter, availableFields, vTensorFlds);
|
||||||
|
|
||||||
|
|
||||||
// Read surface fields and subset
|
// Read surface fields and subset
|
||||||
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
|
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
|
||||||
|
|
||||||
PtrList<surfaceScalarField> surfScalarFlds;
|
PtrList<surfaceScalarField> sScalarFlds;
|
||||||
subsetFields(subsetter, availableFields, surfScalarFlds);
|
subsetFields(subsetter, availableFields, sScalarFlds);
|
||||||
|
|
||||||
PtrList<surfaceVectorField> surfVectorFlds;
|
PtrList<surfaceVectorField> sVectorFlds;
|
||||||
subsetFields(subsetter, availableFields, surfVectorFlds);
|
subsetFields(subsetter, availableFields, sVectorFlds);
|
||||||
|
|
||||||
PtrList<surfaceSphericalTensorField> surfSphTensorFlds;
|
PtrList<surfaceSphericalTensorField> sSphTensorFlds;
|
||||||
subsetFields(subsetter, availableFields, surfSphTensorFlds);
|
subsetFields(subsetter, availableFields, sSphTensorFlds);
|
||||||
|
|
||||||
PtrList<surfaceSymmTensorField> surfSymmTensorFlds;
|
PtrList<surfaceSymmTensorField> sSymmTensorFlds;
|
||||||
subsetFields(subsetter, availableFields, surfSymmTensorFlds);
|
subsetFields(subsetter, availableFields, sSymmTensorFlds);
|
||||||
|
|
||||||
PtrList<surfaceTensorField> surfTensorFlds;
|
PtrList<surfaceTensorField> sTensorFlds;
|
||||||
subsetFields(subsetter, availableFields, surfTensorFlds);
|
subsetFields(subsetter, availableFields, sTensorFlds);
|
||||||
|
|
||||||
|
|
||||||
// Read point fields and subset
|
// Read point fields and subset
|
||||||
@ -525,43 +600,43 @@ int main(int argc, char *argv[])
|
|||||||
|
|
||||||
const pointMesh& pMesh = pointMesh::New(mesh);
|
const pointMesh& pMesh = pointMesh::New(mesh);
|
||||||
|
|
||||||
PtrList<pointScalarField> pointScalarFlds;
|
PtrList<pointScalarField> pScalarFlds;
|
||||||
subsetPointFields(subsetter, pMesh, availableFields, pointScalarFlds);
|
subsetPointFields(subsetter, pMesh, availableFields, pScalarFlds);
|
||||||
|
|
||||||
PtrList<pointVectorField> pointVectorFlds;
|
PtrList<pointVectorField> pVectorFlds;
|
||||||
subsetPointFields(subsetter, pMesh, availableFields, pointVectorFlds);
|
subsetPointFields(subsetter, pMesh, availableFields, pVectorFlds);
|
||||||
|
|
||||||
PtrList<pointSphericalTensorField> pointSphTensorFlds;
|
PtrList<pointSphericalTensorField> pSphTensorFlds;
|
||||||
subsetPointFields(subsetter, pMesh, availableFields, pointSphTensorFlds);
|
subsetPointFields(subsetter, pMesh, availableFields, pSphTensorFlds);
|
||||||
|
|
||||||
PtrList<pointSymmTensorField> pointSymmTensorFlds;
|
PtrList<pointSymmTensorField> pSymmTensorFlds;
|
||||||
subsetPointFields(subsetter, pMesh, availableFields, pointSymmTensorFlds);
|
subsetPointFields(subsetter, pMesh, availableFields, pSymmTensorFlds);
|
||||||
|
|
||||||
PtrList<pointTensorField> pointTensorFlds;
|
PtrList<pointTensorField> pTensorFlds;
|
||||||
subsetPointFields(subsetter, pMesh, availableFields, pointTensorFlds);
|
subsetPointFields(subsetter, pMesh, availableFields, pTensorFlds);
|
||||||
|
|
||||||
|
|
||||||
// Read dimensioned fields and subset
|
// Read dimensioned fields and subset
|
||||||
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
|
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
|
||||||
|
|
||||||
PtrList<volScalarField::Internal> scalarDimFlds;
|
PtrList<volScalarField::Internal> dScalarFlds;
|
||||||
subsetDimensionedFields(subsetter, availableFields, scalarDimFlds);
|
subsetDimensionedFields(subsetter, availableFields, dScalarFlds);
|
||||||
|
|
||||||
PtrList<volVectorField::Internal> vectorDimFlds;
|
PtrList<volVectorField::Internal> dVectorFlds;
|
||||||
subsetDimensionedFields(subsetter, availableFields, vectorDimFlds);
|
subsetDimensionedFields(subsetter, availableFields, dVectorFlds);
|
||||||
|
|
||||||
PtrList<volSphericalTensorField::Internal> sphTensorDimFlds;
|
PtrList<volSphericalTensorField::Internal> dSphTensorFlds;
|
||||||
subsetDimensionedFields(subsetter, availableFields, sphTensorDimFlds);
|
subsetDimensionedFields(subsetter, availableFields, dSphTensorFlds);
|
||||||
|
|
||||||
PtrList<volSymmTensorField::Internal> symmTensorDimFlds;
|
PtrList<volSymmTensorField::Internal> dSymmTensorFlds;
|
||||||
subsetDimensionedFields(subsetter, availableFields, symmTensorDimFlds);
|
subsetDimensionedFields(subsetter, availableFields, dSymmTensorFlds);
|
||||||
|
|
||||||
PtrList<volTensorField::Internal> tensorDimFlds;
|
PtrList<volTensorField::Internal> dTensorFlds;
|
||||||
subsetDimensionedFields(subsetter, availableFields, tensorDimFlds);
|
subsetDimensionedFields(subsetter, availableFields, dTensorFlds);
|
||||||
|
|
||||||
|
|
||||||
// topoSets and subset
|
// Read topoSets and subset
|
||||||
// ~~~~~~~~~~~~~~~~~~~
|
// ~~~~~~~~~~~~~~~~~~~~~~~~
|
||||||
|
|
||||||
PtrList<cellSet> cellSets;
|
PtrList<cellSet> cellSets;
|
||||||
PtrList<faceSet> faceSets;
|
PtrList<faceSet> faceSets;
|
||||||
@ -569,7 +644,10 @@ int main(int argc, char *argv[])
|
|||||||
|
|
||||||
{
|
{
|
||||||
IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
|
IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
|
||||||
objects.remove(currentSet);
|
if (cellSetPtr)
|
||||||
|
{
|
||||||
|
objects.remove(*cellSetPtr);
|
||||||
|
}
|
||||||
subsetTopoSets
|
subsetTopoSets
|
||||||
(
|
(
|
||||||
mesh,
|
mesh,
|
||||||
@ -621,94 +699,34 @@ int main(int argc, char *argv[])
|
|||||||
|
|
||||||
|
|
||||||
// Volume fields
|
// Volume fields
|
||||||
forAll(scalarFlds, i)
|
for (const auto& fld : vScalarFlds) { fld.write(); }
|
||||||
{
|
for (const auto& fld : vVectorFlds) { fld.write(); }
|
||||||
scalarFlds[i].write();
|
for (const auto& fld : vSphTensorFlds) { fld.write(); }
|
||||||
}
|
for (const auto& fld : vSymmTensorFlds) { fld.write(); }
|
||||||
forAll(vectorFlds, i)
|
for (const auto& fld : vTensorFlds) { fld.write(); }
|
||||||
{
|
|
||||||
vectorFlds[i].write();
|
|
||||||
}
|
|
||||||
forAll(sphTensorFlds, i)
|
|
||||||
{
|
|
||||||
sphTensorFlds[i].write();
|
|
||||||
}
|
|
||||||
forAll(symmTensorFlds, i)
|
|
||||||
{
|
|
||||||
symmTensorFlds[i].write();
|
|
||||||
}
|
|
||||||
forAll(tensorFlds, i)
|
|
||||||
{
|
|
||||||
tensorFlds[i].write();
|
|
||||||
}
|
|
||||||
|
|
||||||
// Surface fields.
|
// Surface fields
|
||||||
forAll(surfScalarFlds, i)
|
for (const auto& fld : sScalarFlds) { fld.write(); }
|
||||||
{
|
for (const auto& fld : sVectorFlds) { fld.write(); }
|
||||||
surfScalarFlds[i].write();
|
for (const auto& fld : sSphTensorFlds) { fld.write(); }
|
||||||
}
|
for (const auto& fld : sSymmTensorFlds) { fld.write(); }
|
||||||
forAll(surfVectorFlds, i)
|
for (const auto& fld : sTensorFlds) { fld.write(); }
|
||||||
{
|
|
||||||
surfVectorFlds[i].write();
|
|
||||||
}
|
|
||||||
forAll(surfSphTensorFlds, i)
|
|
||||||
{
|
|
||||||
surfSphTensorFlds[i].write();
|
|
||||||
}
|
|
||||||
forAll(surfSymmTensorFlds, i)
|
|
||||||
{
|
|
||||||
surfSymmTensorFlds[i].write();
|
|
||||||
}
|
|
||||||
forAll(surfTensorFlds, i)
|
|
||||||
{
|
|
||||||
surfTensorFlds[i].write();
|
|
||||||
}
|
|
||||||
|
|
||||||
// Point fields
|
// Point fields
|
||||||
forAll(pointScalarFlds, i)
|
for (const auto& fld : pScalarFlds) { fld.write(); }
|
||||||
{
|
for (const auto& fld : pVectorFlds) { fld.write(); }
|
||||||
pointScalarFlds[i].write();
|
for (const auto& fld : pSphTensorFlds) { fld.write(); }
|
||||||
}
|
for (const auto& fld : pSymmTensorFlds) { fld.write(); }
|
||||||
forAll(pointVectorFlds, i)
|
for (const auto& fld : pTensorFlds) { fld.write(); }
|
||||||
{
|
|
||||||
pointVectorFlds[i].write();
|
|
||||||
}
|
|
||||||
forAll(pointSphTensorFlds, i)
|
|
||||||
{
|
|
||||||
pointSphTensorFlds[i].write();
|
|
||||||
}
|
|
||||||
forAll(pointSymmTensorFlds, i)
|
|
||||||
{
|
|
||||||
pointSymmTensorFlds[i].write();
|
|
||||||
}
|
|
||||||
forAll(pointTensorFlds, i)
|
|
||||||
{
|
|
||||||
pointTensorFlds[i].write();
|
|
||||||
}
|
|
||||||
|
|
||||||
// Dimensioned fields
|
// Dimensioned fields
|
||||||
forAll(scalarDimFlds, i)
|
for (const auto& fld : dScalarFlds) { fld.write(); }
|
||||||
{
|
for (const auto& fld : dVectorFlds) { fld.write(); }
|
||||||
scalarDimFlds[i].write();
|
for (const auto& fld : dSphTensorFlds) { fld.write(); }
|
||||||
}
|
for (const auto& fld : dSymmTensorFlds) { fld.write(); }
|
||||||
forAll(vectorDimFlds, i)
|
for (const auto& fld : dTensorFlds) { fld.write(); }
|
||||||
{
|
|
||||||
vectorDimFlds[i].write();
|
|
||||||
}
|
|
||||||
forAll(sphTensorDimFlds, i)
|
|
||||||
{
|
|
||||||
sphTensorDimFlds[i].write();
|
|
||||||
}
|
|
||||||
forAll(symmTensorDimFlds, i)
|
|
||||||
{
|
|
||||||
symmTensorDimFlds[i].write();
|
|
||||||
}
|
|
||||||
forAll(tensorDimFlds, i)
|
|
||||||
{
|
|
||||||
tensorDimFlds[i].write();
|
|
||||||
}
|
|
||||||
|
|
||||||
Info<< "End\n" << endl;
|
Info<< "\nEnd\n" << endl;
|
||||||
|
|
||||||
return 0;
|
return 0;
|
||||||
}
|
}
|
||||||
|
|||||||
@ -18,7 +18,7 @@ _of_complete_cache_[ansysToFoam]="-case -fileHandler -scale | -noFunctionObjects
|
|||||||
_of_complete_cache_[applyBoundaryLayer]="-case -decomposeParDict -fileHandler -hostRoots -listScalarBCs -listVectorBCs -region -roots -ybl | -listFunctionObjects -listRegisteredSwitches -listSwitches -listTurbulenceModels -listUnsetSwitches -noFunctionObjects -parallel -write-nut -doc -doc-source -help"
|
_of_complete_cache_[applyBoundaryLayer]="-case -decomposeParDict -fileHandler -hostRoots -listScalarBCs -listVectorBCs -region -roots -ybl | -listFunctionObjects -listRegisteredSwitches -listSwitches -listTurbulenceModels -listUnsetSwitches -noFunctionObjects -parallel -write-nut -doc -doc-source -help"
|
||||||
_of_complete_cache_[attachMesh]="-case -fileHandler | -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -overwrite -doc -doc-source -help"
|
_of_complete_cache_[attachMesh]="-case -fileHandler | -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -overwrite -doc -doc-source -help"
|
||||||
_of_complete_cache_[autoPatch]="-case -fileHandler | -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -overwrite -doc -doc-source -help"
|
_of_complete_cache_[autoPatch]="-case -fileHandler | -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -overwrite -doc -doc-source -help"
|
||||||
_of_complete_cache_[blockMesh]="-case -dict -fileHandler -region | -blockTopology -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -noClean -noFunctionObjects -sets -doc -doc-source -help"
|
_of_complete_cache_[blockMesh]="-case -dict -fileHandler -region -time | -blockTopology -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -noClean -noFunctionObjects -sets -doc -doc-source -help"
|
||||||
_of_complete_cache_[boundaryFoam]="-case -fileHandler -listScalarBCs -listVectorBCs | -dry-run -dry-run-write -listFunctionObjects -listFvOptions -listRegisteredSwitches -listSwitches -listTurbulenceModels -listUnsetSwitches -noFunctionObjects -doc -doc-source -help"
|
_of_complete_cache_[boundaryFoam]="-case -fileHandler -listScalarBCs -listVectorBCs | -dry-run -dry-run-write -listFunctionObjects -listFvOptions -listRegisteredSwitches -listSwitches -listTurbulenceModels -listUnsetSwitches -noFunctionObjects -doc -doc-source -help"
|
||||||
_of_complete_cache_[boxTurb]="-case -fileHandler -listScalarBCs -listVectorBCs | -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -noFunctionObjects -doc -doc-source -help"
|
_of_complete_cache_[boxTurb]="-case -fileHandler -listScalarBCs -listVectorBCs | -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -noFunctionObjects -doc -doc-source -help"
|
||||||
_of_complete_cache_[buoyantBoussinesqPimpleFoam]="-case -decomposeParDict -fileHandler -hostRoots -listScalarBCs -listVectorBCs -roots | -dry-run -dry-run-write -listFunctionObjects -listFvOptions -listRegisteredSwitches -listSwitches -listTurbulenceModels -listUnsetSwitches -noFunctionObjects -parallel -postProcess -doc -doc-source -help"
|
_of_complete_cache_[buoyantBoussinesqPimpleFoam]="-case -decomposeParDict -fileHandler -hostRoots -listScalarBCs -listVectorBCs -roots | -dry-run -dry-run-write -listFunctionObjects -listFvOptions -listRegisteredSwitches -listSwitches -listTurbulenceModels -listUnsetSwitches -noFunctionObjects -parallel -postProcess -doc -doc-source -help"
|
||||||
@ -235,7 +235,7 @@ _of_complete_cache_[SRFSimpleFoam]="-case -decomposeParDict -fileHandler -hostRo
|
|||||||
_of_complete_cache_[star4ToFoam]="-case -fileHandler -scale | -ascii -noFunctionObjects -solids -doc -doc-source -help"
|
_of_complete_cache_[star4ToFoam]="-case -fileHandler -scale | -ascii -noFunctionObjects -solids -doc -doc-source -help"
|
||||||
_of_complete_cache_[steadyParticleTracks]="-case -dict -fileHandler -region -time | -constant -latestTime -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -newTimes -noFunctionObjects -noZero -doc -doc-source -help"
|
_of_complete_cache_[steadyParticleTracks]="-case -dict -fileHandler -region -time | -constant -latestTime -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -newTimes -noFunctionObjects -noZero -doc -doc-source -help"
|
||||||
_of_complete_cache_[stitchMesh]="-case -dict -fileHandler -listScalarBCs -listVectorBCs -region -toleranceDict | -integral -intermediate -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -overwrite -partial -perfect -doc -doc-source -help"
|
_of_complete_cache_[stitchMesh]="-case -dict -fileHandler -listScalarBCs -listVectorBCs -region -toleranceDict | -integral -intermediate -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -overwrite -partial -perfect -doc -doc-source -help"
|
||||||
_of_complete_cache_[subsetMesh]="-case -decomposeParDict -fileHandler -hostRoots -listScalarBCs -listVectorBCs -patch -patches -region -resultTime -roots | -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -overwrite -parallel -doc -doc-source -help"
|
_of_complete_cache_[subsetMesh]="-case -decomposeParDict -fileHandler -hostRoots -listScalarBCs -listVectorBCs -patch -patches -region -resultTime -roots | -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -overwrite -parallel -zone -doc -doc-source -help"
|
||||||
_of_complete_cache_[subsetToPatch]="-case -fileHandler | -noFunctionObjects -doc -doc-source -help"
|
_of_complete_cache_[subsetToPatch]="-case -fileHandler | -noFunctionObjects -doc -doc-source -help"
|
||||||
_of_complete_cache_[surfaceAdd]="-case -fileHandler -points -scale | -mergeRegions -noFunctionObjects -doc -doc-source -help"
|
_of_complete_cache_[surfaceAdd]="-case -fileHandler -points -scale | -mergeRegions -noFunctionObjects -doc -doc-source -help"
|
||||||
_of_complete_cache_[surfaceBooleanFeatures]="-case -fileHandler -scale -trim | -invertedSpace -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -noFunctionObjects -perturb -surf1Baffle -surf2Baffle -doc -doc-source -help"
|
_of_complete_cache_[surfaceBooleanFeatures]="-case -fileHandler -scale -trim | -invertedSpace -listFunctionObjects -listRegisteredSwitches -listSwitches -listUnsetSwitches -noFunctionObjects -perturb -surf1Baffle -surf2Baffle -doc -doc-source -help"
|
||||||
|
|||||||
@ -56,6 +56,11 @@ inline void markUsed
|
|||||||
} // End anonymous namespace
|
} // End anonymous namespace
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::word Foam::fvMeshSubset::exposedPatchName("oldInternalFaces");
|
||||||
|
|
||||||
|
|
||||||
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
||||||
|
|
||||||
bool Foam::fvMeshSubset::checkCellSubset() const
|
bool Foam::fvMeshSubset::checkCellSubset() const
|
||||||
@ -571,7 +576,7 @@ void Foam::fvMeshSubset::setCellSubset
|
|||||||
{
|
{
|
||||||
// No explicit patch specified. Put in oldInternalFaces patch.
|
// No explicit patch specified. Put in oldInternalFaces patch.
|
||||||
// Check if patch with this name already exists.
|
// Check if patch with this name already exists.
|
||||||
wantedPatchID = oldPatches.findPatchID("oldInternalFaces");
|
wantedPatchID = oldPatches.findPatchID(exposedPatchName);
|
||||||
}
|
}
|
||||||
else if (wantedPatchID < 0 || wantedPatchID >= oldPatches.size())
|
else if (wantedPatchID < 0 || wantedPatchID >= oldPatches.size())
|
||||||
{
|
{
|
||||||
@ -1078,7 +1083,7 @@ void Foam::fvMeshSubset::setCellSubset
|
|||||||
{
|
{
|
||||||
newBoundary[nNewPatches] = new emptyPolyPatch
|
newBoundary[nNewPatches] = new emptyPolyPatch
|
||||||
(
|
(
|
||||||
"oldInternalFaces",
|
exposedPatchName,
|
||||||
boundaryPatchSizes[oldInternalPatchID],
|
boundaryPatchSizes[oldInternalPatchID],
|
||||||
patchStart,
|
patchStart,
|
||||||
nNewPatches,
|
nNewPatches,
|
||||||
@ -1086,7 +1091,7 @@ void Foam::fvMeshSubset::setCellSubset
|
|||||||
emptyPolyPatch::typeName
|
emptyPolyPatch::typeName
|
||||||
);
|
);
|
||||||
|
|
||||||
//Pout<< " oldInternalFaces : "
|
//Pout<< " " << exposedPatchName << " : "
|
||||||
// << boundaryPatchSizes[oldInternalPatchID] << endl;
|
// << boundaryPatchSizes[oldInternalPatchID] << endl;
|
||||||
|
|
||||||
// The index for the first patch is -1 as it originates from
|
// The index for the first patch is -1 as it originates from
|
||||||
|
|||||||
@ -153,6 +153,12 @@ class fvMeshSubset
|
|||||||
|
|
||||||
public:
|
public:
|
||||||
|
|
||||||
|
// Static Data Members
|
||||||
|
|
||||||
|
//- Name for exposed internal faces (default: oldInternalFaces)
|
||||||
|
static word exposedPatchName;
|
||||||
|
|
||||||
|
|
||||||
// Constructors
|
// Constructors
|
||||||
|
|
||||||
//- Construct given a mesh to subset
|
//- Construct given a mesh to subset
|
||||||
|
|||||||
Reference in New Issue
Block a user