foamVtkCloudAdaptor.C 4.34 KB
Newer Older
1 2 3 4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
OpenFOAM bot's avatar
OpenFOAM bot committed
5
    \\  /    A nd           | www.openfoam.com
6
     \\/     M anipulation  |
OpenFOAM bot's avatar
OpenFOAM bot committed
7
-------------------------------------------------------------------------------
8
    Copyright (C) 2018-2020 OpenCFD Ltd.
9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41
-------------------------------------------------------------------------------
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 "foamVtkCloudAdaptor.H"
#include "Cloud.H"
#include "IOField.H"
#include "predicates.H"

// Only used locally
#include "foamVtkCloudAdaptorTemplates.C"

// VTK includes
#include <vtkSmartPointer.h>
#include <vtkMultiPieceDataSet.h>

// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //

42
namespace Foam
43
{
44 45 46 47
namespace vtk
{
    defineTypeNameAndDebug(cloudAdaptor, 0);
}
48
} // End namespace Foam
49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75


// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //

template<class UnaryMatchPredicate>
vtkSmartPointer<vtkMultiPieceDataSet>
Foam::vtk::cloudAdaptor::getCloudImpl
(
    const objectRegistry& mesh,
    const word& cloudName,
    const UnaryMatchPredicate& matcher
)
{
    // All individual datasets are vtkMultiPieceDataSet for improved
    // handling downstream.

    label rank = 0;
    label nproc = 1;

    if (Pstream::parRun())
    {
        rank  = Pstream::myProcNo();
        nproc = Pstream::nProcs();
    }

    auto output = vtkSmartPointer<vtkMultiPieceDataSet>::New();

76
    const auto* objPtr = mesh.cfindObject<cloud>(cloudName);
77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96
    if (!objPtr)
    {
        return output;
    }

    objectRegistry obrTmp
    (
        IOobject
        (
            "vtk::catalyst::" + cloudName,
            mesh.time().constant(),
            mesh,
            IOobject::NO_READ,
            IOobject::NO_WRITE,
            false
        )
    );

    objPtr->writeObjects(obrTmp);

97
    const auto* pointsPtr = obrTmp.cfindObject<vectorField>("position");
98 99 100 101

    vtkSmartPointer<vtkPolyData> vtkmesh;
    if (pointsPtr)
    {
102
        vtkmesh = vtk::Tools::Vertices(*pointsPtr);
103

104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120
        // Prevent any possible conversion of positions as a field
        obrTmp.filterKeys
        (
            [](const word& k)
            {
                return k.startsWith("position")
                    || k.startsWith("coordinate");
            },
            true  // prune
        );

        // Restrict to specified fields
        obrTmp.filterKeys(matcher);

        convertLagrangianFields<label>(vtkmesh, obrTmp);
        convertLagrangianFields<scalar>(vtkmesh, obrTmp);
        convertLagrangianFields<vector>(vtkmesh, obrTmp);
121 122 123
    }
    else
    {
124
        // This should never occur
125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171
        vtkmesh = vtkSmartPointer<vtkPolyData>::New();
    }

    output->SetNumberOfPieces(nproc);
    output->SetPiece(rank, vtkmesh);

    return output;
}


// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

Foam::vtk::cloudAdaptor::cloudAdaptor(const fvMesh& mesh)
:
    mesh_(mesh)
{}


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

vtkSmartPointer<vtkMultiPieceDataSet>
Foam::vtk::cloudAdaptor::getCloud
(
    const word& cloudName
) const
{
    return getCloudImpl(mesh_, cloudName, predicates::always());
}


vtkSmartPointer<vtkMultiPieceDataSet>
Foam::vtk::cloudAdaptor::getCloud
(
    const word& cloudName,
    const wordRes& selectFields
) const
{
    if (selectFields.empty())
    {
        return getCloud(cloudName);
    }

    return getCloudImpl(mesh_, cloudName, selectFields);
}


// ************************************************************************* //