/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2019 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 .
\*---------------------------------------------------------------------------*/
// OpenFOAM includes
#include "surface.H"
#include "runTimePostProcessing.H"
#include "foamVtkTools.H"
#include "polySurface.H"
// VTK includes
#include "vtkMultiPieceDataSet.h"
#include "vtkPolyData.h"
#include "vtkSmartPointer.h"
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
vtkSmartPointer
Foam::functionObjects::runTimePostPro::surface::gatherSurfacePieces
(
const polySurface* surf
) const
{
auto multiPiece = vtkSmartPointer::New();
multiPiece->SetNumberOfPieces(Pstream::nProcs());
if (!needsCollective())
{
// Simple case (serial-serial, parallel-parallel)
if (surf)
{
multiPiece->SetPiece
(
Pstream::myProcNo(),
Foam::vtk::Tools::Patch::mesh(*surf)
);
}
}
else if (Pstream::master())
{
// Gather pieces on master
if (surf)
{
// Add myself
multiPiece->SetPiece
(
Pstream::myProcNo(),
Foam::vtk::Tools::Patch::mesh(*surf)
);
}
// Receive surfaces
pointField points;
faceList faces;
for
(
int slave=Pstream::firstSlave();
slave<=Pstream::lastSlave();
++slave
)
{
IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
points.clear();
faces.clear();
fromSlave >> points >> faces;
if (points.size())
{
multiPiece->SetPiece
(
slave,
Foam::vtk::Tools::Patch::mesh(points, faces)
);
}
}
}
else
{
// Slave - send surfaces
OPstream toMaster
(
Pstream::commsTypes::scheduled,
Pstream::masterNo()
);
if (surf)
{
toMaster
<< surf->points() << surf->faces();
}
else
{
toMaster
<< pointField() << faceList();
}
}
return multiPiece;
}
vtkSmartPointer
Foam::functionObjects::runTimePostPro::surface::gatherFaceCentres
(
const polySurface* surf
) const
{
auto multiPiece = vtkSmartPointer::New();
multiPiece->SetNumberOfPieces(Pstream::nProcs());
if (!needsCollective())
{
// Simple case
if (surf)
{
auto dataset = vtkSmartPointer::New();
auto geom = vtkSmartPointer::New();
geom->SetPoints(Foam::vtk::Tools::Patch::faceCentres(*surf));
geom->SetVerts(Foam::vtk::Tools::identityVertices(surf->nFaces()));
multiPiece->SetPiece(Pstream::myProcNo(), geom);
}
}
else if (Pstream::master())
{
// Gather pieces (face centres) on master
if (surf)
{
// Add myself
auto geom = vtkSmartPointer::New();
geom->SetPoints(Foam::vtk::Tools::Patch::faceCentres(*surf));
geom->SetVerts(Foam::vtk::Tools::identityVertices(surf->nFaces()));
multiPiece->SetPiece(Pstream::myProcNo(), geom);
}
// Receive points
pointField points;
for
(
int slave=Pstream::firstSlave();
slave<=Pstream::lastSlave();
++slave
)
{
IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
points.clear();
fromSlave >> points;
if (points.size())
{
multiPiece->SetPiece
(
slave,
Foam::vtk::Tools::Vertices(points)
);
}
}
}
else
{
// Slave - send face centres
OPstream toMaster
(
Pstream::commsTypes::scheduled,
Pstream::masterNo()
);
if (surf)
{
toMaster << surf->faceCentres();
}
else
{
toMaster << pointField();
}
}
return multiPiece;
}
// ************************************************************************* //