/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2016 OpenFOAM Foundation
Copyright (C) 2016-2022 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
createBaffles
Group
grpMeshManipulationUtilities
Description
Makes internal faces into boundary faces. Does not duplicate points, unlike
mergeOrSplitBaffles.
Note: if any coupled patch face is selected for baffling the opposite
member has to be selected for baffling as well.
- if the patch already exists will not override it nor its fields
- if the patch does not exist it will be created together with 'calculated'
patchfields unless the field is mentioned in the patchFields section.
- any 0-sized patches (since faces have been moved out) will get removed
\*---------------------------------------------------------------------------*/
#include "argList.H"
#include "Time.H"
#include "polyTopoChange.H"
#include "polyModifyFace.H"
#include "polyAddFace.H"
#include "ReadFields.H"
#include "volFields.H"
#include "surfaceFields.H"
#include "fvMeshMapper.H"
#include "faceSelection.H"
#include "fvMeshTools.H"
#include "topoSet.H"
#include "processorPolyPatch.H"
#include "processorMeshes.H"
using namespace Foam;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
label addPatch
(
fvMesh& mesh,
const word& patchName,
const word& groupName,
const dictionary& patchDict
)
{
const polyBoundaryMesh& pbm = mesh.boundaryMesh();
if (pbm.findPatchID(patchName) == -1)
{
autoPtr ppPtr
(
polyPatch::New
(
patchName,
patchDict,
0,
pbm
)
);
auto& pp = *ppPtr;
if (!groupName.empty())
{
pp.inGroups().appendUniq(groupName);
}
// Add patch, create calculated everywhere
fvMeshTools::addPatch
(
mesh,
pp,
dictionary(), // do not set specialised patchFields
calculatedFvPatchField::typeName,
true // parallel sync'ed addition
);
}
else
{
Info<< "Patch '" << patchName
<< "' already exists. Only "
<< "moving patch faces - type will remain the same"
<< endl;
}
return pbm.findPatchID(patchName);
}
// Filter out the empty patches.
void filterPatches(fvMesh& mesh, const wordHashSet& addedPatchNames)
{
// Remove any zero-sized ones. Assumes
// - processor patches are already only there if needed
// - all other patches are available on all processors
// - but coupled ones might still be needed, even if zero-size
// (e.g. processorCyclic)
// See also logic in createPatch.
const polyBoundaryMesh& pbm = mesh.boundaryMesh();
labelList oldToNew(pbm.size(), -1);
label newPatchi = 0;
forAll(pbm, patchi)
{
const polyPatch& pp = pbm[patchi];
if (!isA(pp))
{
if
(
isA(pp)
|| returnReduceOr(pp.size())
|| addedPatchNames.found(pp.name())
)
{
// Coupled (and unknown size) or uncoupled and used
oldToNew[patchi] = newPatchi++;
}
}
}
forAll(pbm, patchi)
{
const polyPatch& pp = pbm[patchi];
if (isA(pp))
{
oldToNew[patchi] = newPatchi++;
}
}
const label nKeepPatches = newPatchi;
// Shuffle unused ones to end
if (nKeepPatches != pbm.size())
{
Info<< endl
<< "Removing zero-sized patches:" << endl << incrIndent;
forAll(oldToNew, patchi)
{
if (oldToNew[patchi] == -1)
{
Info<< indent << pbm[patchi].name()
<< " type " << pbm[patchi].type()
<< " at position " << patchi << endl;
oldToNew[patchi] = newPatchi++;
}
}
Info<< decrIndent;
fvMeshTools::reorderPatches(mesh, oldToNew, nKeepPatches, true);
Info<< endl;
}
}
void modifyOrAddFace
(
polyTopoChange& meshMod,
const face& f,
const label facei,
const label own,
const bool flipFaceFlux,
const label newPatchi,
const label zoneID,
const bool zoneFlip,
bitSet& modifiedFace
)
{
if (modifiedFace.set(facei))
{
// First usage of face. Modify.
meshMod.setAction
(
polyModifyFace
(
f, // modified face
facei, // label of face
own, // owner
-1, // neighbour
flipFaceFlux, // face flip
newPatchi, // patch for face
false, // remove from zone
zoneID, // zone for face
zoneFlip // face flip in zone
)
);
}
else
{
// Second or more usage of face. Add.
meshMod.setAction
(
polyAddFace
(
f, // modified face
own, // owner
-1, // neighbour
-1, // master point
-1, // master edge
facei, // master face
flipFaceFlux, // face flip
newPatchi, // patch for face
zoneID, // zone for face
zoneFlip // face flip in zone
)
);
}
}
// Create faces for fZone faces. Usually newMasterPatches, newSlavePatches
// only size one but can be more for duplicate baffle sets
void createFaces
(
const bool internalFacesOnly,
const fvMesh& mesh,
const faceZone& fZone,
const labelList& newMasterPatches,
const labelList& newSlavePatches,
polyTopoChange& meshMod,
bitSet& modifiedFace,
label& nModified
)
{
const polyBoundaryMesh& pbm = mesh.boundaryMesh();
forAll(newMasterPatches, i)
{
// Pass 1. Do selected side of zone
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
{
label zoneFacei = fZone.whichFace(facei);
if (zoneFacei != -1)
{
if (!fZone.flipMap()[zoneFacei])
{
// Use owner side of face
modifyOrAddFace
(
meshMod,
mesh.faces()[facei], // modified face
facei, // label of face
mesh.faceOwner()[facei],// owner
false, // face flip
newMasterPatches[i], // patch for face
fZone.index(), // zone for face
false, // face flip in zone
modifiedFace // modify or add status
);
}
else
{
// Use neighbour side of face.
// To keep faceZone pointing out of original neighbour
// we don't need to set faceFlip since that cell
// now becomes the owner
modifyOrAddFace
(
meshMod,
mesh.faces()[facei].reverseFace(), // modified face
facei, // label of face
mesh.faceNeighbour()[facei],// owner
true, // face flip
newMasterPatches[i], // patch for face
fZone.index(), // zone for face
false, // face flip in zone
modifiedFace // modify or add status
);
}
nModified++;
}
}
// Pass 2. Do other side of zone
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
{
label zoneFacei = fZone.whichFace(facei);
if (zoneFacei != -1)
{
if (!fZone.flipMap()[zoneFacei])
{
// Use neighbour side of face
modifyOrAddFace
(
meshMod,
mesh.faces()[facei].reverseFace(), // modified face
facei, // label of face
mesh.faceNeighbour()[facei], // owner
true, // face flip
newSlavePatches[i], // patch for face
fZone.index(), // zone for face
true, // face flip in zone
modifiedFace // modify or add
);
}
else
{
// Use owner side of face
modifyOrAddFace
(
meshMod,
mesh.faces()[facei], // modified face
facei, // label of face
mesh.faceOwner()[facei],// owner
false, // face flip
newSlavePatches[i], // patch for face
fZone.index(), // zone for face
true, // face flip in zone
modifiedFace // modify or add status
);
}
}
}
// Modify any boundary faces
// ~~~~~~~~~~~~~~~~~~~~~~~~~
// Normal boundary:
// - move to new patch. Might already be back-to-back baffle
// you want to add cyclic to. Do warn though.
//
// Processor boundary:
// - do not move to cyclic
// - add normal patches though.
// For warning once per patch.
labelHashSet patchWarned;
forAll(pbm, patchi)
{
const polyPatch& pp = pbm[patchi];
const label newMasterPatchi = newMasterPatches[i];
const label newSlavePatchi = newSlavePatches[i];
if
(
pp.coupled()
&& (
pbm[newMasterPatchi].coupled()
|| pbm[newSlavePatchi].coupled()
)
)
{
// Do not allow coupled faces to be moved to different
// coupled patches.
}
else if (pp.coupled() || !internalFacesOnly)
{
forAll(pp, i)
{
label facei = pp.start()+i;
label zoneFacei = fZone.whichFace(facei);
if (zoneFacei != -1)
{
if (patchWarned.insert(patchi))
{
WarningInFunction
<< "Found boundary face (in patch "
<< pp.name()
<< ") in faceZone " << fZone.name()
<< " to convert to baffle patches "
<< pbm[newMasterPatchi].name() << "/"
<< pbm[newSlavePatchi].name()
<< endl
<< " Set internalFacesOnly to true in the"
<< " createBaffles control dictionary if you"
<< " don't wish to convert boundary faces."
<< endl;
}
modifyOrAddFace
(
meshMod,
mesh.faces()[facei], // modified face
facei, // label of face
mesh.faceOwner()[facei], // owner
false, // face flip
fZone.flipMap()[zoneFacei]
? newSlavePatchi
: newMasterPatchi, // patch for face
fZone.index(), // zone for face
fZone.flipMap()[zoneFacei], // face flip in zone
modifiedFace // modify or add
);
nModified++;
}
}
}
}
}
}
int main(int argc, char *argv[])
{
argList::addNote
(
"Makes internal faces into boundary faces.\n"
"Does not duplicate points."
);
argList::addOption("dict", "file", "Alternative createBafflesDict");
#include "addOverwriteOption.H"
#include "addRegionOption.H"
argList::noFunctionObjects(); // Never use function objects
#include "setRootCase.H"
#include "createTime.H"
#include "createNamedMesh.H"
const bool overwrite = args.found("overwrite");
const word oldInstance = mesh.pointsInstance();
const word dictName("createBafflesDict");
#include "setSystemMeshDictionaryIO.H"
bool internalFacesOnly(false);
bool noFields(false);
PtrList selectors;
{
Info<< "Reading baffle criteria from " << dictIO.name() << nl << endl;
IOdictionary dict(dictIO);
internalFacesOnly = dict.get("internalFacesOnly");
noFields = dict.getOrDefault("noFields", false);
const dictionary& selectionsDict = dict.subDict("baffles");
selectors.resize(selectionsDict.size());
label nselect = 0;
for (const entry& dEntry : selectionsDict)
{
if (dEntry.isDict())
{
selectors.set
(
nselect,
faceSelection::New(dEntry.keyword(), mesh, dEntry.dict())
);
++nselect;
}
}
selectors.resize(nselect);
}
if (internalFacesOnly)
{
Info<< "Not converting faces on non-coupled patches." << nl << endl;
}
// Read objects in time directory
IOobjectList objects(mesh, runTime.timeName());
// Read vol fields.
Info<< "Reading geometric fields" << nl << endl;
PtrList vsFlds;
if (!noFields) ReadFields(mesh, objects, vsFlds);
PtrList vvFlds;
if (!noFields) ReadFields(mesh, objects, vvFlds);
PtrList vstFlds;
if (!noFields) ReadFields(mesh, objects, vstFlds);
PtrList vsymtFlds;
if (!noFields) ReadFields(mesh, objects, vsymtFlds);
PtrList vtFlds;
if (!noFields) ReadFields(mesh, objects, vtFlds);
// Read surface fields.
PtrList ssFlds;
if (!noFields) ReadFields(mesh, objects, ssFlds);
PtrList svFlds;
if (!noFields) ReadFields(mesh, objects, svFlds);
PtrList sstFlds;
if (!noFields) ReadFields(mesh, objects, sstFlds);
PtrList ssymtFlds;
if (!noFields) ReadFields(mesh, objects, ssymtFlds);
PtrList stFlds;
if (!noFields) ReadFields(mesh, objects, stFlds);
// Creating (if necessary) faceZones
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
forAll(selectors, selectorI)
{
const word& name = selectors[selectorI].name();
if (mesh.faceZones().findZoneID(name) == -1)
{
mesh.faceZones().clearAddressing();
const label zoneID = mesh.faceZones().size();
mesh.faceZones().setSize(zoneID+1);
mesh.faceZones().set
(
zoneID,
new faceZone(name, labelList(), false, zoneID, mesh.faceZones())
);
}
}
// Select faces
// ~~~~~~~~~~~~
//- Per face zoneID it is in and flip status.
labelList faceToZoneID(mesh.nFaces(), -1);
boolList faceToFlip(mesh.nFaces(), false);
forAll(selectors, selectorI)
{
const word& name = selectors[selectorI].name();
label zoneID = mesh.faceZones().findZoneID(name);
selectors[selectorI].select(zoneID, faceToZoneID, faceToFlip);
}
// Add faces to faceZones
labelList nFaces(mesh.faceZones().size(), Zero);
forAll(faceToZoneID, facei)
{
label zoneID = faceToZoneID[facei];
if (zoneID != -1)
{
nFaces[zoneID]++;
}
}
forAll(selectors, selectorI)
{
const word& name = selectors[selectorI].name();
label zoneID = mesh.faceZones().findZoneID(name);
label& n = nFaces[zoneID];
labelList addr(n);
boolList flip(n);
n = 0;
forAll(faceToZoneID, facei)
{
label zone = faceToZoneID[facei];
if (zone == zoneID)
{
addr[n] = facei;
flip[n] = faceToFlip[facei];
++n;
}
}
Info<< "Created zone " << name
<< " at index " << zoneID
<< " with " << returnReduce(n, sumOp