feature magic defaults to 1 now
[engrid.git] / src / libengrid / egvtkobject.h
blob660c706aebb284bfc59f00080f6b426a0f313c55
1 //
2 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 // + +
4 // + This file is part of enGrid. +
5 // + +
6 // + Copyright 2008-2013 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 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
22 //
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 EG_SIMPLE_VERTEX 0
48 #define EG_FEATURE_EDGE_VERTEX 1
49 #define EG_BOUNDARY_EDGE_VERTEX 2
50 #define EG_FEATURE_CORNER_VERTEX 3
51 #define EG_FIXED_VERTEX 4
53 class EgVtkObject
56 public: // data-types
58 typedef const QVector<vtkIdType>& l2g_t;
59 typedef const QVector<int>& g2l_t;
60 typedef const QVector<QVector<int> >& l2l_t;
62 private: // methods
64 void addToC2C
66 vtkIdType id_cell,
67 QVector<int> &_cells,
68 QVector<QVector<int> > &c2c,
69 int j,
70 vtkIdList *nds,
71 vtkIdList *cls,
72 vtkUnstructuredGrid *grid
75 void addToN2N
77 QVector<QSet<int> > &n2n,
78 int n1,
79 int n2
82 void createNodeField(vtkUnstructuredGrid *grid, QString field_name, QString type_name, int Nnodes, bool overwrite = false);
83 void createCellField(vtkUnstructuredGrid *grid, QString field_name, QString type_name, int Ncells, bool overwrite = false);
85 protected: // attributes
87 QSet<int> m_BoundaryCodes;
88 static int DebugLevel;
90 protected: // methods
92 /**
93 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
94 * Version for int variables
96 int getSet(QString group, QString key, int value, int& variable);
98 /**
99 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
100 * Version for double variables
102 double getSet(QString group, QString key, double value, double& variable);
105 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
106 * Version for bool variables
108 bool getSet(QString group, QString key, bool value, bool& variable);
111 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
112 * Version for string variables
114 QString getSet(QString group, QString key, QString value, QString& variable);
117 * if key=value pair not found in settings file, write it + read key value from settings file and assign it to variable
118 * Version for string variables.
120 QString getSet(QString group, QString key, QString value, QString& variable, int type);
123 * Update the cell index array.
125 void UpdateCellIndex(vtkUnstructuredGrid *grid);
128 * Update the point index array.
130 void UpdateNodeIndex(vtkUnstructuredGrid *grid);
133 * Compute normal vectors on nodes and cells of a subset of a grid.
134 * The paramters nodes and cells must be consistent; this means the nodes
135 * represent exactly (not more, not less) the nodes forming the cells.
136 * @param cell_normals On return, this will contain the cell normals (same order as cells)
137 * @param node_normals On return, this will contain the cell normals (same order as cells)
138 * @param cells The cells to compute the normals of
139 * @param nodes The nodes to compute the normals of
140 * @param grid The grid to operate on
142 void computeNormals(QVector<vec3_t> &cell_normals, QVector<vec3_t> &node_normals, QVector<vtkIdType> &cells, QVector<vtkIdType> &nodes, vtkUnstructuredGrid *grid);
145 * Create a mapping from global node indices to the indeces of a subset of nodes.
146 * @param nodes The subset of nodes.
147 * @param _nodes On return, this will contain the mapping.
148 * @param grid The grid to operate on.
150 void createNodeMapping(QVector<vtkIdType> &nodes, QVector<int> &_nodes, vtkUnstructuredGrid *grid);
153 * Create a mapping from global cell indices to the indices of a subset of cells.
154 * @param cells The subset of cells.
155 * @param _cells On return, this will contain the mapping.
156 * @param grid The grid to operate on.
158 void createCellMapping(QVector<vtkIdType> &cells, QVector<int> &_cells, vtkUnstructuredGrid *grid);
161 * Create a node to boundary condition ("cell_code") mapping.
162 * Only non-zero boundary conditions will be considered.
163 * @param bcs On return, this will hold the codes of all boundary elements that are
164 * attached to a node.
165 * @param grid The grid to operate on.
167 void createNodeToBcMapping(QVector<QSet<int> > &bcs, vtkUnstructuredGrid *grid);
170 * Create a node to cell structure for a given set of cells and nodes.
171 * This creates a vector of sets which might have performance issues.
172 * @param cells the subset of cells
173 * @param nodes the subset of nodes
174 * @param _nodes the reverse mapping for the nodes
175 * @param n2c On return, this will hold the node to cell structure
176 * @param grid The grid to operate on
178 void createNodeToCell(QVector<vtkIdType> &cells, QVector<vtkIdType> &nodes, QVector<int> &_nodes, QVector<QSet<int> > &n2c, vtkUnstructuredGrid *grid);
181 * Create a node to cell structure for a given set of cells and nodes.
182 * This creates a vector of vectors.
183 * @param cells the subset of cells
184 * @param nodes the subset of nodes
185 * @param _nodes the reverse mapping for the nodes
186 * @param n2c On return, this will hold the node to cell structure
187 * @param grid The grid to operate on
189 void createNodeToCell(QVector<vtkIdType> &cells, QVector<vtkIdType> &nodes, QVector<int> &_nodes, QVector<QVector<int> > &n2c, vtkUnstructuredGrid *grid);
192 * Create a node to node structure for a given set of cells and nodes.
193 * This creates a vector of sets which might have performance issues.
194 * @param cells the subset of cells
195 * @param nodes the subset of nodes
196 * @param _nodes the reverse mapping for the nodes
197 * @param n2n On return, this will hold the node to node structure
198 * @param grid The grid to operate on
200 void createNodeToNode(QVector<vtkIdType> &cells, QVector<vtkIdType> &nodes, QVector<int> &_nodes, QVector<QSet<int> > &n2n, vtkUnstructuredGrid *grid);
203 * Create a node to node structure for a given set of cells and nodes.
204 * This creates a vector of vectors.
205 * @param cells the subset of cells
206 * @param nodes the subset of nodes
207 * @param _nodes the reverse mapping for the nodes
208 * @param n2n On return, this will hold the node to node structure
209 * @param grid The grid to operate on
211 void createNodeToNode(QVector<vtkIdType> &cells, QVector<vtkIdType> &nodes, QVector<int> &_nodes, QVector<QVector<int> > &n2n, vtkUnstructuredGrid *grid);
214 * Extract the nodes which are part of a given set of cells.
215 * @param cells the subset of cells
216 * @param nodes On return, this will contain the nodes that correspond to the subset of cells
217 * @param grid The grid to operate on
219 template <class C>
220 void getNodesFromCells(const C &cells, QVector<vtkIdType> &nodes, vtkUnstructuredGrid *grid);
223 * Check if a cell is a volume cell.
224 * @param cellId The id fof the cell in question
225 * @param grid The grid to operate on
226 * @return true if the cell represents a volume and false if not
228 bool isVolume(vtkIdType id_cell, vtkUnstructuredGrid *grid);
232 * Check if a cell is a surface cell.
233 * @param cellId The id fof the cell in question
234 * @param grid The grid to operate on
235 * @return true if the cell represents a surface and false if not
237 bool isSurface(vtkIdType id_cell, vtkUnstructuredGrid *grid);
240 * Get all volume cells of a grid.
241 * @param cells On return this will hold the Ids of all volume cells.
242 * @param grid The grid to operate on.
244 void getAllVolumeCells(QVector<vtkIdType> &cells, vtkUnstructuredGrid *grid);
247 * Get all cells of a grid.
248 * @param cells On return this will hold the Ids of all cells.
249 * @param grid The grid to operate on.
251 void getAllCells(QVector<vtkIdType> &cells, vtkUnstructuredGrid *grid);
254 * Get all cells of a grid and a specific type.
255 * @param type The type of the cells (e.g. VTK_TETRA, VTK_TRIANGLE, etc.)
256 * @param cells On return this will hold the Ids of all cells.
257 * @param grid The grid to operate on.
259 void getAllCellsOfType(vtkIdType type, QVector<vtkIdType> &cells, vtkUnstructuredGrid *grid);
262 * Get all surface cells of a grid.
263 * @param cells On return this will hold the Ids of all surface cells.
264 * @param grid The grid to operate on.
266 void getAllSurfaceCells(QVector<vtkIdType> &cells, vtkUnstructuredGrid *grid);
269 * Get all surface cells of a grid with a specific boundary condition.
270 * @param bcs The set of boundary conditions
271 * @param cells On return this will hold the Ids of the surface cells.
272 * @param grid The grid to operate on.
274 void getSurfaceCells(QSet<int> &bcs, QVector<vtkIdType> &cells, vtkUnstructuredGrid *grid);
277 * Create a cell neighbourship list for a subset grid.
278 * This has been implemented using VTK's vtkCellLinks structures.
279 * @param cells the subset of cells
280 * @param c2c On return this will hold the neighbourship list
281 * @param grid The grid to operate on.
283 void createCellToCell(QVector<vtkIdType> &cells, QVector<QVector<int> > &c2c, vtkUnstructuredGrid *grid);
286 * Insert a subset of a grid into a vtkPolyData structure.
287 * This is can be used in order to make use of many readily available
288 * operations within VTK; one example is smoothing.
289 * Cell index and node index arrays will be created and passed to the
290 * poly data structure. Thus any purely geometric effect (no topology change)
291 * can be directly reintroduced into the vtkUnstructuredGrid.
292 * @param cells the subset of cells
293 * @param pdata the vtkPolyData to add the nodes and cells to
294 * @param grid The grid to operate on.
296 void addToPolyData(QVector<vtkIdType> &cells, vtkPolyData *pdata, vtkUnstructuredGrid *grid);
299 * Copy the attributes from an input to an output cell.
300 * @param old_grid the input grid
301 * @param oldId the existing input cell
302 * @param new_grid the output grid
303 * @param newId the new output cell
305 void copyCellData(vtkUnstructuredGrid *old_grid, vtkIdType oldId, vtkUnstructuredGrid *new_grid, vtkIdType newId);
308 * Copy the attributes from an input to an output node.
309 * @param old_grid the input grid
310 * @param oldId the existing input node
311 * @param new_grid the output grid
312 * @param newId the new output node
314 void copyNodeData(vtkUnstructuredGrid *old_grid, vtkIdType oldId, vtkUnstructuredGrid *new_grid, vtkIdType newId);
317 * Create the basic 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 Nnodes the number of output nodes
322 * @param overwrite f set to true existing fields will be re-created
324 void createBasicFields(vtkUnstructuredGrid *grid, vtkIdType Ncells, vtkIdType Nnodes, bool overwrite = false);
327 * Create the basic cell fields on a given grid.
328 * Care should be taken with the overwrite parameter; if it is set to <i>false</i>
329 * and the field does not have the correct size it can lead to <i>segmentation faults</i>.
330 * @param Ncells the number of output cells
331 * @param overwrite f set to true existing fields will be re-created
333 void createBasicCellFields(vtkUnstructuredGrid *grid, vtkIdType Ncells, bool overwrite = false);
336 * Create the basic node fields on a given grid.
337 * Care should be taken with the overwrite parameter; if it is set to <i>false</i>
338 * and the field does not have the correct size it can lead to <i>segmentation faults</i>.
339 * @param Nnodes the number of output nodes
340 * @param overwrite f set to true existing fields will be re-created
342 void createBasicNodeFields(vtkUnstructuredGrid *grid, vtkIdType Nnodes, bool overwrite = false);
345 * Allocate memory for a grid. This method will also create the basic
346 * attribute fields (e.g. "cell_code").
347 * @param grid the grid for which to allocate memory
348 * @param Ncells the number of output cells
349 * @param Nnodes the number of output nodes
350 * @param create_fields flag to determine if node and cell data shall be created
352 void allocateGrid(vtkUnstructuredGrid *grid, vtkIdType Ncells, vtkIdType Nnodes, bool create_fields = true);
355 * Get the names of all node (point) attributes (fields) of a VTK grid
356 * @param field_names On return this vector will contain the names of all fields
357 * @param grid the grid
359 void getAllNodeDataNames(QVector<QString> &field_names, vtkUnstructuredGrid *grid);
362 * Get the names of all cell attributes (fields) of a VTK grid
363 * @param field_names On return this vector will contain the names of all fields
364 * @param grid the grid
366 void getAllCellDataNames(QVector<QString> &field_names, vtkUnstructuredGrid *grid);
369 * Compute the intersection of two Q containers.
370 * This will return a set.
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, QSet<typename C1::value_type> &inters);
379 * Compute the intersection of two Q containers.
380 * This will return a vector.
381 * @param set1 the first container
382 * @param set2 the second container
383 * @param inters on return this will hold the intersection
385 template <class C1, class C2>
386 void qcontIntersection(const C1& c1, const C2& c2, QVector<typename C1::value_type> &inters);
389 * Compute the centre of a cell
390 * @param grid the grid to use
391 * @param id_cell the cell of which to compute the centre
392 * @return the centre of the cell
394 vec3_t cellCentre(vtkUnstructuredGrid *grid, vtkIdType id_cell);
397 * Compute the angle between two faces.
398 * @param grid the grid to use
399 * @param id_face1 index of the first face
400 * @param id_face2 index of the second face
401 * @return the angle between the faces (M_PI for a flat surface)
403 double faceAngle(vtkUnstructuredGrid* grid, vtkIdType id_face1, vtkIdType id_face2);
406 * Get the cells of a grid that are not part of a given set of cells.
407 * @param grid the grid to use
408 * @param cells the given set of cells
409 * @param rest_cells on return this will hold the rest of all cells of the grid (not part of cells)
411 void getRestCells(vtkUnstructuredGrid *grid, const QVector<vtkIdType> &cells, QVector<vtkIdType> &rest_cells);
414 * Find the corresponding volume cell of a surface cell
415 * @param grid the grid to use
416 * @param id_surf the id of the surface cell
417 * @param n2n the node to cell structure for this grid
418 * @return the id of the corresponding volume cell (or -1 if not found)
420 vtkIdType findVolumeCell(vtkUnstructuredGrid *grid, vtkIdType id_surf, g2l_t _nodes, l2g_t cells, g2l_t _cells, l2l_t n2c);
423 * Copy "src" grid to "dst" grid. Allocate "dst" so that it fits the data of "src".
424 * @param src a pointer to the source grid
425 * @param dst a pointer to the destination grid
427 void makeCopy(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst, bool copy_data = true);
430 * Copy a part of "src" grid to "dst" grid. Allocate "dst" so that it fits the data to be copied.
431 * @param src a pointer to the source grid
432 * @param dst a pointer to the destination grid
433 * @param cells a container with the cells to be copied
435 template <class C>
436 void makeCopy(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst, const C &cells);
439 * Copy "src" grid to "dst" grid. DO NOT allocate "dst" so that it fits the data of "src".
440 * Allocation is left for the user to do.
441 * @param src a pointer to the source grid
442 * @param dst a pointer to the destination grid
444 void makeCopyNoAlloc(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst);
447 * Change the orientation of a face.
448 * @param grid the grid to use
449 * @param id_face the id of the face to change
451 void reorientateFace(vtkUnstructuredGrid *grid, vtkIdType id_face);
454 * Reset face orientation to original orientation.
455 * @param grid the grid with the faces
457 void resetOrientation(vtkUnstructuredGrid *grid);
459 void createIndices(vtkUnstructuredGrid *grid);
462 * Get the boundary condition of a boundary code.
463 * @param bc the boundary code
464 * @return the boundary condition
466 BoundaryCondition getBC(int bc);
469 * Save the subgrid defined by cls from grid.
470 * @param grid The source grid
471 * @param cls The cells to extract
472 * @param file_name The file to save to
474 template <class C>
475 void writeCells(vtkUnstructuredGrid *grid, const C &cls, QString file_name);
478 * Get the SubGrid defined by cls from grid. The function takes care of allocation for SubGrid.
479 * @param grid The source grid
480 * @param cls The cells to extract
481 * @param SubGrid The SubGrid to create
483 template <class C>
484 void getSubGrid(vtkUnstructuredGrid *grid, const C &cls, vtkUnstructuredGrid *SubGrid);
488 * Save grid to file filename.
489 * @param grid The source grid
490 * @param filename Name of the file to save to
492 void writeGrid(vtkUnstructuredGrid *grid, QString filename);
495 * Get a file name without extension.
496 * @param file_name the full name (with extension)
497 * @return the name without the extension
499 QString stripFromExtension(QString file_name);
502 * Get the extension of a file name
503 * @param file_name the full name (with extension)
504 * @return the extension
506 QString getExtension(QString file_name);
509 * Get a face of a cell
510 * @param grid the unstructured grid
511 * @param id_cell the index of the cell
512 * @param i_face the index of the face within the cell
513 * @param ids on return this will contain the nodes of the face
515 void getFaceOfCell(vtkUnstructuredGrid *grid, vtkIdType id_cell, int i_face, QVector<vtkIdType> &ids);
518 * Get the normal of a face of a volume cell.
519 * @param grid the unstructured grid
520 * @param id_cell the index of the cell
521 * @param i_face the index of the face within the cell
522 * @return the normal vector (absolute value corresponds to the area)
524 vec3_t getNormalOfCell(vtkUnstructuredGrid *grid, vtkIdType id_cell, int i_face);
527 * Get an edge of a face/cell
528 * @param grid the unstructured grid
529 * @param id_cell the index of the cell
530 * @param i_edge the index of the edge within the cell
531 * @param ids on return this will contain the nodes of the edge
533 void getEdgeOfCell(vtkUnstructuredGrid *grid, vtkIdType id_cell, int i_edge, QVector<vtkIdType> &ids);
536 * Get all boundary codes fo a grid.
537 * @param grid the grid to extract the boundaru codes from
538 * @return a set with all boundary codes
540 QSet<int> getAllBoundaryCodes(vtkUnstructuredGrid *grid);
543 * Check if a cell (face) contains a given node
544 * @param grid the unstructured grid
545 * @param id_cell the id of the cell to investigate
546 * @param id_node the id of the required node
547 * @return true if id_cell contains id_node
549 bool cellContainsNode(vtkUnstructuredGrid *grid, vtkIdType id_cell, vtkIdType id_node);
552 * Get all node indices which are shared by two cells.
553 * These cells can be surface or volume cells; also a combination
554 * of a volume and a surface cell is possible.
555 * @param grid the grid to use
556 * @param id_cell1 index of the first cell
557 * @param id_cell2 index of the second cell
558 * @param cont a generic Qt container which will hold the shared node indices on return
560 template <typename C>
561 void sharedNodesOfCells(vtkUnstructuredGrid* grid, vtkIdType id_cell1, vtkIdType id_cell2, C& cont);
563 template <class C> void createPolyData(const C &x, vtkPolyData *poly_data, bool closed_loop = false);
564 void createPolyDataC2C(vtkPolyData *poly_data, QVector<QVector<vtkIdType> > &c2c);
565 void createPolyDataN2C(vtkPolyData *poly_data, QVector<QSet<vtkIdType> > &n2c);
566 void createPolyDataN2N(vtkPolyData *poly_data, QVector<QSet<vtkIdType> > &n2n);
567 template <class C> double convexRatio(const C &x, vec3_t n_plane, bool closed_loop = false);
569 public: // methods
571 EgVtkObject() { DebugLevel = 0; }
573 void setBoundaryCodes(const QSet<int> &bcs);
574 QSet<int> getBoundaryCodes();
575 void setDebugLevel(int a_DebugLevel) { DebugLevel = a_DebugLevel; }
577 bool saveGrid( vtkUnstructuredGrid* a_grid, QString file_name );
579 vtkIdType addGrid(vtkUnstructuredGrid *main_grid, vtkUnstructuredGrid *grid_to_add, vtkIdType offset);
581 private:
582 void addVtkTypeInfo(vtkUnstructuredGrid* a_grid); ///< Add VTK type information to the grid (useful for visualisation with ParaView).
585 //End of class EgVtkObject
587 template <class C1, class C2>
588 void EgVtkObject::qcontIntersection(const C1& c1, const C2& c2, QSet<typename C1::value_type> &inters)
590 inters.clear();
591 foreach (typename C1::value_type t1, c1) {
592 foreach (typename C2::value_type t2, c2) {
593 if (t1 == t2) {
594 inters.insert(t1);
600 template <class C1, class C2>
601 void EgVtkObject::qcontIntersection(const C1& c1, const C2& c2, QVector<typename C1::value_type> &inters)
603 QSet<typename C1::value_type> inters_set;
604 qcontIntersection(c1, c2, inters_set);
605 inters.resize(inters_set.size());
606 qCopy(inters_set.begin(), inters_set.end(), inters.begin());
609 template <class C>
610 void EgVtkObject::getSubGrid(vtkUnstructuredGrid *grid, const C &cls, vtkUnstructuredGrid *SubGrid)
612 createIndices(grid);
613 QVector<vtkIdType> cells;
614 QVector<vtkIdType> nodes;
615 cells.resize(cls.size());
616 qCopy(cls.begin(), cls.end(), cells.begin());
617 getNodesFromCells(cells, nodes, grid);
618 allocateGrid(SubGrid, cells.size(), nodes.size());
619 vtkIdType id_new_node = 0;
620 QVector<vtkIdType> old2new(grid->GetNumberOfPoints(), -1);
621 foreach (vtkIdType id_node, nodes) {
622 vec3_t x;
623 grid->GetPoint(id_node, x.data());
624 SubGrid->GetPoints()->SetPoint(id_new_node, x.data());
625 old2new[id_node] = id_new_node;
626 copyNodeData(grid, id_node, SubGrid, id_new_node);
627 ++id_new_node;
629 foreach (vtkIdType id_cell, cells) {
630 vtkIdType N_pts, *pts;
631 vtkIdType type_cell = grid->GetCellType(id_cell);
632 grid->GetCellPoints(id_cell, N_pts, pts);
633 QVector<vtkIdType> new_pts(N_pts);
634 for (int i_pts = 0; i_pts < N_pts; ++i_pts) {
635 new_pts[i_pts] = old2new[pts[i_pts]];
637 vtkIdType id_new_cell = SubGrid->InsertNextCell(type_cell, N_pts, new_pts.data());
638 copyCellData(grid, id_cell, SubGrid, id_new_cell);
642 template <class C>
643 void EgVtkObject::writeCells(vtkUnstructuredGrid *grid, const C &cls, QString file_name)
645 qDebug()<<"Saving cells from grid as "<<file_name;
647 EG_VTKSP(vtkUnstructuredGrid,SubGrid);
648 getSubGrid(grid,cls,SubGrid);
650 EG_VTKSP(vtkXMLUnstructuredGridWriter,vtu);
651 vtu->SetFileName(qPrintable(file_name));
652 vtu->SetDataModeToBinary();
653 vtu->SetInput(SubGrid);
654 vtu->Write();
657 template <class C>
658 void EgVtkObject::getNodesFromCells(const C& cells, QVector<vtkIdType> &nodes, vtkUnstructuredGrid *grid)
660 QSet<vtkIdType> ex_nodes;
661 vtkIdType id_cell;
662 foreach(id_cell, cells) {
663 vtkIdType *pts;
664 vtkIdType Npts;
665 grid->GetCellPoints(id_cell, Npts, pts);
666 for (int i = 0; i < Npts; ++i) {
667 ex_nodes.insert(pts[i]);
670 nodes.resize(ex_nodes.size());
672 int j = 0;
673 vtkIdType i;
674 foreach(i,ex_nodes) {
675 nodes[j] = i;
676 ++j;
681 template <class C>
682 void EgVtkObject::makeCopy(vtkUnstructuredGrid *src, vtkUnstructuredGrid *dst, const C& cells)
684 QVector<vtkIdType> nodes;
685 getNodesFromCells(cells, nodes, src);
686 allocateGrid(dst, cells.size(), nodes.size());
687 vtkIdType id_new_node = 0;
688 QVector<vtkIdType> old2new(src->GetNumberOfPoints(), -1);
689 foreach (vtkIdType id_node, nodes) {
690 vec3_t x;
691 src->GetPoints()->GetPoint(id_node, x.data());
692 dst->GetPoints()->SetPoint(id_new_node, x.data());
693 copyNodeData(src, id_node, dst, id_new_node);
694 old2new[id_node] = id_new_node;
695 ++id_new_node;
697 foreach (vtkIdType id_cell, cells) {
698 vtkIdType N_pts, *pts;
699 vtkIdType type_cell = src->GetCellType(id_cell);
700 src->GetCellPoints(id_cell, N_pts, pts);
701 QVector<vtkIdType> new_pts(N_pts);
702 for (int i = 0; i < N_pts; ++i) {
703 if (old2new[pts[i]] == -1) {
704 EG_BUG;
706 new_pts[i] = old2new[pts[i]];
708 vtkIdType id_new_cell = dst->InsertNextCell(type_cell, N_pts, new_pts.data());
709 copyCellData(src, id_cell, dst, id_new_cell);
713 template <class C>
714 void EgVtkObject::createPolyData(const C &x, vtkPolyData *poly_data, bool closed_loop)
716 int N = x.size();
717 if (closed_loop) {
718 --N;
720 EG_VTKSP(vtkPoints, points);
721 points->SetNumberOfPoints(N);
722 QVector<vtkIdType> pts(N);
723 for (int i = 0; i < N; ++i) {
724 points->SetPoint(i, x[i][0], x[i][1], x[i][2]);
725 pts[i] = i;
727 poly_data->Allocate(1);
728 poly_data->SetPoints(points);
729 poly_data->InsertNextCell(VTK_POLYGON, N, pts.data());
732 template <class C>
733 double EgVtkObject::convexRatio(const C &x, vec3_t n_plane, bool closed_loop)
735 double L_max = -1e99;
736 double L_min = 1e99;
737 int N = x.size();
738 if (closed_loop) {
739 --N;
741 for (int i = 0; i < N; ++i) {
742 for (int j = 0; j < N; ++j) {
743 int p1 = j;
744 int p2 = j+1;
745 if (j == N - 1 && !closed_loop) {
746 p2 = 0;
748 if (i != p1 && i != p2) {
749 vec3_t n = n_plane.cross(x[p2] - x[p1]);
750 n.normalise();
751 double L = (x[i] - x[j])*n;
752 L_max = max(L, L_max);
753 L_min = min(L, L_min);
757 return L_min/L_max;
760 template <typename C>
761 void EgVtkObject::sharedNodesOfCells(vtkUnstructuredGrid* grid, vtkIdType id_cell1, vtkIdType id_cell2, C& cont)
763 cont.clear();
764 EG_GET_CELL(id_cell1, grid);
765 for (int i = 0; i < num_pts; ++i) {
766 if (cellContainsNode(grid, id_cell2, pts[i])) {
767 cont << pts[i];
774 #endif