mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
When the GeometricBoundaryField template class was originally written it
was a separate class in the Foam namespace rather than a sub-class of
GeometricField as it is now. Without loss of clarity and simplifying
code which access the boundary field of GeometricFields it is better
that GeometricBoundaryField be renamed Boundary for consistency with the
new naming convention for the type of the dimensioned internal field:
Internal, see commit 4a57b9be2e
This is a very simple text substitution change which can be applied to
any code which compiles with the OpenFOAM-dev libraries.
208 lines
5.9 KiB
C
208 lines
5.9 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2011-2016 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/>.
|
|
|
|
Application
|
|
faceAgglomerate
|
|
|
|
Description
|
|
|
|
Agglomerate boundary faces using the pairPatchAgglomeration algorithm.
|
|
It writes a map from the fine to coarse grid.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "argList.H"
|
|
#include "fvMesh.H"
|
|
#include "Time.H"
|
|
#include "volFields.H"
|
|
#include "CompactListList.H"
|
|
#include "unitConversion.H"
|
|
#include "pairPatchAgglomeration.H"
|
|
#include "labelListIOList.H"
|
|
#include "syncTools.H"
|
|
#include "globalIndex.H"
|
|
|
|
using namespace Foam;
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
int main(int argc, char *argv[])
|
|
{
|
|
#include "addRegionOption.H"
|
|
#include "addDictOption.H"
|
|
#include "setRootCase.H"
|
|
#include "createTime.H"
|
|
#include "createNamedMesh.H"
|
|
|
|
const word dictName("viewFactorsDict");
|
|
|
|
#include "setConstantMeshDictionaryIO.H"
|
|
|
|
// Read control dictionary
|
|
const IOdictionary agglomDict(dictIO);
|
|
|
|
bool writeAgglom = readBool(agglomDict.lookup("writeFacesAgglomeration"));
|
|
|
|
|
|
const polyBoundaryMesh& boundary = mesh.boundaryMesh();
|
|
|
|
labelListIOList finalAgglom
|
|
(
|
|
IOobject
|
|
(
|
|
"finalAgglom",
|
|
mesh.facesInstance(),
|
|
mesh,
|
|
IOobject::NO_READ,
|
|
IOobject::NO_WRITE,
|
|
false
|
|
),
|
|
boundary.size()
|
|
);
|
|
|
|
label nCoarseFaces = 0;
|
|
|
|
forAllConstIter(dictionary, agglomDict, iter)
|
|
{
|
|
labelList patchids = boundary.findIndices(iter().keyword());
|
|
forAll(patchids, i)
|
|
{
|
|
label patchi = patchids[i];
|
|
const polyPatch& pp = boundary[patchi];
|
|
if (!pp.coupled())
|
|
{
|
|
Info << "\nAgglomerating patch : " << pp.name() << endl;
|
|
pairPatchAgglomeration agglomObject
|
|
(
|
|
pp,
|
|
agglomDict.subDict(pp.name())
|
|
);
|
|
agglomObject.agglomerate();
|
|
finalAgglom[patchi] =
|
|
agglomObject.restrictTopBottomAddressing();
|
|
|
|
if (finalAgglom[patchi].size())
|
|
{
|
|
nCoarseFaces += max(finalAgglom[patchi] + 1);
|
|
}
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
// - All patches which are not agglomarated are identity for finalAgglom
|
|
forAll(boundary, patchid)
|
|
{
|
|
if (finalAgglom[patchid].size() == 0)
|
|
{
|
|
finalAgglom[patchid] = identity(boundary[patchid].size());
|
|
}
|
|
}
|
|
|
|
// Sync agglomeration across coupled patches
|
|
labelList nbrAgglom(mesh.nFaces() - mesh.nInternalFaces(), -1);
|
|
|
|
forAll(boundary, patchid)
|
|
{
|
|
const polyPatch& pp = boundary[patchid];
|
|
if (pp.coupled())
|
|
{
|
|
finalAgglom[patchid] = identity(pp.size());
|
|
forAll(pp, i)
|
|
{
|
|
nbrAgglom[pp.start() - mesh.nInternalFaces() + i] =
|
|
finalAgglom[patchid][i];
|
|
}
|
|
}
|
|
}
|
|
|
|
syncTools::swapBoundaryFaceList(mesh, nbrAgglom);
|
|
forAll(boundary, patchid)
|
|
{
|
|
const polyPatch& pp = boundary[patchid];
|
|
if (pp.coupled() && !refCast<const coupledPolyPatch>(pp).owner())
|
|
{
|
|
forAll(pp, i)
|
|
{
|
|
finalAgglom[patchid][i] =
|
|
nbrAgglom[pp.start() - mesh.nInternalFaces() + i];
|
|
}
|
|
}
|
|
}
|
|
|
|
finalAgglom.write();
|
|
|
|
if (writeAgglom)
|
|
{
|
|
globalIndex index(nCoarseFaces);
|
|
volScalarField facesAgglomeration
|
|
(
|
|
IOobject
|
|
(
|
|
"facesAgglomeration",
|
|
mesh.time().timeName(),
|
|
mesh,
|
|
IOobject::NO_READ,
|
|
IOobject::NO_WRITE
|
|
),
|
|
mesh,
|
|
dimensionedScalar("facesAgglomeration", dimless, 0)
|
|
);
|
|
|
|
volScalarField::Boundary& facesAgglomerationBf =
|
|
facesAgglomeration.boundaryFieldRef();
|
|
|
|
label coarsePatchIndex = 0;
|
|
forAll(boundary, patchid)
|
|
{
|
|
const polyPatch& pp = boundary[patchid];
|
|
if (pp.size() > 0)
|
|
{
|
|
fvPatchScalarField& bFacesAgglomeration =
|
|
facesAgglomerationBf[patchid];
|
|
|
|
forAll(bFacesAgglomeration, j)
|
|
{
|
|
bFacesAgglomeration[j] =
|
|
index.toGlobal
|
|
(
|
|
Pstream::myProcNo(),
|
|
finalAgglom[patchid][j] + coarsePatchIndex
|
|
);
|
|
}
|
|
|
|
coarsePatchIndex += max(finalAgglom[patchid]) + 1;
|
|
}
|
|
}
|
|
|
|
Info<< "\nWriting facesAgglomeration" << endl;
|
|
facesAgglomeration.write();
|
|
}
|
|
|
|
Info<< "End\n" << endl;
|
|
return 0;
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|