maybe now...
[engrid.git] / src / egvtkobject.h
bloba61601cac9c47bca22d4a8f35ea929469a2d3f97
1 //
2 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 // + +
4 // + This file is part of enGrid. +
5 // + +
6 // + Copyright 2008-2010 enGits GmbH +
7 // + +
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. +
12 // + +
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. +
17 // + +
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/>. +
20 // + +
21 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
23 #ifndef egvtkobject_H
24 #define egvtkobject_H
26 class EgVtkObject;
28 class BezierTriangle;
30 #include "engrid.h"
31 #include "utilities.h"
32 #include "boundarycondition.h"
34 #include <vtkUnstructuredGrid.h>
35 #include <vtkPolyData.h>
36 #include <vtkPointData.h>
37 #include <vtkCellData.h>
38 #include <vtkLongArray.h>
39 #include <vtkDoubleArray.h>
40 #include <vtkXMLUnstructuredGridWriter.h>
41 #include <vtkCellLocator.h>
43 #include <QSettings>
44 #include <QSet>
45 #include <QVector>
47 #define VTK_SIMPLE_VERTEX 0
48 #define VTK_FIXED_VERTEX 1
49 #define VTK_FEATURE_EDGE_VERTEX 2
50 #define VTK_BOUNDARY_EDGE_VERTEX 3
52 class EgVtkObject
55 public: // data-types
57 typedef const QVector<vtkIdType>& l2g_t;
58 typedef const QVector<int>& g2l_t;
59 typedef const QVector<QVector<int> >& l2l_t;
61 private: // methods
63 void addToC2C
65 vtkIdType id_cell,
66 QVector<int> &_cells,
67 QVector<QVector<int> > &c2c,
68 int j,
69 vtkIdList *nds,
70 vtkIdList *cls,
71 vtkUnstructuredGrid *grid
74 void addToN2N
76 QVector<QSet<int> > &n2n,
77 int n1,
78 int n2
81 void createNodeField(vtkUnstructuredGrid *grid, QString field_name, QString type_name, int Nnodes, bool overwrite = false);
82 void createCellField(vtkUnstructuredGrid *grid, QString field_name, QString type_name, int Ncells, bool overwrite = false);
84 protected: // attributes
86 QSet<int> m_BoundaryCodes;
87 static int DebugLevel;
89 protected: // methods
91 /**
92 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
93 * Version for int variables
95 int getSet(QString group, QString key, int value, int& variable);
97 /**
98 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
99 * Version for double variables
101 double getSet(QString group, QString key, double value, double& variable);
104 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
105 * Version for bool variables
107 bool getSet(QString group, QString key, bool value, bool& variable);
110 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
111 * Version for string variables
113 QString getSet(QString group, QString key, QString value, QString& variable);
116 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
117 * Version for string variables.
119 QString getSet(QString group, QString key, QString value, QString& variable, int type);
122 * Update the cell index array.
124 void UpdateCellIndex(vtkUnstructuredGrid *grid);
127 * Update the point index array.
129 void UpdateNodeIndex(vtkUnstructuredGrid *grid);
132 * Compute normal vectors on nodes and cells of a subset of a grid.
133 * The paramters nodes and cells must be consistent; this means the nodes
134 * represent exactly (not more, not less) the nodes forming the cells.
135 * @param cell_normals On return, this will contain the cell normals (same order as cells)
136 * @param node_normals On return, this will contain the cell normals (same order as cells)
137 * @param cells The cells to compute the normals of
138 * @param nodes The nodes to compute the normals of
139 * @param grid The grid to operate on
141 void computeNormals(QVector<vec3_t> &cell_normals, QVector<vec3_t> &node_normals, QVector<vtkIdType> &cells, QVector<vtkIdType> &nodes, vtkUnstructuredGrid *grid);
144 * Create a mapping from global node indices to the indeces of a subset of nodes.
145 * @param nodes The subset of nodes.
146 * @param _nodes On return, this will contain the mapping.
147 * @param grid The grid to operate on.
149 void createNodeMapping(QVector<vtkIdType> &nodes, QVector<int> &_nodes, vtkUnstructuredGrid *grid);
152 * Create a mapping from global cell indices to the indices of a subset of cells.
153 * @param cells The subset of cells.
154 * @param _cells On return, this will contain the mapping.
155 * @param grid The grid to operate on.
157 void createCellMapping(QVector<vtkIdType> &cells, QVector<int> &_cells, vtkUnstructuredGrid *grid);
160 * Create a node to boundary condition ("cell_code") mapping.
161 * Only non-zero boundary conditions will be considered.
162 * @param bcs On return, this will hold the codes of all boundary elements that are
163 * attached to a node.
164 * @param grid The grid to operate on.
166 void createNodeToBcMapping(QVector<QSet<int> > &bcs, vtkUnstructuredGrid *grid);
169 * Create a node to cell structure for a given set of cells and nodes.
170 * This creates a vector of sets which might have performance issues.
171 * @param cells the subset of cells
172 * @param nodes the subset of nodes
173 * @param _nodes the reverse mapping for the nodes
174 * @param n2c On return, this will hold the node to cell structure
175 * @param grid The grid to operate on
177 void createNodeToCell(QVector<vtkIdType> &cells, QVector<vtkIdType> &nodes, QVector<int> &_nodes, QVector<QSet<int> > &n2c, vtkUnstructuredGrid *grid);
180 * Create a node to cell structure for a given set of cells and nodes.
181 * This creates a vector of vectors.
182 * @param cells the subset of cells
183 * @param nodes the subset of nodes
184 * @param _nodes the reverse mapping for the nodes
185 * @param n2c On return, this will hold the node to cell structure
186 * @param grid The grid to operate on
188 void createNodeToCell(QVector<vtkIdType> &cells, QVector<vtkIdType> &nodes, QVector<int> &_nodes, QVector<QVector<int> > &n2c, vtkUnstructuredGrid *grid);
191 * Create a node to node structure for a given set of cells and nodes.
192 * This creates a vector of sets which might have performance issues.
193 * @param cells the subset of cells
194 * @param nodes the subset of nodes
195 * @param _nodes the reverse mapping for the nodes
196 * @param n2n On return, this will hold the node to node structure
197 * @param grid The grid to operate on
199 void createNodeToNode(QVector<vtkIdType> &cells, QVector<vtkIdType> &nodes, QVector<int> &_nodes, QVector<QSet<int> > &n2n, vtkUnstructuredGrid *grid);
202 * Create a node to node structure for a given set of cells and nodes.
203 * This creates a vector of vectors.
204 * @param cells the subset of cells
205 * @param nodes the subset of nodes
206 * @param _nodes the reverse mapping for the nodes
207 * @param n2n On return, this will hold the node to node structure
208 * @param grid The grid to operate on
210 void createNodeToNode(QVector<vtkIdType> &cells, QVector<vtkIdType> &nodes, QVector<int> &_nodes, QVector<QVector<int> > &n2n, vtkUnstructuredGrid *grid);
213 * Extract the nodes which are part of a given set of cells.
214 * @param cells the subset of cells
215 * @param nodes On return, this will contain the nodes that correspond to the subset of cells
216 * @param grid The grid to operate on
218 template <class C>
219 void getNodesFromCells(const C &cells, QVector<vtkIdType> &nodes, vtkUnstructuredGrid *grid);
222 * Check if a cell is a volume cell.
223 * @param cellId The id fof the cell in question
224 * @param grid The grid to operate on
225 * @return true if the cell represents a volume and false if not
227 bool isVolume(vtkIdType id_cell, vtkUnstructuredGrid *grid);
231 * Check if a cell is a surface cell.
232 * @param cellId The id fof the cell in question
233 * @param grid The grid to operate on
234 * @return true if the cell represents a surface and false if not
236 bool isSurface(vtkIdType id_cell, vtkUnstructuredGrid *grid);
239 * Get all volume cells of a grid.
240 * @param cells On return this will hold the Ids of all volume cells.
241 * @param grid The grid to operate on.
243 void getAllVolumeCells(QVector<vtkIdType> &cells, vtkUnstructuredGrid *grid);
246 * Get all cells of a grid.
247 * @param cells On return this will hold the Ids of all cells.
248 * @param grid The grid to operate on.
250 void getAllCells(QVector<vtkIdType> &cells, vtkUnstructuredGrid *grid);
253 * Get all cells of a grid and a specific type.
254 * @param type The type of the cells (e.g. VTK_TETRA, VTK_TRIANGLE, etc.)
255 * @param cells On return this will hold the Ids of all cells.
256 * @param grid The grid to operate on.
258 void getAllCellsOfType(vtkIdType type, QVector<vtkIdType> &cells, vtkUnstructuredGrid *grid);
261 * Get all surface cells of a grid.
262 * @param cells On return this will hold the Ids of all surface cells.
263 * @param grid The grid to operate on.
265 void getAllSurfaceCells(QVector<vtkIdType> &cells, vtkUnstructuredGrid *grid);
268 * Get all surface cells of a grid with a specific boundary condition.
269 * @param bcs The set of boundary conditions
270 * @param cells On return this will hold the Ids of the surface cells.
271 * @param grid The grid to operate on.
273 void getSurfaceCells(QSet<int> &bcs, QVector<vtkIdType> &cells, vtkUnstructuredGrid *grid);
276 * Create a cell neighbourship list for a subset grid.
277 * This has been implemented using VTK's vtkCellLinks structures.
278 * @param cells the subset of cells
279 * @param c2c On return this will hold the neighbourship list
280 * @param grid The grid to operate on.
282 void createCellToCell(QVector<vtkIdType> &cells, QVector<QVector<int> > &c2c, vtkUnstructuredGrid *grid);
285 * Insert a subset of a grid into a vtkPolyData structure.
286 * This is can be used in order to make use of many readily available
287 * operations within VTK; one example is smoothing.
288 * Cell index and node index arrays will be created and passed to the
289 * poly data structure. Thus any purely geometric effect (no topology change)
290 * can be directly reintroduced into the vtkUnstructuredGrid.
291 * @param cells the subset of cells
292 * @param pdata the vtkPolyData to add the nodes and cells to
293 * @param grid The grid to operate on.
295 void addToPolyData(QVector<vtkIdType> &cells, vtkPolyData *pdata, vtkUnstructuredGrid *grid);
298 * Copy the attributes from an input to an output cell.
299 * @param old_grid the input grid
300 * @param oldId the existing input cell
301 * @param new_grid the output grid
302 * @param newId the new output cell
304 void copyCellData(vtkUnstructuredGrid *old_grid, vtkIdType oldId, vtkUnstructuredGrid *new_grid, vtkIdType newId);
307 * Copy the attributes from an input to an output node.
308 * @param old_grid the input grid
309 * @param oldId the existing input node
310 * @param new_grid the output grid
311 * @param newId the new output node
313 void copyNodeData(vtkUnstructuredGrid *old_grid, vtkIdType oldId, vtkUnstructuredGrid *new_grid, vtkIdType newId);
316 * Create the basic fields on a given grid.
317 * Care should be taken with the overwrite parameter; if it is set to <i>false</i>
318 * and the field does not have the correct size it can lead to <i>segmentation faults</i>.
319 * @param Ncells the number of output cells
320 * @param Nnodes the number of output nodes
321 * @param overwrite f set to true existing fields will be re-created
323 void createBasicFields(vtkUnstructuredGrid *grid, vtkIdType Ncells, vtkIdType Nnodes, bool overwrite = false);
326 * Create the basic cell fields on a given grid.
327 * Care should be taken with the overwrite parameter; if it is set to <i>false</i>
328 * and the field does not have the correct size it can lead to <i>segmentation faults</i>.
329 * @param Ncells the number of output cells
330 * @param overwrite f set to true existing fields will be re-created
332 void createBasicCellFields(vtkUnstructuredGrid *grid, vtkIdType Ncells, bool overwrite = false);
335 * Create the basic node fields on a given grid.
336 * Care should be taken with the overwrite parameter; if it is set to <i>false</i>
337 * and the field does not have the correct size it can lead to <i>segmentation faults</i>.
338 * @param Nnodes the number of output nodes
339 * @param overwrite f set to true existing fields will be re-created
341 void createBasicNodeFields(vtkUnstructuredGrid *grid, vtkIdType Nnodes, bool overwrite = false);
344 * Allocate memory for a grid. This method will also create the basic
345 * attribute fields (e.g. "cell_code").
346 * @param grid the grid for which to allocate memory
347 * @param Ncells the number of output cells
348 * @param Nnodes the number of output nodes
349 * @param create_fields flag to determine if node and cell data shall be created
351 void allocateGrid(vtkUnstructuredGrid *grid, vtkIdType Ncells, vtkIdType Nnodes, bool create_fields = true);
354 * Get the names of all node (point) attributes (fields) of a VTK grid
355 * @param field_names On return this vector will contain the names of all fields
356 * @param grid the grid
358 void getAllNodeDataNames(QVector<QString> &field_names, vtkUnstructuredGrid *grid);
361 * Get the names of all cell attributes (fields) of a VTK grid
362 * @param field_names On return this vector will contain the names of all fields
363 * @param grid the grid
365 void getAllCellDataNames(QVector<QString> &field_names, vtkUnstructuredGrid *grid);
368 * Compute the intersection of two Q containers.
369 * This will return a set.
370 * @param set1 the first container
371 * @param set2 the second container
372 * @param inters on return this will hold the intersection
374 template <class C1, class C2>
375 void qcontIntersection(const C1& c1, const C2& c2, QSet<typename C1::value_type> &inters);
378 * Compute the intersection of two Q containers.
379 * This will return a vector.
380 * @param set1 the first container
381 * @param set2 the second container
382 * @param inters on return this will hold the intersection
384 template <class C1, class C2>
385 void qcontIntersection(const C1& c1, const C2& c2, QVector<typename C1::value_type> &inters);
388 * Compute the centre of a cell
389 * @param grid the grid to use
390 * @param id_cell the cell of which to compute the centre
391 * @return the centre of the cell
393 vec3_t cellCentre(vtkUnstructuredGrid *grid, vtkIdType id_cell);
396 * Get the cells of a grid that are not part of a given set of cells.
397 * @param grid the grid to use
398 * @param cells the given set of cells
399 * @param rest_cells on return this will hold the rest of all cells of the grid (not part of cells)
401 void getRestCells(vtkUnstructuredGrid *grid, const QVector<vtkIdType> &cells, QVector<vtkIdType> &rest_cells);
404 * Find the corresponding volume cell of a surface cell
405 * @param grid the grid to use
406 * @param id_surf the id of the surface cell
407 * @param n2n the node to cell structure for this grid
408 * @return the id of the corresponding volume cell (or -1 if not found)
410 vtkIdType findVolumeCell(vtkUnstructuredGrid *grid, vtkIdType id_surf, g2l_t _nodes, l2g_t cells, g2l_t _cells, l2l_t n2c);
413 * Copy "src" grid to "dst" grid. Allocate "dst" so that it fits the data of "src".
414 * @param src a pointer to the source grid
415 * @param dst a pointer to the destination grid
417 void makeCopy(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst);
420 * Copy a part of "src" grid to "dst" grid. Allocate "dst" so that it fits the data to be copied.
421 * @param src a pointer to the source grid
422 * @param dst a pointer to the destination grid
423 * @param cells a container with the cells to be copied
425 template <class C>
426 void makeCopy(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst, const C &cells);
429 * Copy "src" grid to "dst" grid. DO NOT allocate "dst" so that it fits the data of "src".
430 * Allocation is left for the user to do.
431 * @param src a pointer to the source grid
432 * @param dst a pointer to the destination grid
434 void makeCopyNoAlloc(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst);
437 * Change the orientation of a face.
438 * @param grid the grid to use
439 * @param id_face the id of the face to change
441 void reorientateFace(vtkUnstructuredGrid *grid, vtkIdType id_face);
444 * Reset face orientation to original orientation.
445 * @param grid the grid with the faces
447 void resetOrientation(vtkUnstructuredGrid *grid);
449 void createIndices(vtkUnstructuredGrid *grid);
452 * Get the boundary condition of a boundary code.
453 * @param bc the boundary code
454 * @return the boundary condition
456 BoundaryCondition getBC(int bc);
459 * Save the subgrid defined by cls from grid.
460 * @param grid The source grid
461 * @param cls The cells to extract
462 * @param file_name The file to save to
464 template <class C>
465 void writeCells(vtkUnstructuredGrid *grid, const C &cls, QString file_name);
468 * Get the SubGrid defined by cls from grid. The function takes care of allocation for SubGrid.
469 * @param grid The source grid
470 * @param cls The cells to extract
471 * @param SubGrid The SubGrid to create
473 template <class C>
474 void getSubGrid(vtkUnstructuredGrid *grid, const C &cls, vtkUnstructuredGrid *SubGrid);
478 * Save grid to file filename.
479 * @param grid The source grid
480 * @param filename Name of the file to save to
482 void writeGrid(vtkUnstructuredGrid *grid, QString filename);
485 * Get a file name without extension.
486 * @param file_name the full name (with extension)
487 * @return the name without the extension
489 QString stripFromExtension(QString file_name);
492 * Get the extension of a file name
493 * @param file_name the full name (with extension)
494 * @return the extension
496 QString getExtension(QString file_name);
499 * Get a face of a cell
500 * @param grid the unstructured grid
501 * @param id_cell the index of the cell
502 * @param i_face the index of the face within the cell
503 * @param ids on return this will contain the nodes of the face
505 void getFaceOfCell(vtkUnstructuredGrid *grid, vtkIdType id_cell, int i_face, QVector<vtkIdType> &ids);
508 * Get an edge of a face/cell
509 * @param grid the unstructured grid
510 * @param id_cell the index of the cell
511 * @param i_edge the index of the edge within the cell
512 * @param ids on return this will contain the nodes of the edge
514 void getEdgeOfCell(vtkUnstructuredGrid *grid, vtkIdType id_cell, int i_edge, QVector<vtkIdType> &ids);
516 public: // methods
518 EgVtkObject() { DebugLevel = 0; }
520 void setBoundaryCodes(const QSet<int> &bcs);
521 QSet<int> getBoundaryCodes();
522 void setDebugLevel(int a_DebugLevel) { DebugLevel = a_DebugLevel; }
524 bool saveGrid( vtkUnstructuredGrid* a_grid, QString file_name );
526 vtkIdType addGrid(vtkUnstructuredGrid *main_grid, vtkUnstructuredGrid *grid_to_add, vtkIdType offset);
528 private:
529 void addVtkTypeInfo(vtkUnstructuredGrid* a_grid); ///< Add VTK type information to the grid (useful for visualisation with ParaView).
532 //End of class EgVtkObject
534 template <class C1, class C2>
535 void EgVtkObject::qcontIntersection(const C1& c1, const C2& c2, QSet<typename C1::value_type> &inters)
537 inters.clear();
538 foreach (typename C1::value_type t1, c1) {
539 foreach (typename C2::value_type t2, c2) {
540 if (t1 == t2) {
541 inters.insert(t1);
547 template <class C1, class C2>
548 void EgVtkObject::qcontIntersection(const C1& c1, const C2& c2, QVector<typename C1::value_type> &inters)
550 QSet<typename C1::value_type> inters_set;
551 qcontIntersection(c1, c2, inters_set);
552 inters.resize(inters_set.size());
553 qCopy(inters_set.begin(), inters_set.end(), inters.begin());
556 template <class C>
557 void EgVtkObject::getSubGrid(vtkUnstructuredGrid *grid, const C &cls, vtkUnstructuredGrid *SubGrid)
559 createIndices(grid);
560 QVector<vtkIdType> cells;
561 QVector<vtkIdType> nodes;
562 cells.resize(cls.size());
563 qCopy(cls.begin(), cls.end(), cells.begin());
564 getNodesFromCells(cells, nodes, grid);
565 allocateGrid(SubGrid, cells.size(), nodes.size());
566 vtkIdType id_new_node = 0;
567 QVector<vtkIdType> old2new(grid->GetNumberOfPoints(), -1);
568 foreach (vtkIdType id_node, nodes) {
569 vec3_t x;
570 grid->GetPoint(id_node, x.data());
571 SubGrid->GetPoints()->SetPoint(id_new_node, x.data());
572 old2new[id_node] = id_new_node;
573 copyNodeData(grid, id_node, SubGrid, id_new_node);
574 ++id_new_node;
576 foreach (vtkIdType id_cell, cells) {
577 vtkIdType N_pts, *pts;
578 vtkIdType type_cell = grid->GetCellType(id_cell);
579 grid->GetCellPoints(id_cell, N_pts, pts);
580 QVector<vtkIdType> new_pts(N_pts);
581 for (int i_pts = 0; i_pts < N_pts; ++i_pts) {
582 new_pts[i_pts] = old2new[pts[i_pts]];
584 vtkIdType id_new_cell = SubGrid->InsertNextCell(type_cell, N_pts, new_pts.data());
585 copyCellData(grid, id_cell, SubGrid, id_new_cell);
589 template <class C>
590 void EgVtkObject::writeCells(vtkUnstructuredGrid *grid, const C &cls, QString file_name)
592 qDebug()<<"Saving cells from grid as "<<file_name;
594 EG_VTKSP(vtkUnstructuredGrid,SubGrid);
595 getSubGrid(grid,cls,SubGrid);
597 EG_VTKSP(vtkXMLUnstructuredGridWriter,vtu);
598 vtu->SetFileName(qPrintable(file_name));
599 vtu->SetDataModeToBinary();
600 vtu->SetInput(SubGrid);
601 vtu->Write();
604 template <class C>
605 void EgVtkObject::getNodesFromCells(const C& cells, QVector<vtkIdType> &nodes, vtkUnstructuredGrid *grid)
607 QSet<vtkIdType> ex_nodes;
608 vtkIdType id_cell;
609 foreach(id_cell, cells) {
610 vtkIdType *pts;
611 vtkIdType Npts;
612 grid->GetCellPoints(id_cell, Npts, pts);
613 for (int i = 0; i < Npts; ++i) {
614 ex_nodes.insert(pts[i]);
617 nodes.resize(ex_nodes.size());
619 int j = 0;
620 vtkIdType i;
621 foreach(i,ex_nodes) {
622 nodes[j] = i;
623 ++j;
628 template <class C>
629 void EgVtkObject::makeCopy(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst, const C& cells)
631 QVector<vtkIdType> nodes;
632 getNodesFromCells(cells, nodes, src);
633 allocateGrid(dst, cells.size(), nodes.size());
634 vtkIdType id_new_node = 0;
635 foreach (vtkIdType id_node, nodes) {
636 vec3_t x;
637 src->GetPoints()->GetPoint(id_node, x.data());
638 dst->GetPoints()->SetPoint(id_node, x.data());
639 copyNodeData(src, id_node, dst, id_new_node);
640 ++id_new_node;
642 foreach (vtkIdType id_cell, cells) {
643 vtkIdType N_pts, *pts;
644 vtkIdType type_cell = src->GetCellType(id_cell);
645 src->GetCellPoints(id_cell, N_pts, pts);
646 vtkIdType id_new_cell = dst->InsertNextCell(type_cell, N_pts, pts);
647 copyCellData(src, id_cell, dst, id_new_cell);
651 #endif