/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | foam-extend: Open Source CFD \\ / O peration | \\ / A nd | For copyright notice see file Copyright \\/ M anipulation | ------------------------------------------------------------------------------- License This file is part of foam-extend. foam-extend 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. foam-extend 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 foam-extend. If not, see . InClass Foam::vtkFoam \*---------------------------------------------------------------------------*/ #ifndef vtkFoamConvertPointField_H #define vtkFoamConvertPointField_H // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // #include "interpolatePointToCell.H" template void Foam::vtkFoam::convertPointField ( const GeometricField& ptf, const GeometricField& tf ) { vtkUnstructuredGrid *vtkMesh = vtkUnstructuredGrid::SafeDownCast(reader_->GetOutput(0)); vtkFloatArray *pointTypes = vtkFloatArray::New(); pointTypes->SetNumberOfTuples(ptf.size() + addPointCellLabels_.size()); pointTypes->SetNumberOfComponents(Type::nComponents); pointTypes->Allocate(Type::nComponents*ptf.size()); pointTypes->SetName(ptf.name().c_str()); float vec[Type::nComponents]; forAll(ptf, i) { for (direction d=0; dInsertTuple(i, vec); } label i = ptf.size(); if (&tf != &GeometricField::null()) { forAll(addPointCellLabels_, api) { Type t = tf[addPointCellLabels_[api]]; for (direction d=0; dInsertTuple(i++, vec); } } else { forAll(addPointCellLabels_, api) { Type t = interpolatePointToCell(ptf, addPointCellLabels_[api]); for (direction d=0; dInsertTuple(i++, vec); } } vtkMesh->GetPointData()->AddArray(pointTypes); pointTypes->Delete(); } template<> void Foam::vtkFoam::convertPointField ( const GeometricField& psf, const GeometricField& sf ) { vtkUnstructuredGrid *vtkMesh = vtkUnstructuredGrid::SafeDownCast(reader_->GetOutput(0)); vtkFloatArray *pointScalars = vtkFloatArray::New(); pointScalars->SetNumberOfTuples(psf.size() + addPointCellLabels_.size()); pointScalars->SetNumberOfComponents(1); pointScalars->Allocate(psf.size()); pointScalars->SetName(psf.name().c_str()); for (int i=0; iInsertComponent(i, 0, psf[i]); } label i = psf.size(); if (&sf != &GeometricField::null()) { forAll(addPointCellLabels_, api) { pointScalars->InsertComponent ( i++, 0, sf[addPointCellLabels_[api]] ); } } else { forAll(addPointCellLabels_, api) { pointScalars->InsertComponent ( i++, 0, interpolatePointToCell(psf, addPointCellLabels_[api]) ); } } vtkMesh->GetPointData()->AddArray(pointScalars); if (!vtkMesh->GetPointData()->GetScalars()) { vtkMesh->GetPointData()->SetScalars(pointScalars); } pointScalars->Delete(); } // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // #endif // ************************************************************************* //