some gui changes
[engrid.git] / egvtkobject.h
blob5a446b81da2abf17be401a6bf0c912916b25d1d1
1 //
2 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 // + +
4 // + This file is part of enGrid. +
5 // + +
6 // + Copyright 2008,2009 Oliver Gloth +
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 #include "engrid.h"
29 #include "boundarycondition.h"
31 #include <vtkUnstructuredGrid.h>
32 #include <vtkPolyData.h>
33 #include <vtkPointData.h>
34 #include <vtkCellData.h>
35 #include <vtkLongArray.h>
36 #include <vtkDoubleArray.h>
37 #include <vtkXMLUnstructuredGridWriter.h>
38 #include <QSettings>
39 #include <QSet>
40 #include <QVector>
42 #define VTK_SIMPLE_VERTEX 0
43 #define VTK_FIXED_VERTEX 1
44 #define VTK_FEATURE_EDGE_VERTEX 2
45 #define VTK_BOUNDARY_EDGE_VERTEX 3
47 class EgVtkObject
50 private: // methods
52 void addToC2C
54 vtkIdType id_cell,
55 QVector<int> &_cells,
56 QVector<QVector<int> > &c2c,
57 int j,
58 vtkIdList *nds,
59 vtkIdList *cls,
60 vtkUnstructuredGrid *grid
63 void addToN2N
65 QVector<QSet<int> > &n2n,
66 int n1,
67 int n2
70 protected: // attributes
72 QSet<int> boundary_codes;
73 int DebugLevel;
75 protected: // methods
77 /**
78 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
79 * Version for int variables
81 int getSet(QString group, QString key, int value, int& variable);
83 /**
84 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
85 * Version for double variables
87 double getSet(QString group, QString key, double value, double& variable);
89 /**
90 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
91 * Version for bool variables
93 bool getSet(QString group, QString key, bool value, bool& variable);
95 /**
96 * Update the cell index array.
98 void UpdateCellIndex(vtkUnstructuredGrid *grid);
101 * Update the point index array.
103 void UpdateNodeIndex(vtkUnstructuredGrid *grid);
106 * Compute normal vectors on nodes and cells of a subset of a grid.
107 * The paramters nodes and cells must be consistent; this means the nodes
108 * represent exactly (not more, not less) the nodes forming the cells.
109 * @param cell_normals On return, this will contain the cell normals (same order as cells)
110 * @param node_normals On return, this will contain the cell normals (same order as cells)
111 * @param cells The cells to compute the normals of
112 * @param nodes The nodes to compute the normals of
113 * @param grid The grid to operate on
115 void computeNormals
117 QVector<vec3_t> &cell_normals,
118 QVector<vec3_t> &node_normals,
119 QVector<vtkIdType> &cells,
120 QVector<vtkIdType> &nodes,
121 vtkUnstructuredGrid *grid
125 * Create a mapping from global node indices to the indeces of a subset of nodes.
126 * @param nodes The subset of nodes.
127 * @param _nodes On return, this will contain the mapping.
128 * @param grid The grid to operate on.
130 void createNodeMapping
132 QVector<vtkIdType> &nodes,
133 QVector<int> &_nodes,
134 vtkUnstructuredGrid *grid
138 * Create a mapping from global cell indices to the indeces of a subset of cells.
139 * @param cells The subset of cells.
140 * @param _cells On return, this will contain the mapping.
141 * @param grid The grid to operate on.
143 void createCellMapping
145 QVector<vtkIdType> &cells,
146 QVector<int> &_cells,
147 vtkUnstructuredGrid *grid
151 * Create a node to boundary condition ("cell_code") mapping.
152 * Only non-zero boundary conditions will be considered.
153 * @param bcs On return, this will hold the codes of all boundary elements that are
154 * attached to a node.
155 * @param grid The grid to operate on.
157 void createNodeToBcMapping
159 QVector<QSet<int> > &bcs,
160 vtkUnstructuredGrid *grid
164 * Create a node to cell structure for a given set of cells and nodes.
165 * @param cells the subset of cells
166 * @param nodes the subset of nodes
167 * @param _nodes the reverse mapping for the nodes
168 * @param n2c On return, this will hold the node to cell structure
169 * @param grid The grid to operate on
171 void createNodeToCell
173 QVector<vtkIdType> &cells,
174 QVector<vtkIdType> &nodes,
175 QVector<int> &_nodes,
176 QVector<QSet<int> > &n2c,
177 vtkUnstructuredGrid *grid
181 * Create a node to node structure for a given set of cells and nodes.
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 node structure
186 * @param grid The grid to operate on
188 void createNodeToNode
190 QVector<vtkIdType> &cells,
191 QVector<vtkIdType> &nodes,
192 QVector<int> &_nodes,
193 QVector<QSet<int> > &n2n,
194 vtkUnstructuredGrid *grid
198 * Extract the nodes which are part of a given set of cells.
199 * @param cells the subset of cells
200 * @param nodes On return, this will contain the nodes that correspond to the subset of cells
201 * @param grid The grid to operate on
203 void getNodesFromCells
205 QVector<vtkIdType> &cells,
206 QVector<vtkIdType> &nodes,
207 vtkUnstructuredGrid *grid
211 * Check if a cell is a volume cell.
212 * @param cellId The id fof the cell in question
213 * @param grid The grid to operate on
214 * @return true if the cell represents a volume and false if not
216 bool isVolume
218 vtkUnstructuredGrid *grid,
219 vtkIdType cellId
221 bool isVolume(vtkIdType id_cell, vtkUnstructuredGrid *grid) { return isVolume(grid, id_cell); };
225 * Check if a cell is a surface cell.
226 * @param cellId The id fof the cell in question
227 * @param grid The grid to operate on
228 * @return true if the cell represents a surface and false if not
230 bool isSurface
232 vtkUnstructuredGrid *grid,
233 vtkIdType id_cell
235 bool isSurface(vtkIdType id_cell, vtkUnstructuredGrid *grid) { return isSurface(grid, id_cell); };
238 * Get all volume cells of a grid.
239 * @param cells On return this will hold the Ids of all volume cells.
240 * @param grid The grid to operate on.
242 void getAllVolumeCells
244 QVector<vtkIdType> &cells,
245 vtkUnstructuredGrid *grid
249 * Get all cells of a grid.
250 * @param cells On return this will hold the Ids of all cells.
251 * @param grid The grid to operate on.
253 void getAllCells
255 QVector<vtkIdType> &cells,
256 vtkUnstructuredGrid *grid
260 * Get all cells of a grid and a specific type.
261 * @param type The type of the cells (e.g. VTK_TETRA, VTK_TRIANGLE, etc.)
262 * @param cells On return this will hold the Ids of all cells.
263 * @param grid The grid to operate on.
265 void getAllCellsOfType
267 vtkIdType type,
268 QVector<vtkIdType> &cells,
269 vtkUnstructuredGrid *grid
273 * Get all surface cells of a grid.
274 * @param cells On return this will hold the Ids of all surface cells.
275 * @param grid The grid to operate on.
277 void getAllSurfaceCells
279 QVector<vtkIdType> &cells,
280 vtkUnstructuredGrid *grid
284 * Get all surface cells of a grid with a specific boundary condition.
285 * @param bcs The set of boundary conditions
286 * @param cells On return this will hold the Ids of the surface cells.
287 * @param grid The grid to operate on.
289 void getSurfaceCells
291 QSet<int> &bcs,
292 QVector<vtkIdType> &cells,
293 vtkUnstructuredGrid *grid
297 * Get all surface nodes of a grid with a specific boundary condition.
298 * @param bcs The set of boundary conditions
299 * @param nodes On return this will hold the Ids of the surface nodes.
300 * @param grid The grid to operate on.
302 void getSurfaceNodes
304 QSet<int> &bcs,
305 QSet <vtkIdType> &SelectedNodes,
306 vtkUnstructuredGrid *grid
310 * Create a cell neighbourship list for a subset grid.
311 * This has been implemented using VTK's vtkCellLinks structures.
312 * @param cells the subset of cells
313 * @param c2c On return this will hold the neighbourship list
314 * @param grid The grid to operate on.
316 void createCellToCell
318 QVector<vtkIdType> &cells,
319 QVector<QVector<int> > &c2c,
320 vtkUnstructuredGrid *grid
324 * Insert a subset of a grid into a vtkPolyData structure.
325 * This is can be used in order to make use of many readily available
326 * operations within VTK; one example is smoothing.
327 * Cell index and node index arrays will be created and passed to the
328 * poly data structure. Thus any purely geometric effect (no topology change)
329 * can be directly reintroduced into the vtkUnstructuredGrid.
330 * @param cells the subset of cells
331 * @param pdata the vtkPolyData to add the nodes and cells to
332 * @param grid The grid to operate on.
334 void addToPolyData
336 QVector<vtkIdType> &cells,
337 vtkPolyData *pdata,
338 vtkUnstructuredGrid *grid
342 * Copy the attributes from an input to an output cell.
343 * @param old_grid the input grid
344 * @param oldId the existing input cell
345 * @param new_grid the output grid
346 * @param newId the new output cell
348 void copyCellData
350 vtkUnstructuredGrid *old_grid,
351 vtkIdType oldId,
352 vtkUnstructuredGrid *new_grid,
353 vtkIdType newId
357 * Copy the attributes from an input to an output node.
358 * @param old_grid the input grid
359 * @param oldId the existing input node
360 * @param new_grid the output grid
361 * @param newId the new output node
363 void copyNodeData
365 vtkUnstructuredGrid *old_grid,
366 vtkIdType oldId,
367 vtkUnstructuredGrid *new_grid,
368 vtkIdType newId
372 * Create the basic fields on a given grid.
373 * Care should be taken with the overwrite parameter; if it is set to <i>false</i>
374 * and the field does not have the correct size it can lead to <i>segmentation faults</i>.
375 * @param Ncells the number of output cells
376 * @param Nnodes the number of output nodes
377 * @param overwrite f set to true existing fields will be re-created
379 void createBasicFields
381 vtkUnstructuredGrid *grid,
382 vtkIdType Ncells,
383 vtkIdType Nnodes,
384 bool overwrite = true
388 * Create the basic cell fields on a given grid.
389 * Care should be taken with the overwrite parameter; if it is set to <i>false</i>
390 * and the field does not have the correct size it can lead to <i>segmentation faults</i>.
391 * @param Ncells the number of output cells
392 * @param overwrite f set to true existing fields will be re-created
394 void createBasicCellFields
396 vtkUnstructuredGrid *grid,
397 vtkIdType Ncells,
398 bool overwrite = true
402 * Create the basic node fields on a given grid.
403 * Care should be taken with the overwrite parameter; if it is set to <i>false</i>
404 * and the field does not have the correct size it can lead to <i>segmentation faults</i>.
405 * @param Nnodes the number of output nodes
406 * @param overwrite f set to true existing fields will be re-created
408 void createBasicNodeFields
410 vtkUnstructuredGrid *grid,
411 vtkIdType Nnodes,
412 bool overwrite = true
416 * Allocate memory for a grid. This method will also create the basic
417 * attribute fields (e.g. "cell_code").
418 * @param grid the grid for which to allocate memory
419 * @param Ncells the number of output cells
420 * @param Nnodes the number of output nodes
422 void allocateGrid
424 vtkUnstructuredGrid *grid,
425 vtkIdType Ncells,
426 vtkIdType Nnodes
430 * Compute the intersection of two QSets.
431 * @param set1 the first set
432 * @param set2 the second set
433 * @param inters On return this will hold the intersection
435 template <class T>
436 void setIntersection(const QSet<T> &set1,
437 const QSet<T> &set2,
438 QSet<T> &inters);
441 * Compute the intersection of two QVectors.
442 * @param set1 the first vector
443 * @param set2 the second vector
444 * @param inters On return this will hold the intersection
446 template <class T>
447 void vectorIntersection(const QVector<T> &set1,
448 const QVector<T> &set2,
449 QSet<T> &inters);
452 * Compute the centre of a cell
453 * @param grid the grid to use
454 * @param id_cell the cell of which to compute the centre
455 * @return the centre of the cell
457 vec3_t cellCentre(vtkUnstructuredGrid *grid, vtkIdType id_cell);
460 * Get the cells of a grid that are not part of a given set of cells.
461 * @param grid the grid to use
462 * @param cells the given set of cells
463 * @param rest_cells on return this will hold the rest of all cells of the grid (not part of cells)
465 void getRestCells(vtkUnstructuredGrid *grid,
466 const QVector<vtkIdType> &cells,
467 QVector<vtkIdType> &rest_cells);
470 * Find the corresponding volume cell of a surface cell
471 * @param grid the grid to use
472 * @param id_surf the id of the surface cell
473 * @param n2n the node to cell structure for this grid
474 * @return the id of the corresponding volume cell (or -1 if not found)
476 int findVolumeCell
478 vtkUnstructuredGrid *grid,
479 vtkIdType id_surf,
480 const QVector<int> _nodes,
481 const QVector<vtkIdType> cells,
482 const QVector<int> _cells,
483 QVector<QSet<int> > &n2c
487 * Copy "src" grid to "dst" grid. Allocate "dst" so that it fits the data of "src".
489 void makeCopy(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst);
491 * Copy "src" grid to "dst" grid. DO NOT allocate "dst" so that it fits the data of "src".
492 * Allocation is left for the user to do.
494 void makeCopyNoAlloc(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst);
497 * Copy "src" grid to "dst" grid. DO NOT allocate "dst" so that it fits the data of "src".
498 * Allocation is left for the user to do.
499 * Filter is a vector specifying whether a node should be removed or not.
500 * false: don't remove
501 * true: remove
503 void makeCopyNoAllocFiltered(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst, vector <bool> DeadNode);
504 // void makeCopyNoAllocFiltered(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst, vector <bool> DeadNode, QVector <QSet <vtkIdType>> newCells);
506 void createIndices(vtkUnstructuredGrid *grid);
509 * Get the boundary condition of a boundary code.
510 * @param bc the boundary code
511 * @return the boundary condition
513 BoundaryCondition getBC(int bc);
515 template <class T>
516 void writeCells(vtkUnstructuredGrid *grid, const T &cls, QString file_name);
518 public: // methods
520 void setBoundaryCodes(const QSet<int> &bcs);
521 EgVtkObject(){
522 DebugLevel=0;
526 //End of class EgVtkObject
528 template <class T>
529 void EgVtkObject::setIntersection(const QSet<T> &set1,
530 const QSet<T> &set2,
531 QSet<T> &inters)
533 inters.clear();
534 foreach (T t1, set1) {
535 if (set2.contains(t1)) {
536 inters.insert(t1);
541 template <class T>
542 void EgVtkObject::vectorIntersection(const QVector<T> &set1,
543 const QVector<T> &set2,
544 QSet<T> &inters)
546 inters.clear();
547 foreach (T t1, set1) {
548 if (set2.has(t1)) {
549 inters.insert(t1);
554 template <class T>
555 void EgVtkObject::writeCells(vtkUnstructuredGrid *grid, const T &cls, QString file_name)
557 createIndices(grid);
558 EG_VTKSP(vtkUnstructuredGrid,tmp_grid);
559 QVector<vtkIdType> cells;
560 QVector<vtkIdType> nodes;
561 cells.resize(cls.size());
562 qCopy(cls.begin(), cls.end(), cells.begin());
563 getNodesFromCells(cells, nodes, grid);
564 allocateGrid(tmp_grid, cells.size(), nodes.size());
565 vtkIdType id_new_node = 0;
566 QVector<vtkIdType> old2new(grid->GetNumberOfPoints(), -1);
567 foreach (vtkIdType id_node, nodes) {
568 vec3_t x;
569 grid->GetPoint(id_node, x.data());
570 tmp_grid->GetPoints()->SetPoint(id_new_node, x.data());
571 old2new[id_node] = id_new_node;
572 copyNodeData(grid, id_node, tmp_grid, id_new_node);
573 ++id_new_node;
575 foreach (vtkIdType id_cell, cells) {
576 vtkIdType N_pts, *pts;
577 vtkIdType type_cell = grid->GetCellType(id_cell);
578 grid->GetCellPoints(id_cell, N_pts, pts);
579 QVector<vtkIdType> new_pts(N_pts);
580 for (int i_pts = 0; i_pts < N_pts; ++i_pts) {
581 new_pts[i_pts] = old2new[pts[i_pts]];
583 vtkIdType id_new_cell = tmp_grid->InsertNextCell(type_cell, N_pts, new_pts.data());
584 copyCellData(grid, id_cell, tmp_grid, id_new_cell);
586 EG_VTKSP(vtkXMLUnstructuredGridWriter,vtu);
587 vtu->SetFileName(file_name.toAscii().data());
588 vtu->SetDataModeToBinary();
589 vtu->SetInput(tmp_grid);
590 vtu->Write();
594 * Utility function that allows printing selected data from an vtkUnstructuredGrid to any ostream (includes ofstream objects)
595 * @param stream ostream object to print to
596 * @param grid vtkUnstructuredGrid you want to print data from
597 * @param npoints print number of points in the grid
598 * @param ncells print number of cells in the grid
599 * @param points print points in the grid
600 * @param cells print cells in the grid
602 int cout_grid(ostream &stream, vtkUnstructuredGrid *grid, bool npoints=true, bool ncells=true, bool points=false, bool cells=false);
604 ///////////////////////////////////////////
605 int addPoint(vtkUnstructuredGrid* a_grid,vtkIdType index,vec3_t x);
606 int addCell(vtkUnstructuredGrid* a_grid, vtkIdType A, vtkIdType B, vtkIdType C, int bc);
608 int getShortestSide(vtkIdType a_id_cell,vtkUnstructuredGrid* a_grid);
609 int getLongestSide(vtkIdType a_id_cell,vtkUnstructuredGrid* a_grid);
611 int getSide(vtkIdType a_id_cell,vtkUnstructuredGrid* a_grid,vtkIdType a_id_node1,vtkIdType a_id_node2);
613 QSet <int> complementary_bcs(QSet <int> &bcs, vtkUnstructuredGrid *a_grid, QVector <vtkIdType> &a_cells);
614 QString cell2str(vtkIdType id_cell,vtkUnstructuredGrid* grid);
615 Qt::CheckState int2CheckState(int a);
616 int CheckState2int(Qt::CheckState a);
618 ///////////////////////////////////////////
619 template <class T>
620 ostream &operator<<(ostream &out, QVector<T> & vector)
622 int N=vector.size();
623 out<<"[";
624 for (int i = 0; i < N; ++i) {
625 out<<vector.at(i);
626 if(i!=N-1) out<<",";
628 out<<"]";
629 return(out);
632 template <class T>
633 ostream &operator<<(ostream &out, QSet<T> const & set )
635 out << "[ ";
636 foreach (T value, set) out << value << " ";
637 out << "]";
638 return(out);
641 template <class T>
642 ostream &operator<<(ostream &out, QVector<QSet<T> > & vector)
644 int N=vector.size();
645 out<<"[";
646 for (int i = 0; i < N; ++i) {
647 QSet<T> set=vector.at(i);
648 out<<set;
649 if(i!=N-1) out<<",";
651 out<<"]";
652 return(out);
655 template <class T>
656 ostream &operator<<(ostream &out, QVector<QVector<T> > & vector)
658 int N=vector.size();
659 out<<"[";
660 for (int i = 0; i < N; ++i) {
661 QVector<T> subvector=vector.at(i);
662 out<<subvector;
663 if(i!=N-1) out<<",";
665 out<<"]";
666 return(out);
669 template <class T>
670 ostream &operator<<(ostream &out, QMap<T,bool> & map)
672 QMapIterator<T, bool> i(map);
673 out<<"[";
674 while (i.hasNext()) {
675 i.next();
676 out << " [" << i.key() << ": " << i.value() << "]";
678 out<<"]";
679 return(out);
681 ///////////////////////////////////////////
682 // ///////////////////////////////////////////
683 // /* Here is how we we get QTextStreams that look like iostreams */
684 // QTextStream Qcin;
685 // QTextStream Qcout;
686 // QTextStream Qcerr;
687 // ///////////////////////////////////////////
689 pair<vtkIdType,vtkIdType> OrderedPair(vtkIdType a, vtkIdType b);
691 vtkIdType nextcell(vtkIdType a_cell, vtkIdType a_node, QVector< QVector< int > > &a_c2c, vtkUnstructuredGrid *a_grid);
693 const char* VertexType2Str(char T);
694 char Str2VertexType(QString S);
695 const char* vertex_type(char T);
697 #endif