2 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
4 // + This file is part of enGrid. +
6 // + Copyright 2008,2009 Oliver Gloth +
8 // + enGrid is free software: you can redistribute it and/or modify +
9 // + it under the terms of the GNU General Public License as published by +
10 // + the Free Software Foundation, either version 3 of the License, or +
11 // + (at your option) any later version. +
13 // + enGrid is distributed in the hope that it will be useful, +
14 // + but WITHOUT ANY WARRANTY; without even the implied warranty of +
15 // + MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the +
16 // + GNU General Public License for more details. +
18 // + You should have received a copy of the GNU General Public License +
19 // + along with enGrid. If not, see <http://www.gnu.org/licenses/>. +
21 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
29 #include <vtkUnstructuredGrid.h>
30 #include <vtkPolyData.h>
31 #include <vtkPointData.h>
32 #include <vtkCellData.h>
33 #include <vtkLongArray.h>
34 #include <vtkDoubleArray.h>
35 #include <vtkXMLUnstructuredGridWriter.h>
46 QVector
<QVector
<int> > &c2c
,
50 vtkUnstructuredGrid
*grid
55 QVector
<QSet
<int> > &n2n
,
60 protected: // attributes
62 QSet
<int> boundary_codes
;
67 * Update the cell index array.
69 void UpdateCellIndex(vtkUnstructuredGrid
*grid
);
72 * Update the point index array.
74 void UpdateNodeIndex(vtkUnstructuredGrid
*grid
);
77 * Compute normal vectors on nodes and cells of a subset of a grid.
78 * The paramters nodes and cells must be consistent; this means the nodes
79 * represent exactly (not more, not less) the nodes forming the cells.
80 * @param cell_normals On return, this will contain the cell normals (same order as cells)
81 * @param node_normals On return, this will contain the cell normals (same order as cells)
82 * @param cells The cells to compute the normals of
83 * @param nodes The nodes to compute the normals of
84 * @param grid The grid to operate on
88 QVector
<vec3_t
> &cell_normals
,
89 QVector
<vec3_t
> &node_normals
,
90 QVector
<vtkIdType
> &cells
,
91 QVector
<vtkIdType
> &nodes
,
92 vtkUnstructuredGrid
*grid
96 * Create a mapping from global node indices to the indeces of a subset of nodes.
97 * @param nodes The subset of nodes.
98 * @param _nodes On return, this will contain the mapping.
99 * @param grid The grid to operate on.
101 void createNodeMapping
103 QVector
<vtkIdType
> &nodes
,
104 QVector
<int> &_nodes
,
105 vtkUnstructuredGrid
*grid
109 * Create a mapping from global cell indices to the indeces of a subset of cells.
110 * @param cells The subset of cells.
111 * @param _cells On return, this will contain the mapping.
112 * @param grid The grid to operate on.
114 void createCellMapping
116 QVector
<vtkIdType
> &cells
,
117 QVector
<int> &_cells
,
118 vtkUnstructuredGrid
*grid
122 * Create a node to boundary condition ("cell_code") mapping.
123 * Onlu non-zero boundary conditions will be considered.
124 * @param bcs On return, this will hold the codes of all boundary elements that are
125 * attached to a node.
126 * @param grid The grid to operate on.
128 void createNodeToBcMapping
130 QVector
<QSet
<int> > &bcs
,
131 vtkUnstructuredGrid
*grid
135 * Create a node to cell structure for a given set of cells and nodes.
136 * @param cells the subset of cells
137 * @param nodes the subset of nodes
138 * @param _nodes the reverse mapping for the nodes
139 * @param n2c On return, this will hold the node to cell structure
140 * @param grid The grid to operate on
142 void createNodeToCell
144 QVector
<vtkIdType
> &cells
,
145 QVector
<vtkIdType
> &nodes
,
146 QVector
<int> &_nodes
,
147 QVector
<QSet
<int> > &n2c
,
148 vtkUnstructuredGrid
*grid
152 * Create a node to node structure for a given set of cells and nodes.
153 * @param cells the subset of cells
154 * @param nodes the subset of nodes
155 * @param _nodes the reverse mapping for the nodes
156 * @param n2c On return, this will hold the node to node structure
157 * @param grid The grid to operate on
159 void createNodeToNode
161 QVector
<vtkIdType
> &cells
,
162 QVector
<vtkIdType
> &nodes
,
163 QVector
<int> &_nodes
,
164 QVector
<QSet
<int> > &n2n
,
165 vtkUnstructuredGrid
*grid
169 * Extract the nodes which are part of a given set of cells.
170 * @param cells the subset of cells
171 * @param nodes On return, this will contain the nodes that correspond to the subset of cells
172 * @param grid The grid to operate on
174 void getNodesFromCells
176 QVector
<vtkIdType
> &cells
,
177 QVector
<vtkIdType
> &nodes
,
178 vtkUnstructuredGrid
*grid
182 * Check if a cell is a volume cell.
183 * @param cellId The id fof the cell in question
184 * @param grid The grid to operate on
185 * @return true is the cell represents a volume and false if not
189 vtkUnstructuredGrid
*grid
,
192 bool isVolume(vtkIdType id_cell
, vtkUnstructuredGrid
*grid
) { return isVolume(grid
, id_cell
); };
196 * Check if a cell is a surface cell.
197 * @param cellId The id fof the cell in question
198 * @param grid The grid to operate on
199 * @return true is the cell represents a surface and false if not
203 vtkUnstructuredGrid
*grid
,
206 bool isSurface(vtkIdType id_cell
, vtkUnstructuredGrid
*grid
) { return isSurface(grid
, id_cell
); };
209 * Get all volume cells of a grid.
210 * @param cells On return this will hold the Ids of all volume cells.
211 * @param grid The grid to operate on.
213 void getAllVolumeCells
215 QVector
<vtkIdType
> &cells
,
216 vtkUnstructuredGrid
*grid
220 * Get all cells of a grid.
221 * @param cells On return this will hold the Ids of all cells.
222 * @param grid The grid to operate on.
226 QVector
<vtkIdType
> &cells
,
227 vtkUnstructuredGrid
*grid
231 * Get all cells of a grid and a specific type.
232 * @param type The type of the cells (e.g. VTK_TETRA, VTK_TRIANGLE, etc.)
233 * @param cells On return this will hold the Ids of all cells.
234 * @param grid The grid to operate on.
236 void getAllCellsOfType
239 QVector
<vtkIdType
> &cells
,
240 vtkUnstructuredGrid
*grid
244 * Get all surface cells of a grid.
245 * @param cells On return this will hold the Ids of all surface cells.
246 * @param grid The grid to operate on.
248 void getAllSurfaceCells
250 QVector
<vtkIdType
> &cells
,
251 vtkUnstructuredGrid
*grid
255 * Get all surface cells of a grid with a specific boundary condition.
256 * @param bcs The set of boundary conditions
257 * @param cells On return this will hold the Ids of the surface cells.
258 * @param grid The grid to operate on.
263 QVector
<vtkIdType
> &cells
,
264 vtkUnstructuredGrid
*grid
268 * Create a cell neighbourship list for a subset grid.
269 * This has been implemented using VTK's vtkCellLinks structures.
270 * @param cells the subset of cells
271 * @param c2c On return this will hold the neighbourship list
272 * @param grid The grid to operate on.
274 void createCellToCell
276 QVector
<vtkIdType
> &cells
,
277 QVector
<QVector
<int> > &c2c
,
278 vtkUnstructuredGrid
*grid
282 * Insert a subset of subset of a grid into a vtkPolyData structure.
283 * This is can be used in order to make use of many readily available
284 * operations within VTK; one example is smoothing.
285 * Cell index and node index arrays will be created and passed to the
286 * poly data structure. Thus any purely geometric effect (no topology change)
287 * can be directly reintroduced into the vtkUnstructuredGrid.
288 * @param cells the subset of cells
289 * @param pdata the vtkPolyData to add the nodes and cells to
290 * @param grid The grid to operate on.
294 QVector
<vtkIdType
> &cells
,
296 vtkUnstructuredGrid
*grid
300 * Copy the attributes from an input to an output cell.
301 * @param old_grid the input grid
302 * @param oldId the existing input cell
303 * @param new_grid the output grid
304 * @param newId the new output cell
308 vtkUnstructuredGrid
*old_grid
,
310 vtkUnstructuredGrid
*new_grid
,
315 * Copy the attributes from an input to an output node.
316 * @param old_grid the input grid
317 * @param oldId the existing input node
318 * @param new_grid the output grid
319 * @param newId the new output node
323 vtkUnstructuredGrid
*old_grid
,
325 vtkUnstructuredGrid
*new_grid
,
330 * Create the basic fields on a given grid.
331 * Care should be taken with the overwrite parameter; if it is set to <i>false</i>
332 * and the field does not have the correct size it can lead to <i>segmentation faults</i>.
333 * @param Ncells the number of output cells
334 * @param Nnodes the number of output nodes
335 * @param overwrite f set to true existing fields will be re-created
337 void createBasicFields
339 vtkUnstructuredGrid
*grid
,
342 bool overwrite
= true
346 * Create the basic cell fields on a given grid.
347 * Care should be taken with the overwrite parameter; if it is set to <i>false</i>
348 * and the field does not have the correct size it can lead to <i>segmentation faults</i>.
349 * @param Ncells the number of output cells
350 * @param overwrite f set to true existing fields will be re-created
352 void createBasicCellFields
354 vtkUnstructuredGrid
*grid
,
356 bool overwrite
= true
360 * Create the basic node fields on a given grid.
361 * Care should be taken with the overwrite parameter; if it is set to <i>false</i>
362 * and the field does not have the correct size it can lead to <i>segmentation faults</i>.
363 * @param Nnodes the number of output nodes
364 * @param overwrite f set to true existing fields will be re-created
366 void createBasicNodeFields
368 vtkUnstructuredGrid
*grid
,
370 bool overwrite
= true
374 * Allocate memory for a grid. This method will also create the basic
375 * attribute fields (e.g. "cell_code").
376 * @param grid the grid for which to allocate memory
377 * @param Ncells the number of output cells
378 * @param Nnodes the number of output nodes
382 vtkUnstructuredGrid
*grid
,
388 * Compute the intersection of two QSets.
389 * @param set1 the first set
390 * @param set2 the second set
391 * @param inters On return this will hold the intersection
394 void setIntersection(const QSet
<T
> &set1
,
399 * Compute the intersection of two QVectors.
400 * @param set1 the first vector
401 * @param set2 the second vector
402 * @param inters On return this will hold the intersection
405 void vectorIntersection(const QVector
<T
> &set1
,
406 const QVector
<T
> &set2
,
410 * Compute the centre of a cell
411 * @param grid the grid to use
412 * @param id_cell the cell of which to compute the centre
413 * @return the centre of the cell
415 vec3_t
cellCentre(vtkUnstructuredGrid
*grid
, vtkIdType id_cell
);
418 * Get the cells of a grid that are not part of a given set of cells.
419 * @param grid the grid to use
420 * @param cells the given set of cells
421 * @param rest_cells on return this will hold the rest of all cells of the grid (not part of cells)
423 void getRestCells(vtkUnstructuredGrid
*grid
,
424 const QVector
<vtkIdType
> &cells
,
425 QVector
<vtkIdType
> &rest_cells
);
428 * Find the corresponding volume cell of a surface cell
429 * @param grid the grid to use
430 * @param id_surf the id of the surface cell
431 * @param n2n the node to cell structure for this grid
432 * @return the id of the corresponding volume cell (or -1 if not found)
436 vtkUnstructuredGrid
*grid
,
438 const QVector
<int> _nodes
,
439 const QVector
<vtkIdType
> cells
,
440 const QVector
<int> _cells
,
441 QVector
<QSet
<int> > &n2c
444 void makeCopy(vtkUnstructuredGrid
*src
, vtkUnstructuredGrid
*dst
);
445 void createIndices(vtkUnstructuredGrid
*grid
);
448 void writeCells(vtkUnstructuredGrid
*grid
, const T
&cls
, QString file_name
);
452 void setBoundaryCodes(const QSet
<int> &bcs
);
458 void EgVtkObject::setIntersection(const QSet
<T
> &set1
,
463 foreach (T t1
, set1
) {
464 if (set2
.contains(t1
)) {
471 void EgVtkObject::vectorIntersection(const QVector
<T
> &set1
,
472 const QVector
<T
> &set2
,
476 foreach (T t1
, set1
) {
484 void EgVtkObject::writeCells(vtkUnstructuredGrid
*grid
, const T
&cls
, QString file_name
)
487 EG_VTKSP(vtkUnstructuredGrid
,tmp_grid
);
488 QVector
<vtkIdType
> cells
;
489 QVector
<vtkIdType
> nodes
;
490 cells
.resize(cls
.size());
491 qCopy(cls
.begin(), cls
.end(), cells
.begin());
492 getNodesFromCells(cells
, nodes
, grid
);
493 allocateGrid(tmp_grid
, cells
.size(), nodes
.size());
494 vtkIdType id_new_node
= 0;
495 QVector
<vtkIdType
> old2new(grid
->GetNumberOfPoints(), -1);
496 foreach (vtkIdType id_node
, nodes
) {
498 grid
->GetPoint(id_node
, x
.data());
499 tmp_grid
->GetPoints()->SetPoint(id_new_node
, x
.data());
500 old2new
[id_node
] = id_new_node
;
501 copyNodeData(grid
, id_node
, tmp_grid
, id_new_node
);
504 foreach (vtkIdType id_cell
, cells
) {
505 vtkIdType N_pts
, *pts
;
506 vtkIdType type_cell
= grid
->GetCellType(id_cell
);
507 grid
->GetCellPoints(id_cell
, N_pts
, pts
);
508 QVector
<vtkIdType
> new_pts(N_pts
);
509 for (int i_pts
= 0; i_pts
< N_pts
; ++i_pts
) {
510 new_pts
[i_pts
] = old2new
[pts
[i_pts
]];
512 vtkIdType id_new_cell
= tmp_grid
->InsertNextCell(type_cell
, N_pts
, new_pts
.data());
513 copyCellData(grid
, id_cell
, tmp_grid
, id_new_cell
);
515 EG_VTKSP(vtkXMLUnstructuredGridWriter
,vtu
);
516 vtu
->SetFileName(file_name
.toAscii().data());
517 vtu
->SetDataModeToBinary();
518 vtu
->SetInput(tmp_grid
);