/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2016 OpenFOAM Foundation
-------------------------------------------------------------------------------
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
Test-GAMGAgglomeration
Description
Test application for GAMG agglomeration. Hardcoded to expect GAMG on p.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "GAMGAgglomeration.H"
#include "OFstream.H"
#include "meshTools.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// Main program:
int main(int argc, char *argv[])
{
argList::addBoolOption
(
"writeObj",
"write obj files of agglomeration"
);
argList::addBoolOption
(
"normalise",
"normalise agglomeration (0..1)"
);
#include "setRootCase.H"
#include "createTime.H"
bool writeObj = args.found("writeObj");
bool normalise = args.found("normalise");
#include "createMesh.H"
const fvSolution& sol = static_cast(mesh);
const dictionary& pDict = sol.subDict("solvers").subDict("p");
const GAMGAgglomeration& agglom = GAMGAgglomeration::New
(
mesh,
pDict
);
labelList cellToCoarse(identity(mesh.nCells()));
labelListList coarseToCell(invertOneToMany(mesh.nCells(), cellToCoarse));
++runTime;
// Write initial agglomeration
{
volScalarField scalarAgglomeration
(
IOobject
(
"agglomeration",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
mesh,
dimensionedScalar(dimless, Zero)
);
scalarField& fld = scalarAgglomeration.primitiveFieldRef();
forAll(fld, celli)
{
fld[celli] = cellToCoarse[celli];
}
fld /= max(fld);
scalarAgglomeration.correctBoundaryConditions();
scalarAgglomeration.write();
Info<< "Writing initial cell distribution to "
<< runTime.timeName() << endl;
}
for (label level = 0; level < agglom.size(); level++)
{
++runTime;
Info<< "Time = " << runTime.timeName() << nl << endl;
const labelList& addr = agglom.restrictAddressing(level);
label coarseSize = max(addr)+1;
Info<< "Level : " << level << endl
<< returnReduce(addr.size(), sumOp