Files
openfoam/src/functionObjects/graphics/runTimePostProcessing/geometryCloudTemplates.C
Mark Olesen 42fbf6d38c ENH: extended runTimePostProcessing (#1206)
- Extended runTimePostProcessing to include access to "live"
  simulation objects such a geometry patches and sampled surfaces
  stored on the "functionObjectObjects" registry.

- Add 'live' runTimePostProcessing of cloud data.
  Extracts position and fields from the cloud via its objectRegistry writer

- For the "live" simulation objects, there are two new volume filters
  that work directly with the OpenFOAM volume fields:
      * iso-surface
      * cutting planes
  Both use the VTK algorithms directly and support multiple values.
  Eg, can make multiple iso-levels or multiple planes parallel to each
  other.

- When VTK has been compiled with MPI-support, parallel rendering will
  be used.

- Additional title text properties (shadow, italic etc)

- Simplified handling of scalar-bar and visibility switches

- Support multiple text positions. Eg, for adding watermark text.
2019-02-13 11:22:46 +01:00

169 lines
3.9 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / 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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "foamVtkTools.H"
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class Type>
bool Foam::functionObjects::runTimePostPro::geometryCloud::addField
(
vtkDataSet* piece,
const Field<Type>& fld,
const word& fieldName
) const
{
if (!piece) return false;
auto vtkfield = Foam::vtk::Tools::convertFieldToVTK<Type>(fieldName, fld);
// Only has verts
piece->GetPointData()->AddArray(vtkfield);
return true;
}
template<class Type>
bool Foam::functionObjects::runTimePostPro::geometryCloud::addCloudField
(
vtkMultiPieceDataSet* multiPiece,
const IOField<Type>* fldptr,
const word& fieldName
) const
{
if (!multiPiece)
{
return false;
}
if (!needsCollective())
{
// Simple case (serial-serial, parallel-parallel)
return fldptr &&
addField<Type>
(
multiPiece->GetPiece(Pstream::myProcNo()),
*fldptr,
fieldName
);
}
// Gather fields
const bool ok = returnReduce((fldptr != nullptr), orOp<bool>());
if (!ok)
{
return false;
}
if (Pstream::master())
{
if (fldptr)
{
// My field data
addField<Type>
(
multiPiece->GetPiece(Pstream::myProcNo()),
*fldptr,
fieldName
);
}
// Receive field data
Field<Type> recv;
for
(
int slave=Pstream::firstSlave();
slave<=Pstream::lastSlave();
++slave
)
{
IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
recv.clear();
fromSlave
>> recv;
if (recv.size())
{
addField<Type>
(
multiPiece->GetPiece(slave),
recv,
fieldName
);
}
}
}
else
{
// Slave - send field data
OPstream toMaster
(
Pstream::commsTypes::scheduled,
Pstream::masterNo()
);
if (fldptr)
{
toMaster
<< *fldptr;
}
else
{
toMaster
<< List<Type>();
}
}
return ok;
}
template<class Type>
bool Foam::functionObjects::runTimePostPro::geometryCloud::addCloudField
(
vtkMultiPieceDataSet* multiPiece,
const regIOobject* ioptr,
const word& fieldName
) const
{
return addCloudField<Type>
(
multiPiece,
dynamic_cast<const IOField<Type>*>(ioptr),
fieldName
);
}
// ************************************************************************* //