mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
331 lines
8.5 KiB
C
331 lines
8.5 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2009-2009 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 2 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, write to the Free Software Foundation,
|
|
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "faceZonesIntegration.H"
|
|
#include "volFields.H"
|
|
#include "dictionary.H"
|
|
#include "Time.H"
|
|
#include "IOmanip.H"
|
|
#include "ListListOps.H"
|
|
#include "processorPolyPatch.H"
|
|
#include "cyclicPolyPatch.H"
|
|
|
|
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
|
|
|
namespace Foam
|
|
{
|
|
defineTypeNameAndDebug(faceZonesIntegration, 0);
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
|
|
|
Foam::faceZonesIntegration::faceZonesIntegration
|
|
(
|
|
const word& name,
|
|
const objectRegistry& obr,
|
|
const dictionary& dict,
|
|
const bool loadFromFiles
|
|
)
|
|
:
|
|
name_(name),
|
|
obr_(obr),
|
|
active_(true),
|
|
log_(false),
|
|
zoneNames_(),
|
|
fieldNames_(),
|
|
filePtr_(NULL)
|
|
{
|
|
// Check if the available mesh is an fvMesh otherise deactivate
|
|
if (!isA<fvMesh>(obr_))
|
|
{
|
|
active_ = false;
|
|
WarningIn
|
|
(
|
|
"Foam::faceZonesIntegration::faceZonesIntegration"
|
|
"("
|
|
"const word&, "
|
|
"const objectRegistry&, "
|
|
"const dictionary&, "
|
|
"const bool"
|
|
")"
|
|
) << "No fvMesh available, deactivating."
|
|
<< endl;
|
|
}
|
|
|
|
read(dict);
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
|
|
|
Foam::faceZonesIntegration::~faceZonesIntegration()
|
|
{}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
|
|
|
void Foam::faceZonesIntegration::read(const dictionary& dict)
|
|
{
|
|
if (active_)
|
|
{
|
|
log_ = dict.lookupOrDefault<Switch>("log", false);
|
|
|
|
dict.lookup("fields") >> fieldNames_;
|
|
|
|
dict.lookup("faceZones") >> zoneNames_;
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::faceZonesIntegration::makeFile()
|
|
{
|
|
// Create the face Zone file if not already created
|
|
if (filePtr_.empty())
|
|
{
|
|
if (debug)
|
|
{
|
|
Info<< "Creating faceZonesIntegration file." << endl;
|
|
}
|
|
|
|
// File update
|
|
if (Pstream::master())
|
|
{
|
|
fileName faceZonesIntegrationDir;
|
|
if (Pstream::parRun())
|
|
{
|
|
// Put in undecomposed case (Note: gives problems for
|
|
// distributed data running)
|
|
faceZonesIntegrationDir =
|
|
obr_.time().path()/".."/name_/obr_.time().timeName();
|
|
}
|
|
else
|
|
{
|
|
faceZonesIntegrationDir =
|
|
obr_.time().path()/name_/obr_.time().timeName();
|
|
}
|
|
|
|
// Create directory if does not exist.
|
|
mkDir(faceZonesIntegrationDir);
|
|
|
|
// Open new file at start up
|
|
filePtr_.resize(fieldNames_.size());
|
|
|
|
forAll(fieldNames_, fieldI)
|
|
{
|
|
const word& fieldName = fieldNames_[fieldI];
|
|
|
|
OFstream* sPtr = new OFstream
|
|
(
|
|
faceZonesIntegrationDir/fieldName
|
|
);
|
|
|
|
filePtr_.insert(fieldName, sPtr);
|
|
}
|
|
|
|
// Add headers to output data
|
|
writeFileHeader();
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::faceZonesIntegration::writeFileHeader()
|
|
{
|
|
forAllIter(HashPtrTable<OFstream>, filePtr_, iter)
|
|
{
|
|
unsigned int w = IOstream::defaultPrecision() + 7;
|
|
|
|
OFstream& os = *filePtr_[iter.key()];
|
|
|
|
os << "#Time " << setw(w);
|
|
|
|
forAll (zoneNames_, zoneI)
|
|
{
|
|
os << zoneNames_[zoneI] << setw(w);
|
|
}
|
|
|
|
os << nl << endl;
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::faceZonesIntegration::execute()
|
|
{
|
|
// Do nothing - only valid on write
|
|
}
|
|
|
|
|
|
void Foam::faceZonesIntegration::end()
|
|
{
|
|
// Do nothing - only valid on write
|
|
}
|
|
|
|
|
|
void Foam::faceZonesIntegration::write()
|
|
{
|
|
if (active_)
|
|
{
|
|
makeFile();
|
|
|
|
forAll(fieldNames_, fieldI)
|
|
{
|
|
const word& fieldName = fieldNames_[fieldI];
|
|
|
|
const surfaceScalarField& sField =
|
|
obr_.lookupObject<surfaceScalarField>(fieldName);
|
|
|
|
const fvMesh& mesh = sField.mesh();
|
|
|
|
// 1. integrate over all face zones
|
|
|
|
scalarField integralVals(zoneNames_.size());
|
|
|
|
forAll(integralVals, zoneI)
|
|
{
|
|
const word& name = zoneNames_[zoneI];
|
|
|
|
label zoneID = mesh.faceZones().findZoneID(name);
|
|
|
|
const faceZone& fZone = mesh.faceZones()[zoneID];
|
|
|
|
integralVals[zoneI] = returnReduce
|
|
(
|
|
calcIntegral(sField, fZone),
|
|
sumOp<scalar>()
|
|
);
|
|
}
|
|
|
|
|
|
unsigned int w = IOstream::defaultPrecision() + 7;
|
|
|
|
// 2. Write only on master
|
|
|
|
if (Pstream::master() && filePtr_.found(fieldName))
|
|
{
|
|
OFstream& os = *filePtr_(fieldName);
|
|
|
|
os << obr_.time().value();
|
|
|
|
forAll(integralVals, zoneI)
|
|
{
|
|
os << ' ' << setw(w) << integralVals[zoneI];
|
|
|
|
if (log_)
|
|
{
|
|
Info<< "faceZonesIntegration output:" << nl
|
|
<< " Integration[" << zoneI << "] "
|
|
<< integralVals[zoneI] << endl;
|
|
}
|
|
}
|
|
|
|
os << endl;
|
|
}
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
Foam::scalar Foam::faceZonesIntegration::calcIntegral
|
|
(
|
|
const surfaceScalarField& sField,
|
|
const faceZone& fZone
|
|
) const
|
|
{
|
|
scalar sum = 0.0;
|
|
const fvMesh& mesh = sField.mesh();
|
|
|
|
forAll (fZone, i)
|
|
{
|
|
label faceI = fZone[i];
|
|
|
|
if (mesh.isInternalFace(faceI))
|
|
{
|
|
if (fZone.flipMap()[i])
|
|
{
|
|
sum -= sField[faceI];
|
|
}
|
|
else
|
|
{
|
|
sum += sField[faceI];
|
|
}
|
|
}
|
|
else
|
|
{
|
|
label patchI = mesh.boundaryMesh().whichPatch(faceI);
|
|
const polyPatch& pp = mesh.boundaryMesh()[patchI];
|
|
const fvsPatchScalarField& bField = sField.boundaryField()[patchI];
|
|
if (isA<processorPolyPatch>(pp))
|
|
{
|
|
if (refCast<const processorPolyPatch>(pp).owner())
|
|
{
|
|
label patchFaceI = pp.whichFace(faceI);
|
|
if (fZone.flipMap()[i])
|
|
{
|
|
sum -= bField[patchFaceI];
|
|
}
|
|
else
|
|
{
|
|
sum += bField[patchFaceI];
|
|
}
|
|
}
|
|
}
|
|
else if (isA<cyclicPolyPatch>(pp))
|
|
{
|
|
label patchFaceI = faceI - pp.start();
|
|
if (patchFaceI < pp.size()/2)
|
|
{
|
|
if (fZone.flipMap()[i])
|
|
{
|
|
sum -= bField[patchFaceI];
|
|
}
|
|
else
|
|
{
|
|
sum += bField[patchFaceI];
|
|
}
|
|
}
|
|
}
|
|
else if (!isA<emptyPolyPatch>(pp))
|
|
{
|
|
label patchFaceI = faceI - pp.start();
|
|
if (fZone.flipMap()[i])
|
|
{
|
|
sum -= bField[patchFaceI];
|
|
}
|
|
else
|
|
{
|
|
sum += bField[patchFaceI];
|
|
}
|
|
}
|
|
}
|
|
}
|
|
|
|
return sum;
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|