vtkPVblockMeshConvert.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "vtkPVblockMesh.H"
27 #include "vtkPVblockMeshReader.h"
28 #include "vtkOpenFOAMPoints.H"
29 
30 // OpenFOAM includes
31 #include "blockMesh.H"
32 #include "Time.H"
33 
34 // VTK includes
35 #include "vtkCellArray.h"
36 #include "vtkDataArraySelection.h"
37 #include "vtkMultiBlockDataSet.h"
38 #include "vtkPoints.h"
39 #include "vtkPolyData.h"
40 #include "vtkUnstructuredGrid.h"
41 
42 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
43 
44 void Foam::vtkPVblockMesh::convertMeshBlocks
45 (
46  vtkMultiBlockDataSet* output,
47  int& blockNo
48 )
49 {
50  vtkDataArraySelection* selection = reader_->GetBlockSelection();
51  arrayRange& range = arrayRangeBlocks_;
52  range.block(blockNo); // Set output block
53  label datasetNo = 0; // Restart at dataset 0
54 
55  const blockMesh& blkMesh = *meshPtr_;
56  const Foam::pointField& blockPoints = blkMesh.vertices();
57 
58  int blockI = 0;
59  const scalar scaleFactor = blkMesh.scaleFactor();
60 
61  for (int partId = range.start(); partId < range.end(); ++partId, ++blockI)
62  {
63  if (!blockStatus_[partId]) continue;
64 
65  const blockDescriptor& blockDef = blkMesh[blockI];
66 
67  vtkUnstructuredGrid* vtkmesh = vtkUnstructuredGrid::New();
68 
69  // Convert OpenFOAM mesh vertices to VTK
70  vtkPoints *vtkpoints = vtkPoints::New();
71  vtkpoints->Allocate(blockDef.nPoints());
72  const labelList& blockLabels = blockDef.blockShape();
73 
74  vtkmesh->Allocate(1);
75  vtkIdType nodeIds[8];
76 
77  forAll(blockLabels, ptI)
78  {
80  (
81  vtkpoints,
82  blockPoints[blockLabels[ptI]],
83  scaleFactor
84  );
85 
86  nodeIds[ptI] = ptI;
87  }
88 
89  vtkmesh->InsertNextCell
90  (
91  VTK_HEXAHEDRON,
92  8,
93  nodeIds
94  );
95 
96  vtkmesh->SetPoints(vtkpoints);
97  vtkpoints->Delete();
98 
99  AddToBlock
100  (
101  output, vtkmesh, range, datasetNo,
102  selection->GetArrayName(partId)
103  );
104 
105  vtkmesh->Delete();
106  datasetNo++;
107  }
108 
109  // Anything added?
110  if (datasetNo) ++blockNo;
111 }
112 
113 
114 void Foam::vtkPVblockMesh::convertMeshEdges
115 (
116  vtkMultiBlockDataSet* output,
117  int& blockNo
118 )
119 {
120  vtkDataArraySelection* selection = reader_->GetCurvedEdgesSelection();
121  arrayRange& range = arrayRangeEdges_;
122 
123  range.block(blockNo); // Set output block
124  label datasetNo = 0; // Restart at dataset 0
125 
126  const blockMesh& blkMesh = *meshPtr_;
127  const blockEdgeList& edges = blkMesh.edges();
128 
129  int edgeI = 0;
130  const scalar scaleFactor = blkMesh.scaleFactor();
131 
132  for (int partId = range.start(); partId < range.end(); ++partId, ++edgeI)
133  {
134  if (!edgeStatus_[partId]) continue;
135 
136  // Search each block
137  forAll(blkMesh, blockI)
138  {
139  const blockDescriptor& blockDef = blkMesh[blockI];
140 
141  edgeList blkEdges = blockDef.blockShape().edges();
142 
143  // List of edge point and weighting factors
144  pointField edgesPoints[12];
145  scalarList edgesWeights[12];
146  blockDef.edgesPointsWeights(edgesPoints, edgesWeights);
147 
148  // Find the corresponding edge within the block
149  label foundEdgeI = -1;
150  forAll(blkEdges, blkEdgeI)
151  {
152  if (edges[edgeI].compare(blkEdges[blkEdgeI]))
153  {
154  foundEdgeI = blkEdgeI;
155  break;
156  }
157  }
158 
159  if (foundEdgeI != -1)
160  {
161  const List<point>& edgePoints = edgesPoints[foundEdgeI];
162 
163  vtkPolyData* vtkmesh = vtkPolyData::New();
164  vtkPoints* vtkpoints = vtkPoints::New();
165 
166  vtkpoints->Allocate(edgePoints.size());
167  vtkmesh->Allocate(1);
168 
169  List<vtkIdType> pointIds(edgePoints.size());
170  forAll(edgePoints, ptI)
171  {
173  (
174  vtkpoints,
175  edgePoints[ptI],
176  scaleFactor
177  );
178  pointIds[ptI] = ptI;
179  }
180 
181  vtkmesh->InsertNextCell
182  (
183  VTK_POLY_LINE,
184  edgePoints.size(),
185  pointIds.data()
186  );
187 
188  vtkmesh->SetPoints(vtkpoints);
189  vtkpoints->Delete();
190 
191  AddToBlock
192  (
193  output, vtkmesh, range, datasetNo,
194  selection->GetArrayName(partId)
195  );
196 
197  vtkmesh->Delete();
198  datasetNo++;
199 
200  break;
201  }
202  }
203  }
204 
205  // Anything added?
206  if (datasetNo) ++blockNo;
207 }
208 
209 
210 void Foam::vtkPVblockMesh::convertMeshCorners
211 (
212  vtkMultiBlockDataSet* output,
213  int& blockNo
214 )
215 {
216  arrayRange& range = arrayRangeCorners_;
217  range.block(blockNo); // Set output block
218  label datasetNo = 0; // Restart at dataset 0
219 
220  const pointField& blockPoints = meshPtr_->vertices();
221  const scalar& scaleFactor = meshPtr_->scaleFactor();
222 
223  vtkPolyData* vtkmesh = vtkPolyData::New();
224  vtkPoints* vtkpoints = vtkPoints::New();
225  vtkCellArray* vtkcells = vtkCellArray::New();
226 
227  vtkpoints->Allocate(blockPoints.size());
228  vtkcells->Allocate(blockPoints.size());
229 
230  vtkIdType pointId = 0;
231  forAll(blockPoints, ptI)
232  {
234  (
235  vtkpoints,
236  blockPoints[ptI],
237  scaleFactor
238  );
239 
240  vtkcells->InsertNextCell(1, &pointId);
241  pointId++;
242  }
243 
244  vtkmesh->SetPoints(vtkpoints);
245  vtkpoints->Delete();
246 
247  vtkmesh->SetVerts(vtkcells);
248  vtkcells->Delete();
249 
250  AddToBlock
251  (
252  output, vtkmesh, range, datasetNo,
253  arrayRangeCorners_.name()
254  );
255  vtkmesh->Delete();
256 
257  datasetNo++;
258 
259  // Anything added?
260  if (datasetNo) ++blockNo;
261 }
262 
263 
264 // ************************************************************************* //
scalar range
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
PtrList< blockEdge > blockEdgeList
A PtrList of blockEdges.
Definition: blockEdgeList.H:45
List< edge > edgeList
Definition: edgeList.H:38
void vtkInsertNextOpenFOAMPoint(vtkPoints *points, const Foam::point &p)