/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2012 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 .
Application
patchSummary
Description
Writes fields and boundary condition info for each patch at each requested
time instance.
Default action is to write a single entry for patches/patchGroups with the
same boundary conditions. Use the -expand option to print every patch
separately. In case of multiple groups matching it will print only the
first one.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "volFields.H"
#include "IOobjectList.H"
#include "patchSummaryTemplates.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
timeSelector::addOptions();
# include "addRegionOption.H"
argList::addBoolOption
(
"expand",
"Do not combine patches"
);
# include "setRootCase.H"
# include "createTime.H"
instantList timeDirs = timeSelector::select0(runTime, args);
const bool expand = args.optionFound("expand");
# include "createNamedMesh.H"
const polyBoundaryMesh& bm = mesh.boundaryMesh();
forAll(timeDirs, timeI)
{
runTime.setTime(timeDirs[timeI], timeI);
Info<< "Time = " << runTime.timeName() << nl << endl;
// Update the mesh if changed
if (mesh.readUpdate() == polyMesh::TOPO_PATCH_CHANGE)
{
Info<< "Detected changed patches. Recreating patch group table."
<< endl;
}
const IOobjectList fieldObjs(mesh, runTime.timeName());
const wordList objNames = fieldObjs.names();
PtrList vsf(objNames.size());
PtrList vvf(objNames.size());
PtrList vsptf(objNames.size());
PtrList vsytf(objNames.size());
PtrList vtf(objNames.size());
Info<< "Valid fields:" << endl;
forAll(objNames, objI)
{
IOobject obj
(
objNames[objI],
runTime.timeName(),
mesh,
IOobject::MUST_READ
);
if (obj.headerOk())
{
addToFieldList(vsf, obj, objI, mesh);
addToFieldList(vvf, obj, objI, mesh);
addToFieldList(vsptf, obj, objI, mesh);
addToFieldList(vsytf, obj, objI, mesh);
addToFieldList(vtf, obj, objI, mesh);
}
}
Info<< endl;
if (expand)
{
// Print each patch separately
forAll(bm, patchI)
{
Info<< bm[patchI].type() << "\t: " << bm[patchI].name() << nl;
outputFieldList(vsf, patchI);
outputFieldList(vvf, patchI);
outputFieldList(vsptf, patchI);
outputFieldList(vsytf, patchI);
outputFieldList(vtf, patchI);
Info<< endl;
}
}
else
{
// Collect for each patch the bc type per field. Merge similar
// patches.
// Per 'group', the map from fieldname to patchfield type
DynamicList > fieldToTypes(bm.size());
// Per 'group' the patches
DynamicList > groupToPatches(bm.size());
forAll(bm, patchI)
{
HashTable fieldToType;
collectFieldList(vsf, patchI, fieldToType);
collectFieldList(vvf, patchI, fieldToType);
collectFieldList(vsptf, patchI, fieldToType);
collectFieldList(vsytf, patchI, fieldToType);
collectFieldList(vtf, patchI, fieldToType);
label groupI = findIndex(fieldToTypes, fieldToType);
if (groupI == -1)
{
DynamicList