/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2011-2016 OpenFOAM Foundation Copyright (C) 2017-2019 OpenCFD Ltd. ------------------------------------------------------------------------------- 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 . \*---------------------------------------------------------------------------*/ #include "vtkPVFoam.H" // OpenFOAM includes #include "areaFaMesh.H" #include "Cloud.H" #include "IOobjectList.H" #include "vtkPVFoamReader.h" // VTK includes #include "vtkDataArraySelection.h" #include "vtkPolyData.h" #include "vtkUnstructuredGrid.h" // Templates (only needed here) #include "vtkPVFoamFieldTemplates.C" // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // void Foam::vtkPVFoam::convertVolFields() { const fvMesh& mesh = *volMeshPtr_; const bool interpFields = reader_->GetInterpolateVolFields(); wordHashSet selectedFields ( getSelected ( reader_->GetVolFieldSelection() ) ); if (selectedFields.empty()) { return; } // Get objects (fields) for this time - only keep selected fields // the region name is already in the mesh db IOobjectList objects(mesh, dbPtr_().timeName()); objects.filterKeys(selectedFields); if (objects.empty()) { return; } if (debug) { Info<< " " << FUNCTION_NAME << nl; forAllConstIters(objects, iter) { Info<< " " << (*iter)->name() << " == " << (*iter)->objectPath() << nl; } printMemory(); } PtrList interpLst; if (interpFields) { interpLst.setSize(mesh.boundaryMesh().size()); forAll(interpLst, i) { interpLst.set ( i, new PrimitivePatchInterpolation ( mesh.boundaryMesh()[i] ) ); } } convertVolFields(mesh, interpLst, objects); convertVolFields(mesh, interpLst, objects); convertVolFields(mesh, interpLst, objects); convertVolFields(mesh, interpLst, objects); convertVolFields(mesh, interpLst, objects); convertDimFields(mesh, interpLst, objects); convertDimFields(mesh, interpLst, objects); convertDimFields(mesh, interpLst, objects); convertDimFields(mesh, interpLst, objects); convertDimFields(mesh, interpLst, objects); if (debug) { Info<< " " << FUNCTION_NAME << nl; printMemory(); } } void Foam::vtkPVFoam::convertPointFields() { const fvMesh& mesh = *volMeshPtr_; wordHashSet selectedFields ( getSelected ( reader_->GetPointFieldSelection() ) ); if (selectedFields.empty()) { DebugInfo << "no point fields selected" << nl; return; } // Get objects (fields) for this time - only keep selected fields // the region name is already in the mesh db IOobjectList objects(mesh, dbPtr_().timeName()); objects.filterKeys(selectedFields); if (objects.empty()) { return; } if (debug) { Info<< " convert volume -> point fields" << nl; forAllConstIters(objects, iter) { Info<< " " << (*iter)->name() << " == " << (*iter)->objectPath() << nl; } printMemory(); } // Construct interpolation on the raw mesh const pointMesh& pMesh = pointMesh::New(mesh); convertPointFields(pMesh, objects); convertPointFields(pMesh, objects); convertPointFields(pMesh, objects); convertPointFields(pMesh, objects); convertPointFields(pMesh, objects); if (debug) { Info<< " convert volume -> point fields" << nl; printMemory(); } } void Foam::vtkPVFoam::convertAreaFields() { if (!areaMeshPtr_) { return; } const faMesh& mesh = *areaMeshPtr_; vtkDataArraySelection* select = reader_->GetVolFieldSelection(); wordHashSet selectedFields ( getSelected(select) ); if (selectedFields.empty()) { return; } // Get objects (fields) for this time - only keep selected fields // the region name is already in the mesh db IOobjectList objects(mesh.mesh(), dbPtr_().timeName()); objects.filterKeys(selectedFields); if (objects.empty()) { return; } if (debug) { Info<< " " << FUNCTION_NAME << nl; forAllConstIters(objects, iter) { Info<< " " << (*iter)->name() << " == " << (*iter)->objectPath() << nl; } printMemory(); } convertAreaFields(mesh, objects); convertAreaFields(mesh, objects); convertAreaFields(mesh, objects); convertAreaFields(mesh, objects); convertAreaFields(mesh, objects); if (debug) { Info<< " " << FUNCTION_NAME << nl; printMemory(); } } void Foam::vtkPVFoam::convertLagrangianFields() { const List