vtkPVblockMeshConvert.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 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 
29 // OpenFOAM includes
30 #include "blockMesh.H"
31 #include "Time.H"
32 
33 #include "vtkOpenFOAMPoints.H"
34 
35 // VTK includes
36 #include "vtkCellArray.h"
37 #include "vtkDataArraySelection.h"
38 #include "vtkMultiBlockDataSet.h"
39 #include "vtkPoints.h"
40 #include "vtkPolyData.h"
41 #include "vtkUnstructuredGrid.h"
42 
43 
44 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
45 
46 void Foam::vtkPVblockMesh::convertMeshBlocks
47 (
48  vtkMultiBlockDataSet* output,
49  int& blockNo
50 )
51 {
52  vtkDataArraySelection* selection = reader_->GetBlockSelection();
53  arrayRange& range = arrayRangeBlocks_;
54  range.block(blockNo); // set output block
55  label datasetNo = 0; // restart at dataset 0
56 
57  const blockMesh& blkMesh = *meshPtr_;
58  const Foam::pointField& blockPoints = blkMesh.blockPointField();
59 
60  if (debug)
61  {
62  Info<< "<beg> Foam::vtkPVblockMesh::convertMeshBlocks" << endl;
63  }
64 
65  int blockI = 0;
66  const scalar scaleFactor = blkMesh.scaleFactor();
67 
68  for
69  (
70  int partId = range.start();
71  partId < range.end();
72  ++partId, ++blockI
73  )
74  {
75  if (!blockStatus_[partId])
76  {
77  continue;
78  }
79 
80  const blockDescriptor& blockDef = blkMesh[blockI].blockDef();
81 
82  vtkUnstructuredGrid* vtkmesh = vtkUnstructuredGrid::New();
83 
84  // Convert OpenFOAM mesh vertices to VTK
85  vtkPoints *vtkpoints = vtkPoints::New();
86  vtkpoints->Allocate( blockDef.nPoints() );
87  const labelList& blockLabels = blockDef.blockShape();
88 
89  vtkmesh->Allocate(1);
90  vtkIdType nodeIds[8];
91 
92  forAll(blockLabels, ptI)
93  {
95  (
96  vtkpoints,
97  blockPoints[blockLabels[ptI]],
98  scaleFactor
99  );
100 
101  nodeIds[ptI] = ptI;
102  }
103 
104  vtkmesh->InsertNextCell
105  (
106  VTK_HEXAHEDRON,
107  8,
108  nodeIds
109  );
110 
111  vtkmesh->SetPoints(vtkpoints);
112  vtkpoints->Delete();
113 
114  AddToBlock
115  (
116  output, vtkmesh, range, datasetNo,
117  selection->GetArrayName(partId)
118  );
119 
120  vtkmesh->Delete();
121  datasetNo++;
122  }
123 
124 
125  // anything added?
126  if (datasetNo)
127  {
128  ++blockNo;
129  }
130 
131  if (debug)
132  {
133  Info<< "<end> Foam::vtkPVblockMesh::convertMeshBlocks" << endl;
134  }
135 }
136 
137 
138 void Foam::vtkPVblockMesh::convertMeshEdges
139 (
140  vtkMultiBlockDataSet* output,
141  int& blockNo
142 )
143 {
144  vtkDataArraySelection* selection = reader_->GetCurvedEdgesSelection();
145  arrayRange& range = arrayRangeEdges_;
146 
147  range.block(blockNo); // set output block
148  label datasetNo = 0; // restart at dataset 0
149 
150  const blockMesh& blkMesh = *meshPtr_;
151  const curvedEdgeList& edges = blkMesh.edges();
152 
153  int edgeI = 0;
154  const scalar scaleFactor = blkMesh.scaleFactor();
155 
156  for
157  (
158  int partId = range.start();
159  partId < range.end();
160  ++partId, ++edgeI
161  )
162  {
163  if (!edgeStatus_[partId])
164  {
165  continue;
166  }
167 
168  // search each block
169  forAll(blkMesh, blockI)
170  {
171  const blockDescriptor& blockDef = blkMesh[blockI].blockDef();
172 
173  edgeList blkEdges = blockDef.blockShape().edges();
174 
175  // find the corresponding edge within the block
176  label foundEdgeI = -1;
177  forAll(blkEdges, blkEdgeI)
178  {
179  if (edges[edgeI].compare(blkEdges[blkEdgeI]))
180  {
181  foundEdgeI = blkEdgeI;
182  break;
183  }
184  }
185 
186  if (foundEdgeI != -1)
187  {
188  const List<point>& edgePoints =
189  blockDef.blockEdgePoints()[foundEdgeI];
190 
191 
192  vtkPolyData* vtkmesh = vtkPolyData::New();
193  vtkPoints* vtkpoints = vtkPoints::New();
194 
195  vtkpoints->Allocate( edgePoints.size() );
196  vtkmesh->Allocate(1);
197 
198  vtkIdType pointIds[edgePoints.size()];
199  forAll(edgePoints, ptI)
200  {
202  (
203  vtkpoints,
204  edgePoints[ptI],
205  scaleFactor
206  );
207  pointIds[ptI] = ptI;
208  }
209 
210  vtkmesh->InsertNextCell
211  (
212  VTK_POLY_LINE,
213  edgePoints.size(),
214  pointIds
215  );
216 
217  vtkmesh->SetPoints(vtkpoints);
218  vtkpoints->Delete();
219 
220  AddToBlock
221  (
222  output, vtkmesh, range, datasetNo,
223  selection->GetArrayName(partId)
224  );
225 
226  vtkmesh->Delete();
227  datasetNo++;
228 
229  break;
230  }
231  }
232  }
233 
234 
235  // anything added?
236  if (datasetNo)
237  {
238  ++blockNo;
239  }
240 
241  if (debug)
242  {
243  Info<< "<end> Foam::vtkPVblockMesh::convertMeshEdges" << endl;
244  }
245 
246 }
247 
248 
249 void Foam::vtkPVblockMesh::convertMeshCorners
250 (
251  vtkMultiBlockDataSet* output,
252  int& blockNo
253 )
254 {
255  arrayRange& range = arrayRangeCorners_;
256  range.block(blockNo); // set output block
257  label datasetNo = 0; // restart at dataset 0
258 
259  const pointField& blockPoints = meshPtr_->blockPointField();
260  const scalar& scaleFactor = meshPtr_->scaleFactor();
261 
262  if (debug)
263  {
264  Info<< "<beg> Foam::vtkPVblockMesh::convertMeshCorners" << endl;
265  }
266 
267  if (true) // or some flag or other condition
268  {
269  vtkPolyData* vtkmesh = vtkPolyData::New();
270  vtkPoints* vtkpoints = vtkPoints::New();
271  vtkCellArray* vtkcells = vtkCellArray::New();
272 
273  vtkpoints->Allocate( blockPoints.size() );
274  vtkcells->Allocate( blockPoints.size() );
275 
276  vtkIdType pointId = 0;
277  forAll(blockPoints, ptI)
278  {
280  (
281  vtkpoints,
282  blockPoints[ptI],
283  scaleFactor
284  );
285 
286  vtkcells->InsertNextCell(1, &pointId);
287  pointId++;
288  }
289 
290  vtkmesh->SetPoints(vtkpoints);
291  vtkpoints->Delete();
292 
293  vtkmesh->SetVerts(vtkcells);
294  vtkcells->Delete();
295 
296  AddToBlock
297  (
298  output, vtkmesh, range, datasetNo,
299  arrayRangeCorners_.name()
300  );
301  vtkmesh->Delete();
302 
303  datasetNo++;
304  }
305 
306  // anything added?
307  if (datasetNo)
308  {
309  ++blockNo;
310  }
311 
312  if (debug)
313  {
314  Info<< "<end> Foam::vtkPVblockMesh::convertMeshCorners" << endl;
315  }
316 }
317 
318 
319 // ************************************************************************* //
scalar scaleFactor() const
The scaling factor used to convert to metres.
Definition: blockMesh.C:112
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
PtrList< curvedEdge > curvedEdgeList
A PtrList of curvedEdges.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
void vtkInsertNextOpenFOAMPoint(vtkPoints *points, const Foam::point &p)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< edge > edgeList
Definition: edgeList.H:38
const pointField & blockPointField() const
Reference to point field defining the block mesh.
Definition: blockMesh.C:76
List< label > labelList
A List of labels.
Definition: labelList.H:56
messageStream Info