1 // Copyright (C) 2010-2016 CEA/DEN, EDF R&D, OPEN CASCADE
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
11 // Lesser General Public License for more details.
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
20 #include "SMDS_UnstructuredGrid.hxx"
21 #include "SMDS_Mesh.hxx"
22 #include "SMDS_MeshInfo.hxx"
23 #include "SMDS_Downward.hxx"
24 #include "SMDS_MeshVolume.hxx"
26 #include "utilities.h"
29 #include <vtkCellArray.h>
30 #include <vtkCellData.h>
31 #include <vtkCellLinks.h>
32 #include <vtkDoubleArray.h>
33 #include <vtkIdTypeArray.h>
34 #include <vtkUnsignedCharArray.h>
41 SMDS_CellLinks* SMDS_CellLinks::New()
43 return new SMDS_CellLinks();
46 void SMDS_CellLinks::ResizeForPoint(vtkIdType vtkID)
48 if ( vtkID > this->MaxId )
51 if ( vtkID >= this->Size )
52 vtkCellLinks::Resize( vtkID+SMDS_Mesh::chunkSize );
56 void SMDS_CellLinks::BuildLinks(vtkDataSet *data, vtkCellArray *Connectivity, vtkUnsignedCharArray* types)
58 // build links taking into account removed cells
60 vtkIdType numPts = data->GetNumberOfPoints();
61 vtkIdType j, cellId = 0;
62 unsigned short *linkLoc;
65 vtkIdType loc = Connectivity->GetTraversalLocation();
67 // traverse data to determine number of uses of each point
69 for (Connectivity->InitTraversal();
70 Connectivity->GetNextCell(npts,pts); cellId++)
72 if ( types->GetValue( cellId ) != VTK_EMPTY_CELL )
73 for (j=0; j < npts; j++)
75 this->IncrementLinkCount(pts[j]);
79 // now allocate storage for the links
80 this->AllocateLinks(numPts);
81 this->MaxId = numPts - 1;
83 // fill out lists with references to cells
84 linkLoc = new unsigned short[numPts];
85 memset(linkLoc, 0, numPts*sizeof(unsigned short));
88 for (Connectivity->InitTraversal();
89 Connectivity->GetNextCell(npts,pts); cellId++)
91 if ( types->GetValue( cellId ) != VTK_EMPTY_CELL )
92 for (j=0; j < npts; j++)
94 this->InsertCellReference(pts[j], (linkLoc[pts[j]])++, cellId);
98 Connectivity->SetTraversalLocation(loc);
101 SMDS_CellLinks::SMDS_CellLinks() :
106 SMDS_CellLinks::~SMDS_CellLinks()
110 SMDS_UnstructuredGrid* SMDS_UnstructuredGrid::New()
112 return new SMDS_UnstructuredGrid();
115 SMDS_UnstructuredGrid::SMDS_UnstructuredGrid() :
116 vtkUnstructuredGrid()
118 _cellIdToDownId.clear();
124 SMDS_UnstructuredGrid::~SMDS_UnstructuredGrid()
128 unsigned long SMDS_UnstructuredGrid::GetMTime()
130 unsigned long mtime = vtkUnstructuredGrid::GetMTime();
134 vtkPoints* SMDS_UnstructuredGrid::GetPoints()
136 // TODO erreur incomprehensible de la macro vtk GetPoints apparue avec la version paraview de fin aout 2010
140 int SMDS_UnstructuredGrid::InsertNextLinkedCell(int type, int npts, vtkIdType *pts)
145 if (type != VTK_POLYHEDRON)
146 return vtkUnstructuredGrid::InsertNextLinkedCell(type, npts, pts);
148 // --- type = VTK_POLYHEDRON
149 int cellid = this->InsertNextCell(type, npts, pts);
151 set<vtkIdType> setOfNodes;
155 for (int nf = 0; nf < nbfaces; nf++)
157 int nbnodes = pts[i];
159 for (int k = 0; k < nbnodes; k++)
161 setOfNodes.insert(pts[i]);
166 set<vtkIdType>::iterator it = setOfNodes.begin();
167 for (; it != setOfNodes.end(); ++it)
169 this->Links->ResizeCellList(*it, 1);
170 this->Links->AddCellReference(cellid, *it);
176 void SMDS_UnstructuredGrid::setSMDS_mesh(SMDS_Mesh *mesh)
181 void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize,
182 std::vector<int>& idCellsOldToNew, int newCellSize)
184 int alreadyCopied = 0;
188 // --- if newNodeSize, create a new compacted vtkPoints
190 vtkPoints *newPoints = vtkPoints::New();
191 newPoints->SetDataType(VTK_DOUBLE);
192 newPoints->SetNumberOfPoints(newNodeSize);
195 // rnv: to fix bug "21125: EDF 1233 SMESH: Degradation of precision in a test case for quadratic conversion"
196 // using double type for storing coordinates of nodes instead float.
197 int oldNodeSize = idNodesOldToNew.size();
200 while ( i < oldNodeSize )
202 // skip a hole if any
203 while ( i < oldNodeSize && idNodesOldToNew[i] < 0 )
206 // look for a block end
207 while ( i < oldNodeSize && idNodesOldToNew[i] >= 0 )
210 copyNodes(newPoints, idNodesOldToNew, alreadyCopied, startBloc, endBloc);
212 newPoints->Squeeze();
215 if (1/*newNodeSize*/)
217 this->SetPoints(newPoints);
222 // --- create new compacted Connectivity, Locations and Types
224 int oldCellSize = this->Types->GetNumberOfTuples();
226 if ( oldCellSize == newCellSize ) // no holes in elements
228 this->Connectivity->Squeeze();
229 this->Locations->Squeeze();
230 this->Types->Squeeze();
231 if ( this->FaceLocations )
233 this->FaceLocations->Squeeze();
234 this->Faces->Squeeze();
236 for ( int i = 0; i < oldCellSize; ++i )
237 idCellsOldToNew[i] = i;
240 vtkCellArray *newConnectivity = vtkCellArray::New();
241 newConnectivity->Initialize();
242 int oldCellDataSize = this->Connectivity->GetData()->GetSize();
243 newConnectivity->Allocate(oldCellDataSize);
245 vtkUnsignedCharArray *newTypes = vtkUnsignedCharArray::New();
246 newTypes->Initialize();
247 newTypes->SetNumberOfValues(newCellSize);
249 vtkIdTypeArray *newLocations = vtkIdTypeArray::New();
250 newLocations->Initialize();
251 newLocations->SetNumberOfValues(newCellSize);
253 // TODO some polyhedron may be huge (only in some tests)
254 vtkIdType tmpid[NBMAXNODESINCELL];
255 vtkIdType *pointsCell = &tmpid[0]; // --- points id to fill a new cell
259 while ( i < oldCellSize )
261 // skip a hole if any
262 while ( i < oldCellSize && this->Types->GetValue(i) == VTK_EMPTY_CELL )
265 // look for a block end
266 while ( i < oldCellSize && this->Types->GetValue(i) != VTK_EMPTY_CELL )
269 if ( endBloc > startBloc )
271 idCellsOldToNew, idNodesOldToNew,
272 newConnectivity, newLocations,
273 pointsCell, alreadyCopied,
276 newConnectivity->Squeeze();
278 if (vtkDoubleArray* diameters =
279 vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )) // Balls
281 for (int oldCellID = 0; oldCellID < oldCellSize; oldCellID++)
283 if (this->Types->GetValue(oldCellID) == VTK_EMPTY_CELL)
285 int newCellId = idCellsOldToNew[ oldCellID ];
286 if (newTypes->GetValue(newCellId) == VTK_POLY_VERTEX)
287 diameters->SetValue( newCellId, diameters->GetValue( oldCellID ));
291 if (this->FaceLocations)
293 vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New();
294 newFaceLocations->Initialize();
295 newFaceLocations->Allocate(newTypes->GetSize());
296 vtkIdTypeArray *newFaces = vtkIdTypeArray::New();
297 newFaces->Initialize();
298 newFaces->Allocate(this->Faces->GetSize());
299 for (int i = 0; i < oldCellSize; i++)
301 if (this->Types->GetValue(i) == VTK_EMPTY_CELL)
303 int newCellId = idCellsOldToNew[i];
304 if (newTypes->GetValue(newCellId) == VTK_POLYHEDRON)
306 newFaceLocations->InsertNextValue(newFaces->GetMaxId()+1);
307 int oldFaceLoc = this->FaceLocations->GetValue(i);
308 int nCellFaces = this->Faces->GetValue(oldFaceLoc++);
309 newFaces->InsertNextValue(nCellFaces);
310 for (int n=0; n<nCellFaces; n++)
312 int nptsInFace = this->Faces->GetValue(oldFaceLoc++);
313 newFaces->InsertNextValue(nptsInFace);
314 for (int k=0; k<nptsInFace; k++)
316 int oldpt = this->Faces->GetValue(oldFaceLoc++);
317 newFaces->InsertNextValue(idNodesOldToNew[oldpt]);
323 newFaceLocations->InsertNextValue(-1);
326 newFaceLocations->Squeeze();
328 this->SetCells(newTypes, newLocations, newConnectivity, newFaceLocations, newFaces);
329 newFaceLocations->Delete();
334 this->SetCells(newTypes, newLocations, newConnectivity, FaceLocations, Faces);
338 newLocations->Delete();
339 newConnectivity->Delete();
342 void SMDS_UnstructuredGrid::copyNodes(vtkPoints * newPoints,
343 std::vector<int>& idNodesOldToNew,
348 void *target = newPoints->GetVoidPointer(3 * alreadyCopied);
349 void *source = this->Points->GetVoidPointer(3 * start);
350 int nbPoints = end - start;
353 memcpy(target, source, 3 * sizeof(double) * nbPoints);
354 for (int j = start; j < end; j++)
355 idNodesOldToNew[j] = alreadyCopied++; // old vtkId --> new vtkId
359 void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes,
360 std::vector<int>& idCellsOldToNew,
361 std::vector<int>& idNodesOldToNew,
362 vtkCellArray* newConnectivity,
363 vtkIdTypeArray* newLocations,
364 vtkIdType* pointsCell,
369 for (int j = start; j < end; j++)
371 newTypes->SetValue(alreadyCopied, this->Types->GetValue(j));
372 idCellsOldToNew[j] = alreadyCopied; // old vtkId --> new vtkId
373 vtkIdType oldLoc = this->Locations->GetValue(j);
375 vtkIdType *oldPtsCell = 0;
376 this->Connectivity->GetCell(oldLoc, nbpts, oldPtsCell);
377 assert(nbpts < NBMAXNODESINCELL);
378 for (int l = 0; l < nbpts; l++)
380 int oldval = oldPtsCell[l];
381 pointsCell[l] = idNodesOldToNew[oldval];
383 /*int newcnt = */newConnectivity->InsertNextCell(nbpts, pointsCell);
384 int newLoc = newConnectivity->GetInsertLocation(nbpts);
385 newLocations->SetValue(alreadyCopied, newLoc);
390 int SMDS_UnstructuredGrid::CellIdToDownId(int vtkCellId)
392 if ((vtkCellId < 0) || (vtkCellId >= (int)_cellIdToDownId.size()))
396 return _cellIdToDownId[vtkCellId];
399 void SMDS_UnstructuredGrid::setCellIdToDownId(int vtkCellId, int downId)
401 // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
402 _cellIdToDownId[vtkCellId] = downId;
405 void SMDS_UnstructuredGrid::CleanDownwardConnectivity()
407 for (size_t i = 0; i < _downArray.size(); i++)
410 delete _downArray[i];
413 _cellIdToDownId.clear();
416 /*! Build downward connectivity: to do only when needed because heavy memory load.
417 * Downward connectivity is no more valid if vtkUnstructuredGrid is modified.
420 void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
422 MESSAGE("SMDS_UnstructuredGrid::BuildDownwardConnectivity");CHRONO(2);
423 // TODO calcul partiel sans edges
425 // --- erase previous data if any
427 this->CleanDownwardConnectivity();
429 // --- create SMDS_Downward structures (in _downArray vector[vtkCellType])
431 _downArray.resize(VTK_MAXTYPE + 1, 0);
433 _downArray[VTK_LINE] = new SMDS_DownEdge(this);
434 _downArray[VTK_QUADRATIC_EDGE] = new SMDS_DownQuadEdge(this);
435 _downArray[VTK_TRIANGLE] = new SMDS_DownTriangle(this);
436 _downArray[VTK_QUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
437 _downArray[VTK_BIQUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
438 _downArray[VTK_QUAD] = new SMDS_DownQuadrangle(this);
439 _downArray[VTK_QUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
440 _downArray[VTK_BIQUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
441 _downArray[VTK_TETRA] = new SMDS_DownTetra(this);
442 _downArray[VTK_QUADRATIC_TETRA] = new SMDS_DownQuadTetra(this);
443 _downArray[VTK_PYRAMID] = new SMDS_DownPyramid(this);
444 _downArray[VTK_QUADRATIC_PYRAMID] = new SMDS_DownQuadPyramid(this);
445 _downArray[VTK_WEDGE] = new SMDS_DownPenta(this);
446 _downArray[VTK_QUADRATIC_WEDGE] = new SMDS_DownQuadPenta(this);
447 _downArray[VTK_HEXAHEDRON] = new SMDS_DownHexa(this);
448 _downArray[VTK_QUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
449 _downArray[VTK_TRIQUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
450 _downArray[VTK_HEXAGONAL_PRISM] = new SMDS_DownPenta(this);
452 // --- get detailed info of number of cells of each type, allocate SMDS_downward structures
454 const SMDS_MeshInfo &meshInfo = _mesh->GetMeshInfo();
456 int nbLinTetra = meshInfo.NbTetras (ORDER_LINEAR);
457 int nbQuadTetra = meshInfo.NbTetras (ORDER_QUADRATIC);
458 int nbLinPyra = meshInfo.NbPyramids(ORDER_LINEAR);
459 int nbQuadPyra = meshInfo.NbPyramids(ORDER_QUADRATIC);
460 int nbLinPrism = meshInfo.NbPrisms (ORDER_LINEAR);
461 int nbQuadPrism = meshInfo.NbPrisms (ORDER_QUADRATIC);
462 int nbLinHexa = meshInfo.NbHexas (ORDER_LINEAR);
463 int nbQuadHexa = meshInfo.NbHexas (ORDER_QUADRATIC);
464 int nbHexPrism = meshInfo.NbHexPrisms();
466 int nbLineGuess = int((4.0 / 3.0) * nbLinTetra + 2 * nbLinPrism + 2.5 * nbLinPyra + 3 * nbLinHexa);
467 int nbQuadEdgeGuess = int((4.0 / 3.0) * nbQuadTetra + 2 * nbQuadPrism + 2.5 * nbQuadPyra + 3 * nbQuadHexa);
468 int nbLinTriaGuess = 2 * nbLinTetra + nbLinPrism + 2 * nbLinPyra;
469 int nbQuadTriaGuess = 2 * nbQuadTetra + nbQuadPrism + 2 * nbQuadPyra;
470 int nbLinQuadGuess = int((2.0 / 3.0) * nbLinPrism + (1.0 / 2.0) * nbLinPyra + 3 * nbLinHexa);
471 int nbQuadQuadGuess = int((2.0 / 3.0) * nbQuadPrism + (1.0 / 2.0) * nbQuadPyra + 3 * nbQuadHexa);
473 int GuessSize[VTK_MAXTYPE];
474 GuessSize[VTK_LINE] = nbLineGuess;
475 GuessSize[VTK_QUADRATIC_EDGE] = nbQuadEdgeGuess;
476 GuessSize[VTK_TRIANGLE] = nbLinTriaGuess;
477 GuessSize[VTK_QUADRATIC_TRIANGLE] = nbQuadTriaGuess;
478 GuessSize[VTK_BIQUADRATIC_TRIANGLE] = nbQuadTriaGuess;
479 GuessSize[VTK_QUAD] = nbLinQuadGuess;
480 GuessSize[VTK_QUADRATIC_QUAD] = nbQuadQuadGuess;
481 GuessSize[VTK_BIQUADRATIC_QUAD] = nbQuadQuadGuess;
482 GuessSize[VTK_TETRA] = nbLinTetra;
483 GuessSize[VTK_QUADRATIC_TETRA] = nbQuadTetra;
484 GuessSize[VTK_PYRAMID] = nbLinPyra;
485 GuessSize[VTK_QUADRATIC_PYRAMID] = nbQuadPyra;
486 GuessSize[VTK_WEDGE] = nbLinPrism;
487 GuessSize[VTK_QUADRATIC_WEDGE] = nbQuadPrism;
488 GuessSize[VTK_HEXAHEDRON] = nbLinHexa;
489 GuessSize[VTK_QUADRATIC_HEXAHEDRON] = nbQuadHexa;
490 GuessSize[VTK_TRIQUADRATIC_HEXAHEDRON] = nbQuadHexa;
491 GuessSize[VTK_HEXAGONAL_PRISM] = nbHexPrism;
493 _downArray[VTK_LINE] ->allocate(nbLineGuess);
494 _downArray[VTK_QUADRATIC_EDGE] ->allocate(nbQuadEdgeGuess);
495 _downArray[VTK_TRIANGLE] ->allocate(nbLinTriaGuess);
496 _downArray[VTK_QUADRATIC_TRIANGLE] ->allocate(nbQuadTriaGuess);
497 _downArray[VTK_BIQUADRATIC_TRIANGLE] ->allocate(nbQuadTriaGuess);
498 _downArray[VTK_QUAD] ->allocate(nbLinQuadGuess);
499 _downArray[VTK_QUADRATIC_QUAD] ->allocate(nbQuadQuadGuess);
500 _downArray[VTK_BIQUADRATIC_QUAD] ->allocate(nbQuadQuadGuess);
501 _downArray[VTK_TETRA] ->allocate(nbLinTetra);
502 _downArray[VTK_QUADRATIC_TETRA] ->allocate(nbQuadTetra);
503 _downArray[VTK_PYRAMID] ->allocate(nbLinPyra);
504 _downArray[VTK_QUADRATIC_PYRAMID] ->allocate(nbQuadPyra);
505 _downArray[VTK_WEDGE] ->allocate(nbLinPrism);
506 _downArray[VTK_QUADRATIC_WEDGE] ->allocate(nbQuadPrism);
507 _downArray[VTK_HEXAHEDRON] ->allocate(nbLinHexa);
508 _downArray[VTK_QUADRATIC_HEXAHEDRON] ->allocate(nbQuadHexa);
509 _downArray[VTK_TRIQUADRATIC_HEXAHEDRON]->allocate(nbQuadHexa);
510 _downArray[VTK_HEXAGONAL_PRISM] ->allocate(nbHexPrism);
512 // --- iteration on vtkUnstructuredGrid cells, only faces
513 // for each vtk face:
514 // create a downward face entry with its downward id.
515 // compute vtk volumes, create downward volumes entry.
516 // mark face in downward volumes
517 // mark volumes in downward face
519 MESSAGE("--- iteration on vtkUnstructuredGrid cells, only faces");CHRONO(20);
520 int cellSize = this->Types->GetNumberOfTuples();
521 _cellIdToDownId.resize(cellSize, -1);
523 for (int i = 0; i < cellSize; i++)
525 int vtkFaceType = this->GetCellType(i);
526 if (SMDS_Downward::getCellDimension(vtkFaceType) == 2)
529 //ASSERT(_downArray[vtkFaceType]);
530 int connFaceId = _downArray[vtkFaceType]->addCell(vtkFaceId);
531 SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
532 downFace->setTempNodes(connFaceId, vtkFaceId);
533 int vols[2] = { -1, -1 };
534 int nbVolumes = downFace->computeVolumeIds(vtkFaceId, vols);
535 //MESSAGE("nbVolumes="<< nbVolumes);
536 for (int ivol = 0; ivol < nbVolumes; ivol++)
538 int vtkVolId = vols[ivol];
539 int vtkVolType = this->GetCellType(vtkVolId);
540 //ASSERT(_downArray[vtkVolType]);
541 int connVolId = _downArray[vtkVolType]->addCell(vtkVolId);
542 _downArray[vtkVolType]->addDownCell(connVolId, connFaceId, vtkFaceType);
543 _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId, vtkVolType);
544 // MESSAGE("Face " << vtkFaceId << " belongs to volume " << vtkVolId);
549 // --- iteration on vtkUnstructuredGrid cells, only volumes
550 // for each vtk volume:
551 // create downward volumes entry if not already done
552 // build a temporary list of faces described with their nodes
554 // compute the vtk volumes containing this face
555 // check if the face is already listed in the volumes (comparison of ordered list of nodes)
556 // if not, create a downward face entry (resizing of structure required)
557 // (the downward faces store a temporary list of nodes to ease the comparison)
558 // create downward volumes entry if not already done
559 // mark volumes in downward face
560 // mark face in downward volumes
563 MESSAGE("--- iteration on vtkUnstructuredGrid cells, only volumes");CHRONO(21);
565 for (int i = 0; i < cellSize; i++)
567 int vtkType = this->GetCellType(i);
568 if (SMDS_Downward::getCellDimension(vtkType) == 3)
572 // MESSAGE("vtk volume " << vtkVolId);
573 //ASSERT(_downArray[vtkType]);
574 /*int connVolId = */_downArray[vtkType]->addCell(vtkVolId);
576 // --- find all the faces of the volume, describe the faces by their nodes
578 SMDS_Down3D* downVol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
579 ListElemByNodesType facesWithNodes;
580 downVol->computeFacesWithNodes(vtkVolId, facesWithNodes);
581 // MESSAGE("vtk volume " << vtkVolId << " contains " << facesWithNodes.nbElems << " faces");
583 for (int iface = 0; iface < facesWithNodes.nbElems; iface++)
585 // --- find the volumes containing the face
588 int vtkFaceType = facesWithNodes.elems[iface].vtkType;
589 SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
590 int vols[2] = { -1, -1 };
591 int *nodes = &facesWithNodes.elems[iface].nodeIds[0];
592 int lg = facesWithNodes.elems[iface].nbNodes;
593 int nbVolumes = downFace->computeVolumeIdsFromNodesFace(nodes, lg, vols);
594 // MESSAGE("vtk volume " << vtkVolId << " face " << iface << " belongs to " << nbVolumes << " volumes");
596 // --- check if face is registered in the volumes
601 for (int ivol = 0; ivol < nbVolumes; ivol++)
603 int vtkVolId2 = vols[ivol];
604 int vtkVolType = this->GetCellType(vtkVolId2);
605 //ASSERT(_downArray[vtkVolType]);
606 int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
607 SMDS_Down3D* downVol2 = static_cast<SMDS_Down3D*> (_downArray[vtkVolType]);
608 connFaceId = downVol2->FindFaceByNodes(connVolId2, facesWithNodes.elems[iface]);
610 break; // --- face already created
613 // --- if face is not registered in the volumes, create face
618 connFaceId = _downArray[vtkFaceType]->addCell();
619 downFace->setTempNodes(connFaceId, facesWithNodes.elems[iface]);
622 // --- mark volumes in downward face and mark face in downward volumes
625 for (int ivol = 0; ivol < nbVolumes; ivol++)
627 int vtkVolId2 = vols[ivol];
628 int vtkVolType = this->GetCellType(vtkVolId2);
629 //ASSERT(_downArray[vtkVolType]);
630 int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
631 _downArray[vtkVolType]->addDownCell(connVolId2, connFaceId, vtkFaceType);
632 _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId2, vtkVolType);
633 // MESSAGE(" From volume " << vtkVolId << " face " << connFaceId << " belongs to volume " << vtkVolId2);
639 // --- iteration on vtkUnstructuredGrid cells, only edges
640 // for each vtk edge:
641 // create downward edge entry
642 // store the nodes id's in downward edge (redundant with vtkUnstructuredGrid)
643 // find downward faces
644 // (from vtk faces or volumes, get downward faces, they have a temporary list of nodes)
645 // mark edge in downward faces
646 // mark faces in downward edge
649 MESSAGE("--- iteration on vtkUnstructuredGrid cells, only edges");CHRONO(22);
651 for (int i = 0; i < cellSize; i++)
653 int vtkEdgeType = this->GetCellType(i);
654 if (SMDS_Downward::getCellDimension(vtkEdgeType) == 1)
657 //ASSERT(_downArray[vtkEdgeType]);
658 int connEdgeId = _downArray[vtkEdgeType]->addCell(vtkEdgeId);
659 SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
660 downEdge->setNodes(connEdgeId, vtkEdgeId);
662 int nbVtkCells = downEdge->computeVtkCells(connEdgeId, vtkIds);
664 unsigned char downTypes[1000];
665 int nbDownFaces = downEdge->computeFaces(connEdgeId, &vtkIds[0], nbVtkCells, downFaces, downTypes);
666 for (int n = 0; n < nbDownFaces; n++)
668 _downArray[downTypes[n]]->addDownCell(downFaces[n], connEdgeId, vtkEdgeType);
669 _downArray[vtkEdgeType]->addUpCell(connEdgeId, downFaces[n], downTypes[n]);
674 // --- iteration on downward faces (they are all listed now)
675 // for each downward face:
676 // build a temporary list of edges with their ordered list of nodes
678 // find all the vtk cells containing this edge
679 // then identify all the downward faces containing the edge, from the vtk cells
680 // check if the edge is already listed in the faces (comparison of ordered list of nodes)
681 // if not, create a downward edge entry with the node id's
682 // mark edge in downward faces
683 // mark downward faces in edge (size of list unknown, to be allocated)
685 CHRONOSTOP(22);CHRONO(23);
687 for (int vtkFaceType = 0; vtkFaceType < VTK_QUADRATIC_PYRAMID; vtkFaceType++)
689 if (SMDS_Downward::getCellDimension(vtkFaceType) != 2)
692 // --- find all the edges of the face, describe the edges by their nodes
694 SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
695 int maxId = downFace->getMaxId();
696 for (int faceId = 0; faceId < maxId; faceId++)
699 ListElemByNodesType edgesWithNodes;
700 downFace->computeEdgesWithNodes(faceId, edgesWithNodes);
701 // MESSAGE("downward face type " << vtkFaceType << " num " << faceId << " contains " << edgesWithNodes.nbElems << " edges");
704 for (int iedge = 0; iedge < edgesWithNodes.nbElems; iedge++)
707 // --- check if the edge is already registered by exploration of the faces
711 unsigned char vtkEdgeType = edgesWithNodes.elems[iedge].vtkType;
712 int *pts = &edgesWithNodes.elems[iedge].nodeIds[0];
713 SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
714 int nbVtkCells = downEdge->computeVtkCells(pts, vtkIds);
715 //CHRONOSTOP(41);CHRONO(42);
717 unsigned char downTypes[1000];
718 int nbDownFaces = downEdge->computeFaces(pts, &vtkIds[0], nbVtkCells, downFaces, downTypes);
723 for (int idf = 0; idf < nbDownFaces; idf++)
725 int faceId2 = downFaces[idf];
726 int faceType = downTypes[idf];
727 //ASSERT(_downArray[faceType]);
728 SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
729 connEdgeId = downFace2->FindEdgeByNodes(faceId2, edgesWithNodes.elems[iedge]);
731 break; // --- edge already created
734 // --- if edge is not registered in the faces, create edge
739 connEdgeId = _downArray[vtkEdgeType]->addCell();
740 downEdge->setNodes(connEdgeId, edgesWithNodes.elems[iedge].nodeIds);
744 // --- mark faces in downward edge and mark edge in downward faces
747 for (int idf = 0; idf < nbDownFaces; idf++)
749 int faceId2 = downFaces[idf];
750 int faceType = downTypes[idf];
751 //ASSERT(_downArray[faceType]);
752 //SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
753 _downArray[vtkEdgeType]->addUpCell(connEdgeId, faceId2, faceType);
754 _downArray[faceType]->addDownCell(faceId2, connEdgeId, vtkEdgeType);
755 // MESSAGE(" From face t:" << vtkFaceType << " " << faceId <<
756 // " edge " << connEdgeId << " belongs to face t:" << faceType << " " << faceId2);
762 CHRONOSTOP(23);CHRONO(24);
764 // compact downward connectivity structure: adjust downward arrays size, replace vector<vector int>> by a single vector<int>
765 // 3D first then 2D and last 1D to release memory before edge upCells reorganization, (temporary memory use)
767 for (int vtkType = VTK_QUADRATIC_PYRAMID; vtkType >= 0; vtkType--)
769 if (SMDS_Downward *down = _downArray[vtkType])
771 down->compactStorage();
777 for (int vtkType = 0; vtkType <= VTK_QUADRATIC_PYRAMID; vtkType++)
779 if (SMDS_Downward *down = _downArray[vtkType])
781 if (down->getMaxId())
783 MESSAGE("Cells of Type " << vtkType << " : number of entities, est: "
784 << GuessSize[vtkType] << " real: " << down->getMaxId());
787 }CHRONOSTOP(24);CHRONOSTOP(2);
791 /*! Get the neighbors of a cell.
792 * Only the neighbors having the dimension of the cell are taken into account
793 * (neighbors of a volume are the volumes sharing a face with this volume,
794 * neighbors of a face are the faces sharing an edge with this face...).
795 * @param neighborsVtkIds vector of neighbors vtk id's to fill (reserve enough space).
796 * @param downIds downward id's of cells of dimension n-1, to fill (reserve enough space).
797 * @param downTypes vtk types of cells of dimension n-1, to fill (reserve enough space).
798 * @param vtkId the vtk id of the cell
799 * @return number of neighbors
801 int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId, bool getSkin)
803 int vtkType = this->GetCellType(vtkId);
804 int cellDim = SMDS_Downward::getCellDimension(vtkType);
806 return 0; // TODO voisins des edges = edges connectees
807 int cellId = this->_cellIdToDownId[vtkId];
809 int nbCells = _downArray[vtkType]->getNumberOfDownCells(cellId);
810 const int *downCells = _downArray[vtkType]->getDownCells(cellId);
811 const unsigned char* downTyp = _downArray[vtkType]->getDownTypes(cellId);
813 // --- iteration on faces of the 3D cell (or edges on the 2D cell).
816 for (int i = 0; i < nbCells; i++)
818 int downId = downCells[i];
819 int cellType = downTyp[i];
820 int nbUp = _downArray[cellType]->getNumberOfUpCells(downId);
821 const int *upCells = _downArray[cellType]->getUpCells(downId);
822 const unsigned char* upTypes = _downArray[cellType]->getUpTypes(downId);
824 // ---for a volume, max 2 upCells, one is this cell, the other is a neighbor
825 // for a face, number of neighbors (connected faces) not known
827 for (int j = 0; j < nbUp; j++)
829 if ((upCells[j] == cellId) && (upTypes[j] == vtkType))
831 int vtkNeighbor = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
832 neighborsVtkIds[nb] = vtkNeighbor;
833 downIds[nb] = downId;
834 downTypes[nb] = cellType;
836 if (nb >= NBMAXNEIGHBORS)
838 INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
844 if (cellDim == 3 && nbUp == 1) // this face is on the skin of the volume
846 neighborsVtkIds[nb] = _downArray[cellType]->getVtkCellId(downId); // OK if skin present
847 downIds[nb] = downId;
848 downTypes[nb] = cellType;
850 if (nb >= NBMAXNEIGHBORS)
852 INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
861 /*! get the volumes containing a face or an edge of the grid
862 * The edge or face belongs to the vtkUnstructuredGrid
863 * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
864 * @param vtkId vtk id of the face or edge
866 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int vtkId)
868 int vtkType = this->GetCellType(vtkId);
869 int dim = SMDS_Downward::getCellDimension(vtkType);
871 unsigned char cellTypes[1000];
872 int downCellId[1000];
875 int downId = this->CellIdToDownId(vtkId);
878 MESSAGE("Downward structure not up to date: new edge not taken into account");
881 nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
882 const int *upCells = _downArray[vtkType]->getUpCells(downId);
883 const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
884 for (int i=0; i< nbFaces; i++)
886 cellTypes[i] = upTypes[i];
887 downCellId[i] = upCells[i];
893 cellTypes[0] = this->GetCellType(vtkId);
894 int downId = this->CellIdToDownId(vtkId);
897 MESSAGE("Downward structure not up to date: new face not taken into account");
900 downCellId[0] = downId;
904 for (int i=0; i<nbFaces; i++)
906 int vtkTypeFace = cellTypes[i];
907 int downId = downCellId[i];
908 int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
909 const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
910 const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
911 for (int j=0; j<nv; j++)
913 int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
915 volVtkIds[nbvol++] = vtkVolId;
921 /*! get the volumes containing a face or an edge of the downward structure
922 * The edge or face does not necessary belong to the vtkUnstructuredGrid
923 * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
924 * @param downId id in the downward structure
925 * @param downType type of cell
927 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int downId, unsigned char downType)
929 int vtkType = downType;
930 int dim = SMDS_Downward::getCellDimension(vtkType);
932 unsigned char cellTypes[1000];
933 int downCellId[1000];
936 nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
937 const int *upCells = _downArray[vtkType]->getUpCells(downId);
938 const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
939 for (int i=0; i< nbFaces; i++)
941 cellTypes[i] = upTypes[i];
942 downCellId[i] = upCells[i];
948 cellTypes[0] = vtkType;
949 downCellId[0] = downId;
953 for (int i=0; i<nbFaces; i++)
955 int vtkTypeFace = cellTypes[i];
956 int downId = downCellId[i];
957 int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
958 const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
959 const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
960 for (int j=0; j<nv; j++)
962 int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
964 volVtkIds[nbvol++] = vtkVolId;
970 /*! get the node id's of a cell.
971 * The cell is defined by it's downward connectivity id and type.
972 * @param nodeSet set of of vtk node id's to fill.
973 * @param downId downward connectivity id of the cell.
974 * @param downType type of cell.
976 void SMDS_UnstructuredGrid::GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType)
978 _downArray[downType]->getNodeIds(downId, nodeSet);
981 /*! change some nodes in cell without modifying type or internal connectivity.
982 * Nodes inverse connectivity is maintained up to date.
983 * @param vtkVolId vtk id of the cell
984 * @param localClonedNodeIds map old node id to new node id.
986 void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds)
989 vtkIdType *pts; // will refer to the point id's of the face
990 this->GetCellPoints(vtkVolId, npts, pts);
991 for (int i = 0; i < npts; i++)
993 if (localClonedNodeIds.count(pts[i]))
995 vtkIdType oldpt = pts[i];
996 pts[i] = localClonedNodeIds[oldpt];
997 //MESSAGE(oldpt << " --> " << pts[i]);
998 //this->RemoveReferenceToCell(oldpt, vtkVolId);
999 //this->AddReferenceToCell(pts[i], vtkVolId);
1004 /*! reorder the nodes of a face
1005 * @param vtkVolId vtk id of a volume containing the face, to get an orientation for the face.
1006 * @param orderedNodes list of nodes to reorder (in out)
1007 * @return size of the list
1009 int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes)
1011 int vtkType = this->GetCellType(vtkVolId);
1012 dim = SMDS_Downward::getCellDimension(vtkType);
1015 SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
1016 int downVolId = this->_cellIdToDownId[vtkVolId];
1017 downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
1019 // else nothing to do;
1020 return orderedNodes.size();
1023 void SMDS_UnstructuredGrid::BuildLinks()
1025 // Remove the old links if they are already built
1028 this->Links->UnRegister(this);
1031 SMDS_CellLinks* links;
1032 this->Links = links = SMDS_CellLinks::New();
1033 this->Links->Allocate(this->GetNumberOfPoints());
1034 this->Links->Register(this);
1035 links->BuildLinks(this, this->Connectivity,this->GetCellTypesArray() );
1036 this->Links->Delete();
1039 void SMDS_UnstructuredGrid::DeleteLinks()
1041 // Remove the old links if they are already built
1044 this->Links->UnRegister(this);
1048 SMDS_CellLinks* SMDS_UnstructuredGrid::GetLinks()
1052 return static_cast< SMDS_CellLinks* >( this->Links );
1055 /*! Create a volume (prism or hexahedron) by duplication of a face.
1056 * Designed for use in creation of flat elements separating volume domains.
1057 * A face separating two domains is shared by two volume cells.
1058 * All the nodes are already created (for the two faces).
1059 * Each original Node is associated to corresponding nodes in the domains.
1060 * Some nodes may be duplicated for more than two domains, when domain separations intersect.
1061 * In that case, even some of the nodes to use for the original face may be changed.
1062 * @param vtkVolId: vtk id of a volume containing the face, to get an orientation for the face.
1063 * @param domain1: domain of the original face
1064 * @param domain2: domain of the duplicated face
1065 * @param originalNodes: the vtk node ids of the original face
1066 * @param nodeDomains: map(original id --> map(domain --> duplicated node id))
1067 * @return ok if success.
1070 SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
1073 std::set<int>& originalNodes,
1074 std::map<int, std::map<int, int> >& nodeDomains,
1075 std::map<int, std::map<long, int> >& nodeQuadDomains)
1077 //MESSAGE("extrudeVolumeFromFace " << vtkVolId);
1078 vector<vtkIdType> orderedOriginals;
1079 orderedOriginals.clear();
1080 set<int>::const_iterator it = originalNodes.begin();
1081 for (; it != originalNodes.end(); ++it)
1082 orderedOriginals.push_back(*it);
1085 int nbNodes = this->getOrderedNodesOfFace(vtkVolId, dim, orderedOriginals);
1086 vector<vtkIdType> orderedNodes;
1088 bool isQuadratic = false;
1089 switch (orderedOriginals.size())
1100 isQuadratic = false;
1106 long dom1 = domain1;
1107 long dom2 = domain2;
1108 long dom1_2; // for nodeQuadDomains
1109 if (domain1 < domain2)
1110 dom1_2 = dom1 + INT_MAX * dom2;
1112 dom1_2 = dom2 + INT_MAX * dom1;
1113 //cerr << "dom1=" << dom1 << " dom2=" << dom2 << " dom1_2=" << dom1_2 << endl;
1114 int ima = orderedOriginals.size();
1115 int mid = orderedOriginals.size() / 2;
1116 //cerr << "ima=" << ima << " mid=" << mid << endl;
1117 for (int i = 0; i < mid; i++)
1118 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1119 for (int i = 0; i < mid; i++)
1120 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1121 for (int i = mid; i < ima; i++)
1122 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1123 for (int i = mid; i < ima; i++)
1124 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1125 for (int i = 0; i < mid; i++)
1127 int oldId = orderedOriginals[i];
1129 if (nodeQuadDomains.count(oldId) && nodeQuadDomains[oldId].count(dom1_2))
1130 newId = nodeQuadDomains[oldId][dom1_2];
1133 double *coords = this->GetPoint(oldId);
1134 SMDS_MeshNode *newNode = _mesh->AddNode(coords[0], coords[1], coords[2]);
1135 newId = newNode->getVtkId();
1136 if (! nodeQuadDomains.count(oldId))
1138 std::map<long, int> emptyMap;
1139 nodeQuadDomains[oldId] = emptyMap;
1141 nodeQuadDomains[oldId][dom1_2] = newId;
1143 orderedNodes.push_back(newId);
1148 for (int i = 0; i < nbNodes; i++)
1149 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1151 for (int i = 0; i < nbNodes; i++)
1152 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1154 for (int i = nbNodes-1; i >= 0; i--)
1155 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1161 SMDS_MeshVolume *vol = _mesh->AddVolumeFromVtkIds(orderedNodes);
1166 SMDS_MeshFace *face = _mesh->AddFaceFromVtkIds(orderedNodes);
1170 // TODO update sub-shape list of elements and nodes
1174 //================================================================================
1176 * \brief Allocates data array for ball diameters
1177 * \param MaxVtkID - max ID of a ball element
1179 //================================================================================
1181 void SMDS_UnstructuredGrid::AllocateDiameters( vtkIdType MaxVtkID )
1183 SetBallDiameter( MaxVtkID, 0 );
1186 //================================================================================
1188 * \brief Sets diameter of a ball element
1189 * \param vtkID - vtk id of the ball element
1190 * \param diameter - diameter of the ball element
1192 //================================================================================
1194 void SMDS_UnstructuredGrid::SetBallDiameter( vtkIdType vtkID, double diameter )
1196 vtkDoubleArray* array = vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() );
1199 array = vtkDoubleArray::New();
1200 array->SetNumberOfComponents(1);
1201 vtkDataSet::CellData->SetScalars( array );
1203 array->InsertValue( vtkID, diameter );
1206 //================================================================================
1208 * \brief Returns diameter of a ball element
1209 * \param vtkID - vtk id of the ball element
1211 //================================================================================
1213 double SMDS_UnstructuredGrid::GetBallDiameter( vtkIdType vtkID ) const
1215 if ( vtkDataSet::CellData )
1216 return vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )->GetValue( vtkID );