210 lines
5.9 KiB
C
210 lines
5.9 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2011 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/>.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "sampledSurfaces.H"
|
|
#include "volFields.H"
|
|
#include "surfaceFields.H"
|
|
#include "ListListOps.H"
|
|
|
|
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
|
|
|
template<class Type>
|
|
void Foam::sampledSurfaces::writeSurface
|
|
(
|
|
const Field<Type>& values,
|
|
const label surfI,
|
|
const word& fieldName,
|
|
const fileName& outputDir
|
|
)
|
|
{
|
|
const sampledSurface& s = operator[](surfI);
|
|
|
|
if (Pstream::parRun())
|
|
{
|
|
// Collect values from all processors
|
|
List<Field<Type> > gatheredValues(Pstream::nProcs());
|
|
gatheredValues[Pstream::myProcNo()] = values;
|
|
Pstream::gatherList(gatheredValues);
|
|
|
|
if (Pstream::master())
|
|
{
|
|
// Combine values into single field
|
|
Field<Type> allValues
|
|
(
|
|
ListListOps::combine<Field<Type> >
|
|
(
|
|
gatheredValues,
|
|
accessOp<Field<Type> >()
|
|
)
|
|
);
|
|
|
|
// Renumber (point data) to correspond to merged points
|
|
if (mergeList_[surfI].pointsMap.size() == allValues.size())
|
|
{
|
|
inplaceReorder(mergeList_[surfI].pointsMap, allValues);
|
|
allValues.setSize(mergeList_[surfI].points.size());
|
|
}
|
|
|
|
// Write to time directory under outputPath_
|
|
// skip surface without faces (eg, a failed cut-plane)
|
|
if (mergeList_[surfI].faces.size())
|
|
{
|
|
formatter_->write
|
|
(
|
|
outputDir,
|
|
s.name(),
|
|
mergeList_[surfI].points,
|
|
mergeList_[surfI].faces,
|
|
fieldName,
|
|
allValues,
|
|
s.interpolate()
|
|
);
|
|
}
|
|
}
|
|
}
|
|
else
|
|
{
|
|
// Write to time directory under outputPath_
|
|
// skip surface without faces (eg, a failed cut-plane)
|
|
if (s.faces().size())
|
|
{
|
|
formatter_->write
|
|
(
|
|
outputDir,
|
|
s.name(),
|
|
s.points(),
|
|
s.faces(),
|
|
fieldName,
|
|
values,
|
|
s.interpolate()
|
|
);
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
template<class Type>
|
|
void Foam::sampledSurfaces::sampleAndWrite
|
|
(
|
|
const GeometricField<Type, fvPatchField, volMesh>& vField
|
|
)
|
|
{
|
|
// interpolator for this field
|
|
autoPtr< interpolation<Type> > interpolator;
|
|
|
|
const word& fieldName = vField.name();
|
|
const fileName outputDir = outputPath_/vField.time().timeName();
|
|
|
|
forAll(*this, surfI)
|
|
{
|
|
const sampledSurface& s = operator[](surfI);
|
|
|
|
Field<Type> values;
|
|
|
|
if (s.interpolate())
|
|
{
|
|
if (interpolator.empty())
|
|
{
|
|
interpolator = interpolation<Type>::New
|
|
(
|
|
interpolationScheme_,
|
|
vField
|
|
);
|
|
}
|
|
|
|
values = s.interpolate(interpolator());
|
|
}
|
|
else
|
|
{
|
|
values = s.sample(vField);
|
|
}
|
|
|
|
writeSurface<Type>(values, surfI, fieldName, outputDir);
|
|
}
|
|
}
|
|
|
|
|
|
|
|
template<class Type>
|
|
void Foam::sampledSurfaces::sampleAndWrite
|
|
(
|
|
const GeometricField<Type, fvsPatchField, surfaceMesh>& sField
|
|
)
|
|
{
|
|
const word& fieldName = sField.name();
|
|
const fileName outputDir = outputPath_/sField.time().timeName();
|
|
|
|
forAll(*this, surfI)
|
|
{
|
|
const sampledSurface& s = operator[](surfI);
|
|
Field<Type> values = s.sample(sField);
|
|
writeSurface<Type>(values, surfI, fieldName, outputDir);
|
|
}
|
|
}
|
|
|
|
|
|
template<class GeoField>
|
|
void Foam::sampledSurfaces::sampleAndWrite(const IOobjectList& allObjects)
|
|
{
|
|
IOobjectList fields = allObjects.lookupClass(GeoField::typeName);
|
|
forAllConstIter(IOobjectList, fields, fieldIter)
|
|
{
|
|
forAll (fieldSelection_, fieldI)
|
|
{
|
|
const wordRe field = fieldSelection_[fieldI];
|
|
if (field.match(fieldIter()->name()))
|
|
{
|
|
if (Pstream::master() && verbose_)
|
|
{
|
|
Pout<< "sampleAndWrite: " << field << endl;
|
|
}
|
|
|
|
if (loadFromFiles_)
|
|
{
|
|
fieldIter()->readOpt() = IOobject::MUST_READ;
|
|
sampleAndWrite
|
|
(
|
|
GeoField
|
|
(
|
|
*fieldIter(),
|
|
mesh_
|
|
)
|
|
);
|
|
}
|
|
else
|
|
{
|
|
sampleAndWrite
|
|
(
|
|
mesh_.lookupObject<GeoField>(fieldIter()->name())
|
|
);
|
|
}
|
|
}
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|