replaced abort() with EG_BUG
[engrid.git] / src / egvtkobject.h
blob690631c0890cf60083a46612d44a8cf612179bd9
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 <vtkCellLocator.h>
40 #include <QSettings>
41 #include <QSet>
42 #include <QVector>
44 #define VTK_SIMPLE_VERTEX 0
45 #define VTK_FIXED_VERTEX 1
46 #define VTK_FEATURE_EDGE_VERTEX 2
47 #define VTK_BOUNDARY_EDGE_VERTEX 3
49 class EgVtkObject
52 public: // data-types
54 typedef const QVector<vtkIdType>& l2g_t;
55 typedef const QVector<int>& g2l_t;
56 typedef const QVector<QVector<int> >& l2l_t;
58 private: // methods
60 void addToC2C
62 vtkIdType id_cell,
63 QVector<int> &_cells,
64 QVector<QVector<int> > &c2c,
65 int j,
66 vtkIdList *nds,
67 vtkIdList *cls,
68 vtkUnstructuredGrid *grid
71 void addToN2N
73 QVector<QSet<int> > &n2n,
74 int n1,
75 int n2
78 void createNodeField(vtkUnstructuredGrid *grid, QString field_name, QString type_name, int Nnodes, bool overwrite = false);
79 void createCellField(vtkUnstructuredGrid *grid, QString field_name, QString type_name, int Ncells, bool overwrite = false);
81 protected: // attributes
83 QSet<int> m_BoundaryCodes;
84 static int DebugLevel;
86 protected: // methods
88 /**
89 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
90 * Version for int variables
92 int getSet(QString group, QString key, int value, int& variable);
94 /**
95 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
96 * Version for double variables
98 double getSet(QString group, QString key, double value, double& variable);
101 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
102 * Version for bool variables
104 bool getSet(QString group, QString key, bool value, bool& variable);
107 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
108 * Version for string variables
110 QString getSet(QString group, QString key, QString value, QString& variable);
113 * Update the cell index array.
115 void UpdateCellIndex(vtkUnstructuredGrid *grid);
118 * Update the point index array.
120 void UpdateNodeIndex(vtkUnstructuredGrid *grid);
123 * Compute normal vectors on nodes and cells of a subset of a grid.
124 * The paramters nodes and cells must be consistent; this means the nodes
125 * represent exactly (not more, not less) the nodes forming the cells.
126 * @param cell_normals On return, this will contain the cell normals (same order as cells)
127 * @param node_normals On return, this will contain the cell normals (same order as cells)
128 * @param cells The cells to compute the normals of
129 * @param nodes The nodes to compute the normals of
130 * @param grid The grid to operate on
132 void computeNormals(QVector<vec3_t> &cell_normals, QVector<vec3_t> &node_normals, QVector<vtkIdType> &cells, QVector<vtkIdType> &nodes, vtkUnstructuredGrid *grid);
135 * Create a mapping from global node indices to the indeces of a subset of nodes.
136 * @param nodes The subset of nodes.
137 * @param _nodes On return, this will contain the mapping.
138 * @param grid The grid to operate on.
140 void createNodeMapping(QVector<vtkIdType> &nodes, QVector<int> &_nodes, vtkUnstructuredGrid *grid);
143 * Create a mapping from global cell indices to the indices of a subset of cells.
144 * @param cells The subset of cells.
145 * @param _cells On return, this will contain the mapping.
146 * @param grid The grid to operate on.
148 void createCellMapping(QVector<vtkIdType> &cells, QVector<int> &_cells, 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(QVector<QSet<int> > &bcs, vtkUnstructuredGrid *grid);
160 * Create a node to cell structure for a given set of cells and nodes.
161 * This creates a vector of sets which might have performance issues.
162 * @param cells the subset of cells
163 * @param nodes the subset of nodes
164 * @param _nodes the reverse mapping for the nodes
165 * @param n2c On return, this will hold the node to cell structure
166 * @param grid The grid to operate on
168 void createNodeToCell(QVector<vtkIdType> &cells, QVector<vtkIdType> &nodes, QVector<int> &_nodes, QVector<QSet<int> > &n2c, vtkUnstructuredGrid *grid);
171 * Create a node to cell structure for a given set of cells and nodes.
172 * This creates a vector of vectors.
173 * @param cells the subset of cells
174 * @param nodes the subset of nodes
175 * @param _nodes the reverse mapping for the nodes
176 * @param n2c On return, this will hold the node to cell structure
177 * @param grid The grid to operate on
179 void createNodeToCell(QVector<vtkIdType> &cells, QVector<vtkIdType> &nodes, QVector<int> &_nodes, QVector<QVector<int> > &n2c, vtkUnstructuredGrid *grid);
182 * Create a node to node structure for a given set of cells and nodes.
183 * This creates a vector of sets which might have performance issues.
184 * @param cells the subset of cells
185 * @param nodes the subset of nodes
186 * @param _nodes the reverse mapping for the nodes
187 * @param n2n On return, this will hold the node to node structure
188 * @param grid The grid to operate on
190 void createNodeToNode(QVector<vtkIdType> &cells, QVector<vtkIdType> &nodes, QVector<int> &_nodes, QVector<QSet<int> > &n2n, vtkUnstructuredGrid *grid);
193 * Create a node to node structure for a given set of cells and nodes.
194 * This creates a vector of vectors.
195 * @param cells the subset of cells
196 * @param nodes the subset of nodes
197 * @param _nodes the reverse mapping for the nodes
198 * @param n2n On return, this will hold the node to node structure
199 * @param grid The grid to operate on
201 void createNodeToNode(QVector<vtkIdType> &cells, QVector<vtkIdType> &nodes, QVector<int> &_nodes, QVector<QVector<int> > &n2n, vtkUnstructuredGrid *grid);
204 * Extract the nodes which are part of a given set of cells.
205 * @param cells the subset of cells
206 * @param nodes On return, this will contain the nodes that correspond to the subset of cells
207 * @param grid The grid to operate on
209 template <class C>
210 void getNodesFromCells(const C &cells, QVector<vtkIdType> &nodes, vtkUnstructuredGrid *grid);
213 * Check if a cell is a volume cell.
214 * @param cellId The id fof the cell in question
215 * @param grid The grid to operate on
216 * @return true if the cell represents a volume and false if not
218 bool isVolume(vtkIdType id_cell, vtkUnstructuredGrid *grid);
222 * Check if a cell is a surface 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 surface and false if not
227 bool isSurface(vtkIdType id_cell, vtkUnstructuredGrid *grid);
230 * Get all volume cells of a grid.
231 * @param cells On return this will hold the Ids of all volume cells.
232 * @param grid The grid to operate on.
234 void getAllVolumeCells(QVector<vtkIdType> &cells, vtkUnstructuredGrid *grid);
237 * Get all cells of a grid.
238 * @param cells On return this will hold the Ids of all cells.
239 * @param grid The grid to operate on.
241 void getAllCells(QVector<vtkIdType> &cells, vtkUnstructuredGrid *grid);
244 * Get all cells of a grid and a specific type.
245 * @param type The type of the cells (e.g. VTK_TETRA, VTK_TRIANGLE, etc.)
246 * @param cells On return this will hold the Ids of all cells.
247 * @param grid The grid to operate on.
249 void getAllCellsOfType(vtkIdType type, QVector<vtkIdType> &cells, vtkUnstructuredGrid *grid);
252 * Get all surface cells of a grid.
253 * @param cells On return this will hold the Ids of all surface cells.
254 * @param grid The grid to operate on.
256 void getAllSurfaceCells(QVector<vtkIdType> &cells, vtkUnstructuredGrid *grid);
259 * Get all surface cells of a grid with a specific boundary condition.
260 * @param bcs The set of boundary conditions
261 * @param cells On return this will hold the Ids of the surface cells.
262 * @param grid The grid to operate on.
264 void getSurfaceCells(QSet<int> &bcs, QVector<vtkIdType> &cells, vtkUnstructuredGrid *grid);
267 * Create a cell neighbourship list for a subset grid.
268 * This has been implemented using VTK's vtkCellLinks structures.
269 * @param cells the subset of cells
270 * @param c2c On return this will hold the neighbourship list
271 * @param grid The grid to operate on.
273 void createCellToCell(QVector<vtkIdType> &cells, QVector<QVector<int> > &c2c, vtkUnstructuredGrid *grid);
276 * Insert a subset of a grid into a vtkPolyData structure.
277 * This is can be used in order to make use of many readily available
278 * operations within VTK; one example is smoothing.
279 * Cell index and node index arrays will be created and passed to the
280 * poly data structure. Thus any purely geometric effect (no topology change)
281 * can be directly reintroduced into the vtkUnstructuredGrid.
282 * @param cells the subset of cells
283 * @param pdata the vtkPolyData to add the nodes and cells to
284 * @param grid The grid to operate on.
286 void addToPolyData(QVector<vtkIdType> &cells, vtkPolyData *pdata, vtkUnstructuredGrid *grid);
289 * Copy the attributes from an input to an output cell.
290 * @param old_grid the input grid
291 * @param oldId the existing input cell
292 * @param new_grid the output grid
293 * @param newId the new output cell
295 void copyCellData(vtkUnstructuredGrid *old_grid, vtkIdType oldId, vtkUnstructuredGrid *new_grid, vtkIdType newId);
298 * Copy the attributes from an input to an output node.
299 * @param old_grid the input grid
300 * @param oldId the existing input node
301 * @param new_grid the output grid
302 * @param newId the new output node
304 void copyNodeData(vtkUnstructuredGrid *old_grid, vtkIdType oldId, vtkUnstructuredGrid *new_grid, vtkIdType newId);
307 * Create the basic fields on a given grid.
308 * Care should be taken with the overwrite parameter; if it is set to <i>false</i>
309 * and the field does not have the correct size it can lead to <i>segmentation faults</i>.
310 * @param Ncells the number of output cells
311 * @param Nnodes the number of output nodes
312 * @param overwrite f set to true existing fields will be re-created
314 void createBasicFields(vtkUnstructuredGrid *grid, vtkIdType Ncells, vtkIdType Nnodes, bool overwrite = false);
317 * Create the basic cell fields on a given grid.
318 * Care should be taken with the overwrite parameter; if it is set to <i>false</i>
319 * and the field does not have the correct size it can lead to <i>segmentation faults</i>.
320 * @param Ncells the number of output cells
321 * @param overwrite f set to true existing fields will be re-created
323 void createBasicCellFields(vtkUnstructuredGrid *grid, vtkIdType Ncells, bool overwrite = false);
326 * Create the basic node 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 Nnodes the number of output nodes
330 * @param overwrite f set to true existing fields will be re-created
332 void createBasicNodeFields(vtkUnstructuredGrid *grid, vtkIdType Nnodes, bool overwrite = false);
335 * Allocate memory for a grid. This method will also create the basic
336 * attribute fields (e.g. "cell_code").
337 * @param grid the grid for which to allocate memory
338 * @param Ncells the number of output cells
339 * @param Nnodes the number of output nodes
340 * @param create_fields flag to determine if node and cell data shall be created
342 void allocateGrid(vtkUnstructuredGrid *grid, vtkIdType Ncells, vtkIdType Nnodes, bool create_fields = true);
345 * Get the names of all node (point) attributes (fields) of a VTK grid
346 * @param field_names On return this vector will contain the names of all fields
347 * @param grid the grid
349 void getAllNodeDataNames(QVector<QString> &field_names, vtkUnstructuredGrid *grid);
352 * Get the names of all cell attributes (fields) of a VTK grid
353 * @param field_names On return this vector will contain the names of all fields
354 * @param grid the grid
356 void getAllCellDataNames(QVector<QString> &field_names, vtkUnstructuredGrid *grid);
359 * Compute the intersection of two Q containers.
360 * This will return a set.
361 * @param set1 the first container
362 * @param set2 the second container
363 * @param inters on return this will hold the intersection
365 template <class C1, class C2>
366 void qcontIntersection(const C1& c1, const C2& c2, QSet<typename C1::value_type> &inters);
369 * Compute the intersection of two Q containers.
370 * This will return a vector.
371 * @param set1 the first container
372 * @param set2 the second container
373 * @param inters on return this will hold the intersection
375 template <class C1, class C2>
376 void qcontIntersection(const C1& c1, const C2& c2, QVector<typename C1::value_type> &inters);
379 * Compute the centre of a cell
380 * @param grid the grid to use
381 * @param id_cell the cell of which to compute the centre
382 * @return the centre of the cell
384 vec3_t cellCentre(vtkUnstructuredGrid *grid, vtkIdType id_cell);
387 * Get the cells of a grid that are not part of a given set of cells.
388 * @param grid the grid to use
389 * @param cells the given set of cells
390 * @param rest_cells on return this will hold the rest of all cells of the grid (not part of cells)
392 void getRestCells(vtkUnstructuredGrid *grid, const QVector<vtkIdType> &cells, QVector<vtkIdType> &rest_cells);
395 * Find the corresponding volume cell of a surface cell
396 * @param grid the grid to use
397 * @param id_surf the id of the surface cell
398 * @param n2n the node to cell structure for this grid
399 * @return the id of the corresponding volume cell (or -1 if not found)
401 vtkIdType findVolumeCell(vtkUnstructuredGrid *grid, vtkIdType id_surf, g2l_t _nodes, l2g_t cells, g2l_t _cells, l2l_t n2c);
404 * Copy "src" grid to "dst" grid. Allocate "dst" so that it fits the data of "src".
405 * @param src a pointer to the source grid
406 * @param dst a pointer to the destination grid
408 void makeCopy(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst);
411 * Copy a part of "src" grid to "dst" grid. Allocate "dst" so that it fits the data to be copied.
412 * @param src a pointer to the source grid
413 * @param dst a pointer to the destination grid
414 * @param cells a container with the cells to be copied
416 template <class C>
417 void makeCopy(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst, const C &cells);
420 * Copy "src" grid to "dst" grid. DO NOT allocate "dst" so that it fits the data of "src".
421 * Allocation is left for the user to do.
422 * @param src a pointer to the source grid
423 * @param dst a pointer to the destination grid
425 void makeCopyNoAlloc(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst);
428 * Change the orientation of a face.
429 * @param grid the grid to use
430 * @param id_face the id of the face to change
432 void reorientateFace(vtkUnstructuredGrid *grid, vtkIdType id_face);
435 * Reset face orientation to original orientation.
436 * @param grid the grid with the faces
438 void resetOrientation(vtkUnstructuredGrid *grid);
440 void createIndices(vtkUnstructuredGrid *grid);
443 * Get the boundary condition of a boundary code.
444 * @param bc the boundary code
445 * @return the boundary condition
447 BoundaryCondition getBC(int bc);
450 * Save the subgrid defined by cls from grid.
451 * @param grid The source grid
452 * @param cls The cells to extract
453 * @param file_name The file to save to
455 template <class C>
456 void writeCells(vtkUnstructuredGrid *grid, const C &cls, QString file_name);
459 * Get the SubGrid defined by cls from grid. The function takes care of allocation for SubGrid.
460 * @param grid The source grid
461 * @param cls The cells to extract
462 * @param SubGrid The SubGrid to create
464 template <class C>
465 void getSubGrid(vtkUnstructuredGrid *grid, const C &cls, vtkUnstructuredGrid *SubGrid);
467 void writeGrid(vtkUnstructuredGrid *grid, QString name);
470 * Get a file name without extension.
471 * @param file_name the full name (with extension)
472 * @return the name without the extension
474 QString stripFromExtension(QString file_name);
477 * Get the extension of a file name
478 * @param file_name the full name (with extension)
479 * @return the extension
481 QString getExtension(QString file_name);
483 public: // methods
485 EgVtkObject() { DebugLevel = 0; }
487 void setBoundaryCodes(const QSet<int> &bcs);
488 void setDebugLevel(int a_DebugLevel) { DebugLevel = a_DebugLevel; }
492 //End of class EgVtkObject
494 template <class C1, class C2>
495 void EgVtkObject::qcontIntersection(const C1& c1, const C2& c2, QSet<typename C1::value_type> &inters)
497 inters.clear();
498 foreach (typename C1::value_type t1, c1) {
499 foreach (typename C2::value_type t2, c2) {
500 if (t1 == t2) {
501 inters.insert(t1);
507 template <class C1, class C2>
508 void EgVtkObject::qcontIntersection(const C1& c1, const C2& c2, QVector<typename C1::value_type> &inters)
510 QSet<typename C1::value_type> inters_set;
511 qcontIntersection(c1, c2, inters_set);
512 inters.resize(inters_set.size());
513 qCopy(inters_set.begin(), inters_set.end(), inters.begin());
516 template <class C>
517 void EgVtkObject::getSubGrid(vtkUnstructuredGrid *grid, const C &cls, vtkUnstructuredGrid *SubGrid)
519 createIndices(grid);
520 QVector<vtkIdType> cells;
521 QVector<vtkIdType> nodes;
522 cells.resize(cls.size());
523 qCopy(cls.begin(), cls.end(), cells.begin());
524 getNodesFromCells(cells, nodes, grid);
525 allocateGrid(SubGrid, cells.size(), nodes.size());
526 vtkIdType id_new_node = 0;
527 QVector<vtkIdType> old2new(grid->GetNumberOfPoints(), -1);
528 foreach (vtkIdType id_node, nodes) {
529 vec3_t x;
530 grid->GetPoint(id_node, x.data());
531 SubGrid->GetPoints()->SetPoint(id_new_node, x.data());
532 old2new[id_node] = id_new_node;
533 copyNodeData(grid, id_node, SubGrid, id_new_node);
534 ++id_new_node;
536 foreach (vtkIdType id_cell, cells) {
537 vtkIdType N_pts, *pts;
538 vtkIdType type_cell = grid->GetCellType(id_cell);
539 grid->GetCellPoints(id_cell, N_pts, pts);
540 QVector<vtkIdType> new_pts(N_pts);
541 for (int i_pts = 0; i_pts < N_pts; ++i_pts) {
542 new_pts[i_pts] = old2new[pts[i_pts]];
544 vtkIdType id_new_cell = SubGrid->InsertNextCell(type_cell, N_pts, new_pts.data());
545 copyCellData(grid, id_cell, SubGrid, id_new_cell);
549 template <class C>
550 void EgVtkObject::writeCells(vtkUnstructuredGrid *grid, const C &cls, QString file_name)
552 EG_VTKSP(vtkUnstructuredGrid,SubGrid);
553 getSubGrid(grid,cls,SubGrid);
555 EG_VTKSP(vtkXMLUnstructuredGridWriter,vtu);
556 vtu->SetFileName(qPrintable(file_name));
557 vtu->SetDataModeToBinary();
558 vtu->SetInput(SubGrid);
559 vtu->Write();
562 template <class C>
563 void EgVtkObject::getNodesFromCells(const C& cells, QVector<vtkIdType> &nodes, vtkUnstructuredGrid *grid)
565 QSet<vtkIdType> ex_nodes;
566 vtkIdType id_cell;
567 foreach(id_cell, cells) {
568 vtkIdType *pts;
569 vtkIdType Npts;
570 grid->GetCellPoints(id_cell, Npts, pts);
571 for (int i = 0; i < Npts; ++i) {
572 ex_nodes.insert(pts[i]);
575 nodes.resize(ex_nodes.size());
577 int j = 0;
578 vtkIdType i;
579 foreach(i,ex_nodes) {
580 nodes[j] = i;
581 ++j;
586 template <class C>
587 void EgVtkObject::makeCopy(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst, const C& cells)
589 QVector<vtkIdType> nodes;
590 getNodesFromCells(cells, nodes, src);
591 allocateGrid(dst, cells.size(), nodes.size());
592 vtkIdType id_new_node = 0;
593 foreach (vtkIdType id_node, nodes) {
594 vec3_t x;
595 src->GetPoints()->GetPoint(id_node, x.data());
596 dst->GetPoints()->SetPoint(id_node, x.data());
597 copyNodeData(src, id_node, dst, id_new_node);
598 ++id_new_node;
600 foreach (vtkIdType id_cell, cells) {
601 vtkIdType N_pts, *pts;
602 vtkIdType type_cell = src->GetCellType(id_cell);
603 src->GetCellPoints(id_cell, N_pts, pts);
604 vtkIdType id_new_cell = dst->InsertNextCell(type_cell, N_pts, pts);
605 copyCellData(src, id_cell, dst, id_new_cell);
610 * Utility function that allows printing selected data from an vtkUnstructuredGrid to any ostream (includes ofstream objects)
611 * @param stream ostream object to print to
612 * @param grid vtkUnstructuredGrid you want to print data from
613 * @param npoints print number of points in the grid
614 * @param ncells print number of cells in the grid
615 * @param points print points in the grid
616 * @param cells print cells in the grid
618 int cout_grid(ostream &stream, vtkUnstructuredGrid *grid, bool npoints=true, bool ncells=true, bool points=false, bool cells=false);
620 ///////////////////////////////////////////
621 int addCell(vtkUnstructuredGrid* a_grid, vtkIdType A, vtkIdType B, vtkIdType C, int bc);
623 ///get number of the shortest side of the cell
624 int getShortestSide(vtkIdType a_id_cell,vtkUnstructuredGrid* a_grid);
625 ///get number of the longest side of the cell
626 int getLongestSide(vtkIdType a_id_cell,vtkUnstructuredGrid* a_grid);
627 ///sort sides by length
628 //QVector <vtkIdType> sortSidesByLength(vtkIdType a_id_cell,vtkUnstructuredGrid* a_grid);
630 ///get number of the edge corresponding to node1-node2
631 int getSide(vtkIdType a_id_cell,vtkUnstructuredGrid* a_grid,vtkIdType a_id_node1,vtkIdType a_id_node2);
633 QSet <int> complementary_bcs(QSet <int> &bcs, vtkUnstructuredGrid *a_grid, QVector <vtkIdType> &a_cells);
634 QString cell2str(vtkIdType id_cell,vtkUnstructuredGrid* grid);
635 Qt::CheckState int2CheckState(int a);
636 int CheckState2int(Qt::CheckState a);
638 ///////////////////////////////////////////
640 template <class T>
641 ostream &operator<<(ostream &out, QVector<T> const & vector)
643 int N=vector.size();
644 out<<"[";
645 for (int i = 0; i < N; ++i) {
646 out<<vector.at(i);
647 if(i!=N-1) out<<",";
649 out<<"]";
650 return(out);
653 template <class T>
654 ostream &operator<<(ostream &out, QSet<T> const & set )
656 out << "[ ";
657 foreach (T value, set) out << value << " ";
658 out << "]";
659 return(out);
662 template <class T>
663 ostream &operator<<(ostream &out, QVector<QSet<T> > & vector)
665 int N=vector.size();
666 out<<"[";
667 for (int i = 0; i < N; ++i) {
668 QSet<T> set=vector.at(i);
669 out<<set;
670 if(i!=N-1) out<<",";
672 out<<"]";
673 return(out);
676 template <class T>
677 ostream &operator<<(ostream &out, QVector<QVector<T> > & vector)
679 int N=vector.size();
680 out<<"[";
681 for (int i = 0; i < N; ++i) {
682 QVector<T> subvector=vector.at(i);
683 out<<subvector;
684 if(i!=N-1) out<<",";
686 out<<"]";
687 return(out);
690 template <class T1, class T2>
691 ostream &operator<<(ostream &out, QMap<T1,T2> & map)
693 QMapIterator<T1, T2> i(map);
694 out<<"[";
695 while (i.hasNext()) {
696 i.next();
697 out << " [" << i.key() << ": " << i.value() << "]";
699 out<<"]";
700 return(out);
703 template <class T1, class T2>
704 ostream &operator<<(ostream &out, QVector < pair<T1,T2> > & vector)
706 int N=vector.size();
707 out<<"[";
708 for (int i = 0; i < N; ++i) {
709 out<<"<";
710 out<<vector.at(i).first;
711 out<<",";
712 out<<vector.at(i).second;
713 out<<">";
714 if(i!=N-1) out<<",";
716 out<<"]";
717 return(out);
720 template <class T>
721 QVector <T> Set2Vector(QSet <T> a_set, bool a_sort)
723 QVector <T> l_vector(a_set.size());
724 qCopy(a_set.begin(),a_set.end(),l_vector.begin());
725 if(a_sort) qSort(l_vector.begin(),l_vector.end());
726 return(l_vector);
729 template <class T>
730 QSet <T> Vector2Set(QVector <T> a_vector, bool a_sort)
732 QSet <T> l_set;
733 foreach(T element, a_vector) l_set.insert(element);
734 if(a_sort) qSort(l_set.begin(),l_set.end());
735 return(l_set);
738 template <class T>
739 bool duplicates(QVector <T> a_vector)
741 QSet <T> l_set;
742 foreach(T element, a_vector) l_set.insert(element);
743 return l_set.size()!=a_vector.size();
746 ///////////////////////////////////////////
747 pair<vtkIdType,vtkIdType> OrderedPair(vtkIdType a, vtkIdType b);
749 const char* VertexType2Str(char T);
750 char Str2VertexType(QString S);
752 #endif