1 // Copyright (C) 2010-2019 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>
39 SMDS_CellLinks* SMDS_CellLinks::New()
41 return new SMDS_CellLinks();
44 void SMDS_CellLinks::ResizeForPoint(vtkIdType vtkID)
46 if ( vtkID > this->MaxId )
49 if ( vtkID >= this->Size )
50 vtkCellLinks::Resize( vtkID+SMDS_Mesh::chunkSize );
54 void SMDS_CellLinks::BuildLinks(vtkDataSet *data, vtkCellArray *Connectivity, vtkUnsignedCharArray* types)
56 // build links taking into account removed cells
58 vtkIdType numPts = data->GetNumberOfPoints();
59 vtkIdType j, cellId = 0;
60 unsigned short *linkLoc;
62 vtkIdType const *pts(nullptr);
63 vtkIdType loc = Connectivity->GetTraversalLocation();
65 // traverse data to determine number of uses of each point
67 for (Connectivity->InitTraversal();
68 Connectivity->GetNextCell(npts,pts); cellId++)
70 if ( types->GetValue( cellId ) != VTK_EMPTY_CELL )
71 for (j=0; j < npts; j++)
73 this->IncrementLinkCount(pts[j]);
77 // now allocate storage for the links
78 this->AllocateLinks(numPts);
79 this->MaxId = numPts - 1;
81 // fill out lists with references to cells
82 linkLoc = new unsigned short[numPts];
83 memset(linkLoc, 0, numPts*sizeof(unsigned short));
86 for (Connectivity->InitTraversal();
87 Connectivity->GetNextCell(npts,pts); cellId++)
89 if ( types->GetValue( cellId ) != VTK_EMPTY_CELL )
90 for (j=0; j < npts; j++)
92 this->InsertCellReference(pts[j], (linkLoc[pts[j]])++, cellId);
96 Connectivity->SetTraversalLocation(loc);
99 SMDS_CellLinks::SMDS_CellLinks() :
104 SMDS_CellLinks::~SMDS_CellLinks()
108 SMDS_UnstructuredGrid* SMDS_UnstructuredGrid::New()
110 return new SMDS_UnstructuredGrid();
113 SMDS_UnstructuredGrid::SMDS_UnstructuredGrid() :
114 vtkUnstructuredGrid()
116 _cellIdToDownId.clear();
122 SMDS_UnstructuredGrid::~SMDS_UnstructuredGrid()
126 vtkMTimeType SMDS_UnstructuredGrid::GetMTime()
128 vtkMTimeType mtime = vtkUnstructuredGrid::GetMTime();
132 vtkPoints* SMDS_UnstructuredGrid::GetPoints()
134 // TODO erreur incomprehensible de la macro vtk GetPoints apparue avec la version paraview de fin aout 2010
138 int SMDS_UnstructuredGrid::InsertNextLinkedCell(int type, int npts, vtkIdType *pts)
140 if ( !this->Links ) // don't create Links until they are needed
142 return this->InsertNextCell(type, npts, 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 std::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 if ( setOfNodes.insert( pts[i] ).second )
163 (static_cast< vtkCellLinks * >(this->Links.GetPointer()))->ResizeCellList( pts[i], 1 );
164 (static_cast< vtkCellLinks * >(this->Links.GetPointer()))->AddCellReference( cellid, pts[i] );
173 void SMDS_UnstructuredGrid::setSMDS_mesh(SMDS_Mesh *mesh)
178 void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize,
179 std::vector<int>& idCellsNewToOld, int newCellSize)
183 // IDs of VTK nodes always correspond to SMDS IDs but there can be "holes" in SMDS numeration.
184 // We compact only if there were holes
186 int oldNodeSize = this->GetNumberOfPoints();
187 bool updateNodes = ( oldNodeSize > newNodeSize );
188 if ( true /*updateNodes*/ )
190 // 21125: EDF 1233 SMESH: Degradation of precision in a test case for quadratic conversion
191 // Use double type for storing coordinates of nodes instead float.
192 vtkPoints *newPoints = vtkPoints::New();
193 newPoints->SetDataType( VTK_DOUBLE );
194 newPoints->SetNumberOfPoints( newNodeSize );
196 int i = 0, alreadyCopied = 0;
197 while ( i < oldNodeSize )
199 // skip a hole if any
200 while ( i < oldNodeSize && idNodesOldToNew[i] < 0 )
203 // look for a block end
204 while ( i < oldNodeSize && idNodesOldToNew[i] >= 0 )
207 copyNodes(newPoints, idNodesOldToNew, alreadyCopied, startBloc, endBloc);
209 this->SetPoints(newPoints);
214 this->Points->Squeeze();
215 this->Points->Modified();
218 // Compact cells if VTK IDs do not correspond to SMDS IDs or nodes compacted
220 int oldCellSize = this->Types->GetNumberOfTuples();
221 bool updateCells = ( updateNodes || newCellSize != oldCellSize );
222 for ( int newID = 0, nbIDs = idCellsNewToOld.size(); newID < nbIDs && !updateCells; ++newID )
223 updateCells = ( idCellsNewToOld[ newID ] != newID );
225 if ( false /*!updateCells*/ ) // no holes in elements
227 this->Connectivity->Squeeze();
228 this->CellLocations->Squeeze();
229 this->Types->Squeeze();
230 if ( this->FaceLocations )
232 this->FaceLocations->Squeeze();
233 this->Faces->Squeeze();
235 this->Connectivity->Modified();
239 if ((int) idNodesOldToNew.size() < oldNodeSize )
241 idNodesOldToNew.reserve( oldNodeSize );
242 for ( int i = idNodesOldToNew.size(); i < oldNodeSize; ++i )
243 idNodesOldToNew.push_back( i );
246 // --- create new compacted Connectivity, Locations and Types
248 int newConnectivitySize = this->Connectivity->GetNumberOfConnectivityEntries();
249 if ( newCellSize != oldCellSize )
250 for ( int i = 0; i < oldCellSize - 1; ++i )
251 if ( this->Types->GetValue( i ) == VTK_EMPTY_CELL )
252 newConnectivitySize -= this->CellLocations->GetValue( i+1 ) - this->CellLocations->GetValue( i );
254 vtkCellArray *newConnectivity = vtkCellArray::New();
255 newConnectivity->Initialize();
256 newConnectivity->Allocate( newConnectivitySize );
258 vtkUnsignedCharArray *newTypes = vtkUnsignedCharArray::New();
259 newTypes->Initialize();
260 newTypes->SetNumberOfValues(newCellSize);
262 vtkIdTypeArray *newLocations = vtkIdTypeArray::New();
263 newLocations->Initialize();
264 newLocations->SetNumberOfValues(newCellSize);
266 std::vector< vtkIdType > pointsCell(1024); // --- points id to fill a new cell
268 copyBloc(newTypes, idCellsNewToOld, idNodesOldToNew,
269 newConnectivity, newLocations, pointsCell );
271 if (vtkDoubleArray* diameters =
272 vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )) // Balls
274 vtkDoubleArray* newDiameters = vtkDoubleArray::New();
275 newDiameters->SetNumberOfComponents(1);
276 for ( int newCellID = 0; newCellID < newCellSize; newCellID++ )
278 if ( newTypes->GetValue( newCellID ) == VTK_POLY_VERTEX )
280 int oldCellID = idCellsNewToOld[ newCellID ];
281 newDiameters->InsertValue( newCellID, diameters->GetValue( oldCellID ));
283 vtkDataSet::CellData->SetScalars( newDiameters );
287 if (this->FaceLocations)
289 vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New();
290 newFaceLocations->Initialize();
291 newFaceLocations->Allocate(newTypes->GetSize());
292 vtkIdTypeArray *newFaces = vtkIdTypeArray::New();
293 newFaces->Initialize();
294 newFaces->Allocate(this->Faces->GetSize());
295 for ( int newCellID = 0; newCellID < newCellSize; newCellID++ )
297 if ( newTypes->GetValue( newCellID ) == VTK_POLYHEDRON )
299 int oldCellId = idCellsNewToOld[ newCellID ];
300 newFaceLocations->InsertNextValue( newFaces->GetMaxId()+1 );
301 int oldFaceLoc = this->FaceLocations->GetValue( oldCellId );
302 int nCellFaces = this->Faces->GetValue( oldFaceLoc++ );
303 newFaces->InsertNextValue( nCellFaces );
304 for ( int n = 0; n < nCellFaces; n++ )
306 int nptsInFace = this->Faces->GetValue( oldFaceLoc++ );
307 newFaces->InsertNextValue( nptsInFace );
308 for ( int k = 0; k < nptsInFace; k++ )
310 int oldpt = this->Faces->GetValue( oldFaceLoc++ );
311 newFaces->InsertNextValue( idNodesOldToNew[ oldpt ]);
317 newFaceLocations->InsertNextValue(-1);
320 newFaceLocations->Squeeze();
322 this->SetCells( newTypes, newLocations, newConnectivity, newFaceLocations, newFaces );
323 newFaceLocations->Delete();
328 this->SetCells( newTypes, newLocations, newConnectivity, FaceLocations, Faces );
332 newLocations->Delete();
333 newConnectivity->Delete();
336 void SMDS_UnstructuredGrid::copyNodes(vtkPoints * newPoints,
337 std::vector<int>& idNodesOldToNew,
342 void *target = newPoints->GetVoidPointer(3 * alreadyCopied);
343 void *source = this->Points->GetVoidPointer(3 * start);
344 int nbPoints = end - start;
347 memcpy(target, source, 3 * sizeof(double) * nbPoints);
348 alreadyCopied += nbPoints;
352 void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray * newTypes,
353 const std::vector<int>& idCellsNewToOld,
354 const std::vector<int>& idNodesOldToNew,
355 vtkCellArray* newConnectivity,
356 vtkIdTypeArray* newLocations,
357 std::vector<vtkIdType>& pointsCell)
359 for ( size_t iNew = 0; iNew < idCellsNewToOld.size(); iNew++ )
361 int iOld = idCellsNewToOld[ iNew ];
362 newTypes->SetValue( iNew, this->Types->GetValue( iOld ));
363 vtkIdType oldLoc = this->CellLocations->GetValue( iOld );
365 vtkIdType const *oldPtsCell(nullptr);
366 this->Connectivity->GetCell( oldLoc, nbpts, oldPtsCell );
367 if ((vtkIdType) pointsCell.size() < nbpts )
368 pointsCell.resize( nbpts );
369 for ( int l = 0; l < nbpts; l++ )
371 int oldval = oldPtsCell[l];
372 pointsCell[l] = idNodesOldToNew[oldval];
374 /*int newcnt = */newConnectivity->InsertNextCell( nbpts, pointsCell.data() );
375 int newLoc = newConnectivity->GetInsertLocation( nbpts );
376 newLocations->SetValue( iNew, newLoc );
380 int SMDS_UnstructuredGrid::CellIdToDownId(int vtkCellId)
382 if ((vtkCellId < 0) || (vtkCellId >= (int)_cellIdToDownId.size()))
386 return _cellIdToDownId[vtkCellId];
389 void SMDS_UnstructuredGrid::setCellIdToDownId(int vtkCellId, int downId)
391 // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
392 _cellIdToDownId[vtkCellId] = downId;
395 void SMDS_UnstructuredGrid::CleanDownwardConnectivity()
397 for (size_t i = 0; i < _downArray.size(); i++)
400 delete _downArray[i];
403 _cellIdToDownId.clear();
406 /*! Build downward connectivity: to do only when needed because heavy memory load.
407 * Downward connectivity is no more valid if vtkUnstructuredGrid is modified.
410 void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
412 MESSAGE("SMDS_UnstructuredGrid::BuildDownwardConnectivity");CHRONO(2);
413 // TODO calcul partiel sans edges
415 // --- erase previous data if any
417 this->CleanDownwardConnectivity();
419 // --- create SMDS_Downward structures (in _downArray vector[vtkCellType])
421 _downArray.resize(VTK_MAXTYPE + 1, 0);
423 _downArray[VTK_LINE] = new SMDS_DownEdge(this);
424 _downArray[VTK_QUADRATIC_EDGE] = new SMDS_DownQuadEdge(this);
425 _downArray[VTK_TRIANGLE] = new SMDS_DownTriangle(this);
426 _downArray[VTK_QUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
427 _downArray[VTK_BIQUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
428 _downArray[VTK_QUAD] = new SMDS_DownQuadrangle(this);
429 _downArray[VTK_QUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
430 _downArray[VTK_BIQUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
431 _downArray[VTK_TETRA] = new SMDS_DownTetra(this);
432 _downArray[VTK_QUADRATIC_TETRA] = new SMDS_DownQuadTetra(this);
433 _downArray[VTK_PYRAMID] = new SMDS_DownPyramid(this);
434 _downArray[VTK_QUADRATIC_PYRAMID] = new SMDS_DownQuadPyramid(this);
435 _downArray[VTK_WEDGE] = new SMDS_DownPenta(this);
436 _downArray[VTK_QUADRATIC_WEDGE] = new SMDS_DownQuadPenta(this);
437 _downArray[VTK_HEXAHEDRON] = new SMDS_DownHexa(this);
438 _downArray[VTK_QUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
439 _downArray[VTK_TRIQUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
440 _downArray[VTK_HEXAGONAL_PRISM] = new SMDS_DownPenta(this);
442 // --- get detailed info of number of cells of each type, allocate SMDS_downward structures
444 const SMDS_MeshInfo &meshInfo = _mesh->GetMeshInfo();
446 int nbLinTetra = meshInfo.NbTetras (ORDER_LINEAR);
447 int nbQuadTetra = meshInfo.NbTetras (ORDER_QUADRATIC);
448 int nbLinPyra = meshInfo.NbPyramids(ORDER_LINEAR);
449 int nbQuadPyra = meshInfo.NbPyramids(ORDER_QUADRATIC);
450 int nbLinPrism = meshInfo.NbPrisms (ORDER_LINEAR);
451 int nbQuadPrism = meshInfo.NbPrisms (ORDER_QUADRATIC);
452 int nbLinHexa = meshInfo.NbHexas (ORDER_LINEAR);
453 int nbQuadHexa = meshInfo.NbHexas (ORDER_QUADRATIC);
454 int nbHexPrism = meshInfo.NbHexPrisms();
456 int nbLineGuess = int((4.0 / 3.0) * nbLinTetra + 2 * nbLinPrism + 2.5 * nbLinPyra + 3 * nbLinHexa);
457 int nbQuadEdgeGuess = int((4.0 / 3.0) * nbQuadTetra + 2 * nbQuadPrism + 2.5 * nbQuadPyra + 3 * nbQuadHexa);
458 int nbLinTriaGuess = 2 * nbLinTetra + nbLinPrism + 2 * nbLinPyra;
459 int nbQuadTriaGuess = 2 * nbQuadTetra + nbQuadPrism + 2 * nbQuadPyra;
460 int nbLinQuadGuess = int((2.0 / 3.0) * nbLinPrism + (1.0 / 2.0) * nbLinPyra + 3 * nbLinHexa);
461 int nbQuadQuadGuess = int((2.0 / 3.0) * nbQuadPrism + (1.0 / 2.0) * nbQuadPyra + 3 * nbQuadHexa);
463 int GuessSize[VTK_MAXTYPE];
464 GuessSize[VTK_LINE] = nbLineGuess;
465 GuessSize[VTK_QUADRATIC_EDGE] = nbQuadEdgeGuess;
466 GuessSize[VTK_TRIANGLE] = nbLinTriaGuess;
467 GuessSize[VTK_QUADRATIC_TRIANGLE] = nbQuadTriaGuess;
468 GuessSize[VTK_BIQUADRATIC_TRIANGLE] = nbQuadTriaGuess;
469 GuessSize[VTK_QUAD] = nbLinQuadGuess;
470 GuessSize[VTK_QUADRATIC_QUAD] = nbQuadQuadGuess;
471 GuessSize[VTK_BIQUADRATIC_QUAD] = nbQuadQuadGuess;
472 GuessSize[VTK_TETRA] = nbLinTetra;
473 GuessSize[VTK_QUADRATIC_TETRA] = nbQuadTetra;
474 GuessSize[VTK_PYRAMID] = nbLinPyra;
475 GuessSize[VTK_QUADRATIC_PYRAMID] = nbQuadPyra;
476 GuessSize[VTK_WEDGE] = nbLinPrism;
477 GuessSize[VTK_QUADRATIC_WEDGE] = nbQuadPrism;
478 GuessSize[VTK_HEXAHEDRON] = nbLinHexa;
479 GuessSize[VTK_QUADRATIC_HEXAHEDRON] = nbQuadHexa;
480 GuessSize[VTK_TRIQUADRATIC_HEXAHEDRON] = nbQuadHexa;
481 GuessSize[VTK_HEXAGONAL_PRISM] = nbHexPrism;
483 _downArray[VTK_LINE] ->allocate(nbLineGuess);
484 _downArray[VTK_QUADRATIC_EDGE] ->allocate(nbQuadEdgeGuess);
485 _downArray[VTK_TRIANGLE] ->allocate(nbLinTriaGuess);
486 _downArray[VTK_QUADRATIC_TRIANGLE] ->allocate(nbQuadTriaGuess);
487 _downArray[VTK_BIQUADRATIC_TRIANGLE] ->allocate(nbQuadTriaGuess);
488 _downArray[VTK_QUAD] ->allocate(nbLinQuadGuess);
489 _downArray[VTK_QUADRATIC_QUAD] ->allocate(nbQuadQuadGuess);
490 _downArray[VTK_BIQUADRATIC_QUAD] ->allocate(nbQuadQuadGuess);
491 _downArray[VTK_TETRA] ->allocate(nbLinTetra);
492 _downArray[VTK_QUADRATIC_TETRA] ->allocate(nbQuadTetra);
493 _downArray[VTK_PYRAMID] ->allocate(nbLinPyra);
494 _downArray[VTK_QUADRATIC_PYRAMID] ->allocate(nbQuadPyra);
495 _downArray[VTK_WEDGE] ->allocate(nbLinPrism);
496 _downArray[VTK_QUADRATIC_WEDGE] ->allocate(nbQuadPrism);
497 _downArray[VTK_HEXAHEDRON] ->allocate(nbLinHexa);
498 _downArray[VTK_QUADRATIC_HEXAHEDRON] ->allocate(nbQuadHexa);
499 _downArray[VTK_TRIQUADRATIC_HEXAHEDRON]->allocate(nbQuadHexa);
500 _downArray[VTK_HEXAGONAL_PRISM] ->allocate(nbHexPrism);
502 // --- iteration on vtkUnstructuredGrid cells, only faces
503 // for each vtk face:
504 // create a downward face entry with its downward id.
505 // compute vtk volumes, create downward volumes entry.
506 // mark face in downward volumes
507 // mark volumes in downward face
509 MESSAGE("--- iteration on vtkUnstructuredGrid cells, only faces");CHRONO(20);
510 int cellSize = this->Types->GetNumberOfTuples();
511 _cellIdToDownId.resize(cellSize, -1);
513 for (int i = 0; i < cellSize; i++)
515 int vtkFaceType = this->GetCellType(i);
516 if (SMDS_Downward::getCellDimension(vtkFaceType) == 2)
519 //ASSERT(_downArray[vtkFaceType]);
520 int connFaceId = _downArray[vtkFaceType]->addCell(vtkFaceId);
521 SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
522 downFace->setTempNodes(connFaceId, vtkFaceId);
523 int vols[2] = { -1, -1 };
524 int nbVolumes = downFace->computeVolumeIds(vtkFaceId, vols);
525 //MESSAGE("nbVolumes="<< nbVolumes);
526 for (int ivol = 0; ivol < nbVolumes; ivol++)
528 int vtkVolId = vols[ivol];
529 int vtkVolType = this->GetCellType(vtkVolId);
530 //ASSERT(_downArray[vtkVolType]);
531 int connVolId = _downArray[vtkVolType]->addCell(vtkVolId);
532 _downArray[vtkVolType]->addDownCell(connVolId, connFaceId, vtkFaceType);
533 _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId, vtkVolType);
534 // MESSAGE("Face " << vtkFaceId << " belongs to volume " << vtkVolId);
539 // --- iteration on vtkUnstructuredGrid cells, only volumes
540 // for each vtk volume:
541 // create downward volumes entry if not already done
542 // build a temporary list of faces described with their nodes
544 // compute the vtk volumes containing this face
545 // check if the face is already listed in the volumes (comparison of ordered list of nodes)
546 // if not, create a downward face entry (resizing of structure required)
547 // (the downward faces store a temporary list of nodes to ease the comparison)
548 // create downward volumes entry if not already done
549 // mark volumes in downward face
550 // mark face in downward volumes
553 MESSAGE("--- iteration on vtkUnstructuredGrid cells, only volumes");CHRONO(21);
555 for (int i = 0; i < cellSize; i++)
557 int vtkType = this->GetCellType(i);
558 if (SMDS_Downward::getCellDimension(vtkType) == 3)
562 // MESSAGE("vtk volume " << vtkVolId);
563 //ASSERT(_downArray[vtkType]);
564 /*int connVolId = */_downArray[vtkType]->addCell(vtkVolId);
566 // --- find all the faces of the volume, describe the faces by their nodes
568 SMDS_Down3D* downVol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
569 ListElemByNodesType facesWithNodes;
570 downVol->computeFacesWithNodes(vtkVolId, facesWithNodes);
571 // MESSAGE("vtk volume " << vtkVolId << " contains " << facesWithNodes.nbElems << " faces");
573 for (int iface = 0; iface < facesWithNodes.nbElems; iface++)
575 // --- find the volumes containing the face
578 int vtkFaceType = facesWithNodes.elems[iface].vtkType;
579 SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
580 int vols[2] = { -1, -1 };
581 int *nodes = &facesWithNodes.elems[iface].nodeIds[0];
582 int lg = facesWithNodes.elems[iface].nbNodes;
583 int nbVolumes = downFace->computeVolumeIdsFromNodesFace(nodes, lg, vols);
584 // MESSAGE("vtk volume " << vtkVolId << " face " << iface << " belongs to " << nbVolumes << " volumes");
586 // --- check if face is registered in the volumes
591 for (int ivol = 0; ivol < nbVolumes; ivol++)
593 int vtkVolId2 = vols[ivol];
594 int vtkVolType = this->GetCellType(vtkVolId2);
595 //ASSERT(_downArray[vtkVolType]);
596 int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
597 SMDS_Down3D* downVol2 = static_cast<SMDS_Down3D*> (_downArray[vtkVolType]);
598 connFaceId = downVol2->FindFaceByNodes(connVolId2, facesWithNodes.elems[iface]);
600 break; // --- face already created
603 // --- if face is not registered in the volumes, create face
608 connFaceId = _downArray[vtkFaceType]->addCell();
609 downFace->setTempNodes(connFaceId, facesWithNodes.elems[iface]);
612 // --- mark volumes in downward face and mark face in downward volumes
615 for (int ivol = 0; ivol < nbVolumes; ivol++)
617 int vtkVolId2 = vols[ivol];
618 int vtkVolType = this->GetCellType(vtkVolId2);
619 //ASSERT(_downArray[vtkVolType]);
620 int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
621 _downArray[vtkVolType]->addDownCell(connVolId2, connFaceId, vtkFaceType);
622 _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId2, vtkVolType);
623 // MESSAGE(" From volume " << vtkVolId << " face " << connFaceId << " belongs to volume " << vtkVolId2);
629 // --- iteration on vtkUnstructuredGrid cells, only edges
630 // for each vtk edge:
631 // create downward edge entry
632 // store the nodes id's in downward edge (redundant with vtkUnstructuredGrid)
633 // find downward faces
634 // (from vtk faces or volumes, get downward faces, they have a temporary list of nodes)
635 // mark edge in downward faces
636 // mark faces in downward edge
639 MESSAGE("--- iteration on vtkUnstructuredGrid cells, only edges");CHRONO(22);
641 for (int i = 0; i < cellSize; i++)
643 int vtkEdgeType = this->GetCellType(i);
644 if (SMDS_Downward::getCellDimension(vtkEdgeType) == 1)
647 //ASSERT(_downArray[vtkEdgeType]);
648 int connEdgeId = _downArray[vtkEdgeType]->addCell(vtkEdgeId);
649 SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
650 downEdge->setNodes(connEdgeId, vtkEdgeId);
651 std::vector<int> vtkIds;
652 int nbVtkCells = downEdge->computeVtkCells(connEdgeId, vtkIds);
654 unsigned char downTypes[1000];
655 int nbDownFaces = downEdge->computeFaces(connEdgeId, &vtkIds[0], nbVtkCells, downFaces, downTypes);
656 for (int n = 0; n < nbDownFaces; n++)
658 _downArray[downTypes[n]]->addDownCell(downFaces[n], connEdgeId, vtkEdgeType);
659 _downArray[vtkEdgeType]->addUpCell(connEdgeId, downFaces[n], downTypes[n]);
664 // --- iteration on downward faces (they are all listed now)
665 // for each downward face:
666 // build a temporary list of edges with their ordered list of nodes
668 // find all the vtk cells containing this edge
669 // then identify all the downward faces containing the edge, from the vtk cells
670 // check if the edge is already listed in the faces (comparison of ordered list of nodes)
671 // if not, create a downward edge entry with the node id's
672 // mark edge in downward faces
673 // mark downward faces in edge (size of list unknown, to be allocated)
675 CHRONOSTOP(22);CHRONO(23);
677 for (int vtkFaceType = 0; vtkFaceType < VTK_QUADRATIC_PYRAMID; vtkFaceType++)
679 if (SMDS_Downward::getCellDimension(vtkFaceType) != 2)
682 // --- find all the edges of the face, describe the edges by their nodes
684 SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
685 int maxId = downFace->getMaxId();
686 for (int faceId = 0; faceId < maxId; faceId++)
689 ListElemByNodesType edgesWithNodes;
690 downFace->computeEdgesWithNodes(faceId, edgesWithNodes);
691 // MESSAGE("downward face type " << vtkFaceType << " num " << faceId << " contains " << edgesWithNodes.nbElems << " edges");
694 for (int iedge = 0; iedge < edgesWithNodes.nbElems; iedge++)
697 // --- check if the edge is already registered by exploration of the faces
700 std::vector<int> vtkIds;
701 unsigned char vtkEdgeType = edgesWithNodes.elems[iedge].vtkType;
702 int *pts = &edgesWithNodes.elems[iedge].nodeIds[0];
703 SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
704 int nbVtkCells = downEdge->computeVtkCells(pts, vtkIds);
705 //CHRONOSTOP(41);CHRONO(42);
707 unsigned char downTypes[1000];
708 int nbDownFaces = downEdge->computeFaces(pts, &vtkIds[0], nbVtkCells, downFaces, downTypes);
713 for (int idf = 0; idf < nbDownFaces; idf++)
715 int faceId2 = downFaces[idf];
716 int faceType = downTypes[idf];
717 //ASSERT(_downArray[faceType]);
718 SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
719 connEdgeId = downFace2->FindEdgeByNodes(faceId2, edgesWithNodes.elems[iedge]);
721 break; // --- edge already created
724 // --- if edge is not registered in the faces, create edge
729 connEdgeId = _downArray[vtkEdgeType]->addCell();
730 downEdge->setNodes(connEdgeId, edgesWithNodes.elems[iedge].nodeIds);
734 // --- mark faces in downward edge and mark edge in downward faces
737 for (int idf = 0; idf < nbDownFaces; idf++)
739 int faceId2 = downFaces[idf];
740 int faceType = downTypes[idf];
741 //ASSERT(_downArray[faceType]);
742 //SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
743 _downArray[vtkEdgeType]->addUpCell(connEdgeId, faceId2, faceType);
744 _downArray[faceType]->addDownCell(faceId2, connEdgeId, vtkEdgeType);
745 // MESSAGE(" From face t:" << vtkFaceType << " " << faceId <<
746 // " edge " << connEdgeId << " belongs to face t:" << faceType << " " << faceId2);
752 CHRONOSTOP(23);CHRONO(24);
754 // compact downward connectivity structure: adjust downward arrays size, replace std::vector<vector int>> by a single std::vector<int>
755 // 3D first then 2D and last 1D to release memory before edge upCells reorganization, (temporary memory use)
757 for (int vtkType = VTK_QUADRATIC_PYRAMID; vtkType >= 0; vtkType--)
759 if (SMDS_Downward *down = _downArray[vtkType])
761 down->compactStorage();
767 for (int vtkType = 0; vtkType <= VTK_QUADRATIC_PYRAMID; vtkType++)
769 if (SMDS_Downward *down = _downArray[vtkType])
771 if (down->getMaxId())
773 MESSAGE("Cells of Type " << vtkType << " : number of entities, est: "
774 << GuessSize[vtkType] << " real: " << down->getMaxId());
777 }CHRONOSTOP(24);CHRONOSTOP(2);
781 /*! Get the neighbors of a cell.
782 * Only the neighbors having the dimension of the cell are taken into account
783 * (neighbors of a volume are the volumes sharing a face with this volume,
784 * neighbors of a face are the faces sharing an edge with this face...).
785 * @param neighborsVtkIds vector of neighbors vtk id's to fill (reserve enough space).
786 * @param downIds downward id's of cells of dimension n-1, to fill (reserve enough space).
787 * @param downTypes vtk types of cells of dimension n-1, to fill (reserve enough space).
788 * @param vtkId the vtk id of the cell
789 * @return number of neighbors
791 int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId, bool getSkin)
793 int vtkType = this->GetCellType(vtkId);
794 int cellDim = SMDS_Downward::getCellDimension(vtkType);
796 return 0; // TODO voisins des edges = edges connectees
797 int cellId = this->_cellIdToDownId[vtkId];
799 int nbCells = _downArray[vtkType]->getNumberOfDownCells(cellId);
800 const int *downCells = _downArray[vtkType]->getDownCells(cellId);
801 const unsigned char* downTyp = _downArray[vtkType]->getDownTypes(cellId);
803 // --- iteration on faces of the 3D cell (or edges on the 2D cell).
806 for (int i = 0; i < nbCells; i++)
808 int downId = downCells[i];
809 int cellType = downTyp[i];
810 int nbUp = _downArray[cellType]->getNumberOfUpCells(downId);
811 const int *upCells = _downArray[cellType]->getUpCells(downId);
812 const unsigned char* upTypes = _downArray[cellType]->getUpTypes(downId);
814 // ---for a volume, max 2 upCells, one is this cell, the other is a neighbor
815 // for a face, number of neighbors (connected faces) not known
817 for (int j = 0; j < nbUp; j++)
819 if ((upCells[j] == cellId) && (upTypes[j] == vtkType))
821 int vtkNeighbor = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
822 neighborsVtkIds[nb] = vtkNeighbor;
823 downIds[nb] = downId;
824 downTypes[nb] = cellType;
826 if (nb >= NBMAXNEIGHBORS)
828 INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
834 if (cellDim == 3 && nbUp == 1) // this face is on the skin of the volume
836 neighborsVtkIds[nb] = _downArray[cellType]->getVtkCellId(downId); // OK if skin present
837 downIds[nb] = downId;
838 downTypes[nb] = cellType;
840 if (nb >= NBMAXNEIGHBORS)
842 INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
851 /*! get the volumes containing a face or an edge of the grid
852 * The edge or face belongs to the vtkUnstructuredGrid
853 * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
854 * @param vtkId vtk id of the face or edge
856 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int vtkId)
858 int vtkType = this->GetCellType(vtkId);
859 int dim = SMDS_Downward::getCellDimension(vtkType);
861 unsigned char cellTypes[1000];
862 int downCellId[1000];
865 int downId = this->CellIdToDownId(vtkId);
868 MESSAGE("Downward structure not up to date: new edge not taken into account");
871 nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
872 const int *upCells = _downArray[vtkType]->getUpCells(downId);
873 const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
874 for (int i=0; i< nbFaces; i++)
876 cellTypes[i] = upTypes[i];
877 downCellId[i] = upCells[i];
883 cellTypes[0] = this->GetCellType(vtkId);
884 int downId = this->CellIdToDownId(vtkId);
887 MESSAGE("Downward structure not up to date: new face not taken into account");
890 downCellId[0] = downId;
894 for (int i=0; i<nbFaces; i++)
896 int vtkTypeFace = cellTypes[i];
897 int downId = downCellId[i];
898 int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
899 const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
900 const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
901 for (int j=0; j<nv; j++)
903 int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
905 volVtkIds[nbvol++] = vtkVolId;
911 /*! get the volumes containing a face or an edge of the downward structure
912 * The edge or face does not necessary belong to the vtkUnstructuredGrid
913 * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
914 * @param downId id in the downward structure
915 * @param downType type of cell
917 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int downId, unsigned char downType)
919 int vtkType = downType;
920 int dim = SMDS_Downward::getCellDimension(vtkType);
922 unsigned char cellTypes[1000];
923 int downCellId[1000];
926 nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
927 const int *upCells = _downArray[vtkType]->getUpCells(downId);
928 const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
929 for (int i=0; i< nbFaces; i++)
931 cellTypes[i] = upTypes[i];
932 downCellId[i] = upCells[i];
938 cellTypes[0] = vtkType;
939 downCellId[0] = downId;
943 for (int i=0; i<nbFaces; i++)
945 int vtkTypeFace = cellTypes[i];
946 int downId = downCellId[i];
947 int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
948 const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
949 const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
950 for (int j=0; j<nv; j++)
952 int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
954 volVtkIds[nbvol++] = vtkVolId;
960 /*! get the node id's of a cell.
961 * The cell is defined by it's downward connectivity id and type.
962 * @param nodeSet set of of vtk node id's to fill.
963 * @param downId downward connectivity id of the cell.
964 * @param downType type of cell.
966 void SMDS_UnstructuredGrid::GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType)
968 _downArray[downType]->getNodeIds(downId, nodeSet);
971 /*! change some nodes in cell without modifying type or internal connectivity.
972 * Nodes inverse connectivity is maintained up to date.
973 * @param vtkVolId vtk id of the cell
974 * @param localClonedNodeIds map old node id to new node id.
976 void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds)
979 vtkIdType const *tmp(nullptr); // will refer to the point id's of the face
980 vtkIdType *pts; // will refer to the point id's of the face
981 this->GetCellPoints(vtkVolId, npts, tmp);
982 std::copy(tmp, tmp+npts, pts);
983 for (int i = 0; i < npts; i++)
985 if (localClonedNodeIds.count(pts[i]))
987 vtkIdType oldpt = pts[i];
988 pts[i] = localClonedNodeIds[oldpt];
989 //MESSAGE(oldpt << " --> " << pts[i]);
990 //this->RemoveReferenceToCell(oldpt, vtkVolId);
991 //this->AddReferenceToCell(pts[i], vtkVolId);
996 /*! reorder the nodes of a face
997 * @param vtkVolId vtk id of a volume containing the face, to get an orientation for the face.
998 * @param orderedNodes list of nodes to reorder (in out)
999 * @return size of the list
1001 int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes)
1003 int vtkType = this->GetCellType(vtkVolId);
1004 dim = SMDS_Downward::getCellDimension(vtkType);
1007 SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
1008 int downVolId = this->_cellIdToDownId[vtkVolId];
1009 downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
1011 // else nothing to do;
1012 return orderedNodes.size();
1015 void SMDS_UnstructuredGrid::BuildLinks()
1017 // Remove the old links if they are already built
1020 this->Links->UnRegister(this);
1023 SMDS_CellLinks* links;
1024 this->Links = links = SMDS_CellLinks::New();
1025 (static_cast< vtkCellLinks *>(this->Links.GetPointer()))->Allocate(this->GetNumberOfPoints());
1026 this->Links->Register(this);
1027 links->BuildLinks(this, this->Connectivity,this->GetCellTypesArray() );
1028 this->Links->Delete();
1031 void SMDS_UnstructuredGrid::DeleteLinks()
1033 // Remove the old links if they are already built
1036 this->Links->UnRegister(this);
1040 SMDS_CellLinks* SMDS_UnstructuredGrid::GetLinks()
1044 return static_cast< SMDS_CellLinks* >( this->Links.GetPointer() );
1047 /*! Create a volume (prism or hexahedron) by duplication of a face.
1048 * Designed for use in creation of flat elements separating volume domains.
1049 * A face separating two domains is shared by two volume cells.
1050 * All the nodes are already created (for the two faces).
1051 * Each original Node is associated to corresponding nodes in the domains.
1052 * Some nodes may be duplicated for more than two domains, when domain separations intersect.
1053 * In that case, even some of the nodes to use for the original face may be changed.
1054 * @param vtkVolId: vtk id of a volume containing the face, to get an orientation for the face.
1055 * @param domain1: domain of the original face
1056 * @param domain2: domain of the duplicated face
1057 * @param originalNodes: the vtk node ids of the original face
1058 * @param nodeDomains: map(original id --> map(domain --> duplicated node id))
1059 * @return ok if success.
1062 SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
1065 std::set<int>& originalNodes,
1066 std::map<int, std::map<int, int> >& nodeDomains,
1067 std::map<int, std::map<long, int> >& nodeQuadDomains)
1069 //MESSAGE("extrudeVolumeFromFace " << vtkVolId);
1070 std::vector<vtkIdType> orderedOriginals( originalNodes.begin(), originalNodes.end() );
1073 int nbNodes = this->getOrderedNodesOfFace(vtkVolId, dim, orderedOriginals);
1074 std::vector<vtkIdType> orderedNodes;
1076 bool isQuadratic = false;
1077 switch (orderedOriginals.size())
1088 isQuadratic = false;
1094 long dom1 = domain1;
1095 long dom2 = domain2;
1096 long dom1_2; // for nodeQuadDomains
1097 if (domain1 < domain2)
1098 dom1_2 = dom1 + INT_MAX * dom2;
1100 dom1_2 = dom2 + INT_MAX * dom1;
1101 //cerr << "dom1=" << dom1 << " dom2=" << dom2 << " dom1_2=" << dom1_2 << endl;
1102 int ima = orderedOriginals.size();
1103 int mid = orderedOriginals.size() / 2;
1104 //cerr << "ima=" << ima << " mid=" << mid << endl;
1105 for (int i = 0; i < mid; i++)
1106 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1107 for (int i = 0; i < mid; i++)
1108 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1109 for (int i = mid; i < ima; i++)
1110 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1111 for (int i = mid; i < ima; i++)
1112 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1113 for (int i = 0; i < mid; i++)
1115 int oldId = orderedOriginals[i];
1117 if (nodeQuadDomains.count(oldId) && nodeQuadDomains[oldId].count(dom1_2))
1118 newId = nodeQuadDomains[oldId][dom1_2];
1121 double *coords = this->GetPoint(oldId);
1122 SMDS_MeshNode *newNode = _mesh->AddNode(coords[0], coords[1], coords[2]);
1123 newId = newNode->GetVtkID();
1124 if (! nodeQuadDomains.count(oldId))
1126 std::map<long, int> emptyMap;
1127 nodeQuadDomains[oldId] = emptyMap;
1129 nodeQuadDomains[oldId][dom1_2] = newId;
1131 orderedNodes.push_back(newId);
1136 for (int i = 0; i < nbNodes; i++)
1137 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1139 for (int i = 0; i < nbNodes; i++)
1140 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1142 for (int i = nbNodes-1; i >= 0; i--)
1143 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1149 SMDS_MeshVolume *vol = _mesh->AddVolumeFromVtkIds(orderedNodes);
1154 SMDS_MeshFace *face = _mesh->AddFaceFromVtkIds(orderedNodes);
1158 // TODO update sub-shape list of elements and nodes
1162 //================================================================================
1164 * \brief Allocates data array for ball diameters
1165 * \param MaxVtkID - max ID of a ball element
1167 //================================================================================
1169 void SMDS_UnstructuredGrid::AllocateDiameters( vtkIdType MaxVtkID )
1171 SetBallDiameter( MaxVtkID, 0 );
1174 //================================================================================
1176 * \brief Sets diameter of a ball element
1177 * \param vtkID - vtk id of the ball element
1178 * \param diameter - diameter of the ball element
1180 //================================================================================
1182 void SMDS_UnstructuredGrid::SetBallDiameter( vtkIdType vtkID, double diameter )
1184 vtkDoubleArray* array = vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() );
1187 array = vtkDoubleArray::New();
1188 array->SetNumberOfComponents(1);
1189 vtkDataSet::CellData->SetScalars( array );
1191 array->InsertValue( vtkID, diameter );
1194 //================================================================================
1196 * \brief Returns diameter of a ball element
1197 * \param vtkID - vtk id of the ball element
1199 //================================================================================
1201 double SMDS_UnstructuredGrid::GetBallDiameter( vtkIdType vtkID ) const
1203 if ( vtkDataSet::CellData )
1204 return vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )->GetValue( vtkID );