mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
193 lines
5.1 KiB
C
193 lines
5.1 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2011-2011 OpenCFD Ltd.
|
|
\\/ 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"
|
|
|
|
using namespace Foam;
|
|
|
|
// Main program:
|
|
|
|
int main(int argc, char *argv[])
|
|
{
|
|
|
|
#include "addRegionOption.H"
|
|
#include "setRootCase.H"
|
|
#include "createTime.H"
|
|
#include "createNamedMesh.H"
|
|
|
|
const polyBoundaryMesh& patches = mesh.boundaryMesh();
|
|
|
|
labelListIOList finalAgglom
|
|
(
|
|
IOobject
|
|
(
|
|
"finalAgglom",
|
|
mesh.facesInstance(),
|
|
mesh,
|
|
IOobject::NO_READ,
|
|
IOobject::NO_WRITE,
|
|
false
|
|
),
|
|
patches.size()
|
|
);
|
|
|
|
|
|
// Read view factor dictionary
|
|
IOdictionary viewFactorDict
|
|
(
|
|
IOobject
|
|
(
|
|
"viewFactorsDict",
|
|
runTime.constant(),
|
|
mesh,
|
|
IOobject::MUST_READ,
|
|
IOobject::NO_WRITE
|
|
)
|
|
);
|
|
|
|
bool writeAgglo =
|
|
readBool(viewFactorDict.lookup("writeFacesAgglomeration"));
|
|
|
|
const polyBoundaryMesh& boundary = mesh.boundaryMesh();
|
|
|
|
forAll(boundary, patchId)
|
|
{
|
|
const polyPatch& pp = boundary[patchId];
|
|
|
|
label patchI = pp.index();
|
|
finalAgglom[patchI].setSize(pp.size(), 0);
|
|
|
|
if (pp.size() > 0 && !pp.coupled())
|
|
{
|
|
if (viewFactorDict.found(pp.name()))
|
|
{
|
|
Info << "\nAgglomerating name : " << pp.name() << endl;
|
|
pairPatchAgglomeration agglomObject
|
|
(
|
|
pp,
|
|
viewFactorDict.subDict(pp.name())
|
|
);
|
|
agglomObject.agglomerate();
|
|
finalAgglom[patchI] =
|
|
agglomObject.restrictTopBottomAddressing();
|
|
}
|
|
else
|
|
{
|
|
FatalErrorIn
|
|
(
|
|
"main(int argc, char *argv[])"
|
|
) << pp.name()
|
|
<< " not found in dictionary : "
|
|
<< viewFactorDict.name()
|
|
<< exit(FatalError);
|
|
}
|
|
}
|
|
}
|
|
|
|
// 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 (writeAgglo)
|
|
{
|
|
volScalarField facesAgglomeration
|
|
(
|
|
IOobject
|
|
(
|
|
"facesAgglomeration",
|
|
mesh.time().timeName(),
|
|
mesh,
|
|
IOobject::NO_READ,
|
|
IOobject::NO_WRITE
|
|
),
|
|
mesh,
|
|
dimensionedScalar("facesAgglomeration", dimless, 0)
|
|
);
|
|
|
|
forAll(boundary, patchId)
|
|
{
|
|
|
|
fvPatchScalarField& bFacesAgglomeration =
|
|
facesAgglomeration.boundaryField()[patchId];
|
|
|
|
forAll(bFacesAgglomeration, j)
|
|
{
|
|
bFacesAgglomeration[j] = finalAgglom[patchId][j];
|
|
}
|
|
}
|
|
|
|
Info << "\nWriting facesAgglomeration..." << endl;
|
|
facesAgglomeration.write();
|
|
}
|
|
|
|
Info<< "End\n" << endl;
|
|
return 0;
|
|
}
|