Loading...
Searching...
No Matches
foamVtkVtuAdaptor.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::vtk::vtuAdaptor
28
29Description
30 A low-level backend adaptor for converting OpenFOAM volume meshes/fields
31 to VTK internal representation.
32 The output is a single piece vtkUnstructuredGrid dataset.
33 Multiblock composition is done by the caller.
34
35Note
36 This was originally part of the catalystFvMesh function object backend,
37 which was prototyped in the ParaView reader module.
38 This adaptor unifies many common elements of both.
39
40SourceFiles
41 foamVtkVtuAdaptorI.H
42 foamVtkVtuAdaptor.txx
43
44\*---------------------------------------------------------------------------*/
45
46#ifndef Foam_vtk_vtuAdaptor_H
47#define Foam_vtk_vtuAdaptor_H
48
49#include "fvMesh.H"
50#include "volFieldsFwd.H"
51#include "foamVtkTools.H"
52#include "foamVtkMeshMaps.H"
53#include "foamVtuSizing.H"
54
55#include "vtkUnstructuredGrid.h"
56#include "vtkMultiBlockDataSet.h"
57
58// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59
60// Forward Declarations
61class vtkDataSet;
62class vtkFloatArray;
63class vtkDoubleArray;
64
65namespace Foam
66{
67namespace vtk
68{
69
70/*---------------------------------------------------------------------------*\
71 Class vtk::vtuAdaptor Declaration
72\*---------------------------------------------------------------------------*/
73
74//- Bookkeeping for vtkUnstructuredGrid
75class vtuAdaptor
76:
77 public vtk::Caching<vtkUnstructuredGrid>,
78 public foamVtkMeshMaps
79{
80public:
81
82 // Member Functions
83
84 // Mesh Conversion
85
86 //- Allocate vtkPoints
87 inline vtkSmartPointer<vtkPoints> NewPoints
88 (
89 vtkIdType count
90 ) const;
91
92 //- The vtk points for the mesh (and decomposition)
93 inline vtkSmartPointer<vtkPoints> points
94 (
95 const polyMesh& mesh
96 ) const;
97
98 //- The vtk points for the mesh (and decomposition)
99 //- using the provided pointMap
100 inline vtkSmartPointer<vtkPoints> points
101 (
102 const polyMesh& mesh,
103 const labelUList& pointMap
104 ) const;
105
106 //- Internal mesh as vtkUnstructuredGrid.
107 // Calling this automatically populates the sizing and vtk::Caching
108 inline vtkSmartPointer<vtkUnstructuredGrid> internal
109 (
110 const fvMesh& mesh,
112 const bool decompPoly = false
113 );
114
115
116 // Field Conversion
117
118 //- Convert internal volume field (CellData)
119 template<class Type, class DataArrayType = vtkFloatArray>
120 vtkSmartPointer<DataArrayType>
122 (
124 ) const;
125
126 //- Convert internal volume field (CellData)
127 template<class Type, class DataArrayType = vtkFloatArray>
128 vtkSmartPointer<DataArrayType>
130 (
132 ) const;
133};
134
135
136// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137
138} // End namespace vtk
139} // End namespace Foam
140
141// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
142
143#include "foamVtkVtuAdaptorI.H"
144
145// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146
147#ifdef NoRepository
148 #include "foamVtkVtuAdaptor.txx"
149#endif
150
151// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152
153#endif
154
155// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
foamVtkMeshMaps()=default
Default construct: zero-sized, no reserved sizes.
const labelList & pointMap() const noexcept
Point labels for subsetted meshes.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Bookkeeping for vtkUnstructuredGrid.
vtkSmartPointer< vtkPoints > NewPoints(vtkIdType count) const
Allocate vtkPoints.
vtkSmartPointer< DataArrayType > convertField(const GeometricField< Type, fvPatchField, volMesh > &fld) const
Convert internal volume field (CellData).
vtkSmartPointer< DataArrayType > convertField(const DimensionedField< Type, volMesh > &fld) const
Convert internal volume field (CellData).
vtkSmartPointer< vtkUnstructuredGrid > internal(const fvMesh &mesh, const bool decompPoly=false)
Internal mesh as vtkUnstructuredGrid.
dynamicFvMesh & mesh
const pointField & points
Namespace for handling VTK output. Contains classes and functions for writing VTK file content.
Namespace for OpenFOAM.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
Bookkeeping for internal caching.
Forwards and collection of common volume field types.