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();
133 // OUV_PORTING_VTK6: seems to be useless
135 void SMDS_UnstructuredGrid::Update()
137 return vtkUnstructuredGrid::Update();
140 void SMDS_UnstructuredGrid::UpdateInformation()
142 return vtkUnstructuredGrid::UpdateInformation();
145 vtkPoints* SMDS_UnstructuredGrid::GetPoints()
147 // TODO erreur incomprehensible de la macro vtk GetPoints apparue avec la version paraview de fin aout 2010
151 //#ifdef VTK_HAVE_POLYHEDRON
152 int SMDS_UnstructuredGrid::InsertNextLinkedCell(int type, int npts, vtkIdType *pts)
154 if (type != VTK_POLYHEDRON)
155 return vtkUnstructuredGrid::InsertNextLinkedCell(type, npts, pts);
157 // --- type = VTK_POLYHEDRON
158 int cellid = this->InsertNextCell(type, npts, pts);
160 set<vtkIdType> setOfNodes;
164 for (int nf = 0; nf < nbfaces; nf++)
166 int nbnodes = pts[i];
168 for (int k = 0; k < nbnodes; k++)
170 setOfNodes.insert(pts[i]);
175 set<vtkIdType>::iterator it = setOfNodes.begin();
176 for (; it != setOfNodes.end(); ++it)
178 this->Links->ResizeCellList(*it, 1);
179 this->Links->AddCellReference(cellid, *it);
186 void SMDS_UnstructuredGrid::setSMDS_mesh(SMDS_Mesh *mesh)
191 void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize,
192 std::vector<int>& idCellsOldToNew, int newCellSize)
194 //MESSAGE("------------------------- compactGrid " << newNodeSize << " " << newCellSize);//CHRONO(1);
195 int alreadyCopied = 0;
199 // --- if newNodeSize, create a new compacted vtkPoints
201 vtkPoints *newPoints = vtkPoints::New();
202 newPoints->SetDataType(VTK_DOUBLE);
203 newPoints->SetNumberOfPoints(newNodeSize);
206 // rnv: to fix bug "21125: EDF 1233 SMESH: Degradation of precision in a test case for quadratic conversion"
207 // using double type for storing coordinates of nodes instead float.
208 int oldNodeSize = idNodesOldToNew.size();
211 while ( i < oldNodeSize )
213 // skip a hole if any
214 while ( i < oldNodeSize && idNodesOldToNew[i] < 0 )
217 // look for a block end
218 while ( i < oldNodeSize && idNodesOldToNew[i] >= 0 )
221 copyNodes(newPoints, idNodesOldToNew, alreadyCopied, startBloc, endBloc);
223 newPoints->Squeeze();
226 if (1/*newNodeSize*/)
228 this->SetPoints(newPoints);
233 // --- create new compacted Connectivity, Locations and Types
235 int oldCellSize = this->Types->GetNumberOfTuples();
237 if ( oldCellSize == newCellSize )
239 for ( int i = 0; i < oldCellSize; ++i )
240 idCellsOldToNew[i] = i;
243 vtkCellArray *newConnectivity = vtkCellArray::New();
244 newConnectivity->Initialize();
245 int oldCellDataSize = this->Connectivity->GetData()->GetSize();
246 newConnectivity->Allocate(oldCellDataSize);
248 vtkUnsignedCharArray *newTypes = vtkUnsignedCharArray::New();
249 newTypes->Initialize();
250 newTypes->SetNumberOfValues(newCellSize);
252 vtkIdTypeArray *newLocations = vtkIdTypeArray::New();
253 newLocations->Initialize();
254 newLocations->SetNumberOfValues(newCellSize);
256 // TODO some polyhedron may be huge (only in some tests)
257 vtkIdType tmpid[NBMAXNODESINCELL];
258 vtkIdType *pointsCell = &tmpid[0]; // --- points id to fill a new cell
262 while ( i < oldCellSize )
264 // skip a hole if any
265 while ( i < oldCellSize && this->Types->GetValue(i) == VTK_EMPTY_CELL )
268 // look for a block end
269 while ( i < oldCellSize && this->Types->GetValue(i) != VTK_EMPTY_CELL )
272 if ( endBloc > startBloc )
274 idCellsOldToNew, idNodesOldToNew,
275 newConnectivity, newLocations,
276 pointsCell, alreadyCopied,
279 newConnectivity->Squeeze();
281 if (vtkDoubleArray* diameters =
282 vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )) // Balls
284 for (int oldCellID = 0; oldCellID < oldCellSize; oldCellID++)
286 if (this->Types->GetValue(oldCellID) == VTK_EMPTY_CELL)
288 int newCellId = idCellsOldToNew[ oldCellID ];
289 if (newTypes->GetValue(newCellId) == VTK_POLY_VERTEX)
290 diameters->SetValue( newCellId, diameters->GetValue( oldCellID ));
294 if (this->FaceLocations)
296 vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New();
297 newFaceLocations->Initialize();
298 newFaceLocations->Allocate(newTypes->GetSize());
299 vtkIdTypeArray *newFaces = vtkIdTypeArray::New();
300 newFaces->Initialize();
301 newFaces->Allocate(this->Faces->GetSize());
302 for (int i = 0; i < oldCellSize; i++)
304 if (this->Types->GetValue(i) == VTK_EMPTY_CELL)
306 int newCellId = idCellsOldToNew[i];
307 if (newTypes->GetValue(newCellId) == VTK_POLYHEDRON)
309 newFaceLocations->InsertNextValue(newFaces->GetMaxId()+1);
310 int oldFaceLoc = this->FaceLocations->GetValue(i);
311 int nCellFaces = this->Faces->GetValue(oldFaceLoc++);
312 newFaces->InsertNextValue(nCellFaces);
313 for (int n=0; n<nCellFaces; n++)
315 int nptsInFace = this->Faces->GetValue(oldFaceLoc++);
316 newFaces->InsertNextValue(nptsInFace);
317 for (int k=0; k<nptsInFace; k++)
319 int oldpt = this->Faces->GetValue(oldFaceLoc++);
320 newFaces->InsertNextValue(idNodesOldToNew[oldpt]);
326 newFaceLocations->InsertNextValue(-1);
329 newFaceLocations->Squeeze();
331 this->SetCells(newTypes, newLocations, newConnectivity, newFaceLocations, newFaces);
332 newFaceLocations->Delete();
337 this->SetCells(newTypes, newLocations, newConnectivity, FaceLocations, Faces);
341 newLocations->Delete();
342 newConnectivity->Delete();
345 void SMDS_UnstructuredGrid::copyNodes(vtkPoints * newPoints,
346 std::vector<int>& idNodesOldToNew,
351 void *target = newPoints->GetVoidPointer(3 * alreadyCopied);
352 void *source = this->Points->GetVoidPointer(3 * start);
353 int nbPoints = end - start;
356 memcpy(target, source, 3 * sizeof(double) * nbPoints);
357 for (int j = start; j < end; j++)
358 idNodesOldToNew[j] = alreadyCopied++; // old vtkId --> new vtkId
362 void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes,
363 std::vector<int>& idCellsOldToNew,
364 std::vector<int>& idNodesOldToNew,
365 vtkCellArray* newConnectivity,
366 vtkIdTypeArray* newLocations,
367 vtkIdType* pointsCell,
372 for (int j = start; j < end; j++)
374 newTypes->SetValue(alreadyCopied, this->Types->GetValue(j));
375 idCellsOldToNew[j] = alreadyCopied; // old vtkId --> new vtkId
376 vtkIdType oldLoc = this->Locations->GetValue(j);
378 vtkIdType *oldPtsCell = 0;
379 this->Connectivity->GetCell(oldLoc, nbpts, oldPtsCell);
380 assert(nbpts < NBMAXNODESINCELL);
381 for (int l = 0; l < nbpts; l++)
383 int oldval = oldPtsCell[l];
384 pointsCell[l] = idNodesOldToNew[oldval];
386 /*int newcnt = */newConnectivity->InsertNextCell(nbpts, pointsCell);
387 int newLoc = newConnectivity->GetInsertLocation(nbpts);
388 newLocations->SetValue(alreadyCopied, newLoc);
393 int SMDS_UnstructuredGrid::CellIdToDownId(int vtkCellId)
395 if ((vtkCellId < 0) || (vtkCellId >= (int)_cellIdToDownId.size()))
399 return _cellIdToDownId[vtkCellId];
402 void SMDS_UnstructuredGrid::setCellIdToDownId(int vtkCellId, int downId)
404 // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
405 _cellIdToDownId[vtkCellId] = downId;
408 void SMDS_UnstructuredGrid::CleanDownwardConnectivity()
410 for (size_t i = 0; i < _downArray.size(); i++)
413 delete _downArray[i];
416 _cellIdToDownId.clear();
419 /*! Build downward connectivity: to do only when needed because heavy memory load.
420 * Downward connectivity is no more valid if vtkUnstructuredGrid is modified.
423 void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
425 MESSAGE("SMDS_UnstructuredGrid::BuildDownwardConnectivity");CHRONO(2);
426 // TODO calcul partiel sans edges
428 // --- erase previous data if any
430 this->CleanDownwardConnectivity();
432 // --- create SMDS_Downward structures (in _downArray vector[vtkCellType])
434 _downArray.resize(VTK_MAXTYPE + 1, 0);
436 _downArray[VTK_LINE] = new SMDS_DownEdge(this);
437 _downArray[VTK_QUADRATIC_EDGE] = new SMDS_DownQuadEdge(this);
438 _downArray[VTK_TRIANGLE] = new SMDS_DownTriangle(this);
439 _downArray[VTK_QUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
440 _downArray[VTK_BIQUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
441 _downArray[VTK_QUAD] = new SMDS_DownQuadrangle(this);
442 _downArray[VTK_QUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
443 _downArray[VTK_BIQUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
444 _downArray[VTK_TETRA] = new SMDS_DownTetra(this);
445 _downArray[VTK_QUADRATIC_TETRA] = new SMDS_DownQuadTetra(this);
446 _downArray[VTK_PYRAMID] = new SMDS_DownPyramid(this);
447 _downArray[VTK_QUADRATIC_PYRAMID] = new SMDS_DownQuadPyramid(this);
448 _downArray[VTK_WEDGE] = new SMDS_DownPenta(this);
449 _downArray[VTK_QUADRATIC_WEDGE] = new SMDS_DownQuadPenta(this);
450 _downArray[VTK_HEXAHEDRON] = new SMDS_DownHexa(this);
451 _downArray[VTK_QUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
452 _downArray[VTK_TRIQUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
453 _downArray[VTK_HEXAGONAL_PRISM] = new SMDS_DownPenta(this);
455 // --- get detailed info of number of cells of each type, allocate SMDS_downward structures
457 const SMDS_MeshInfo &meshInfo = _mesh->GetMeshInfo();
459 int nbLinTetra = meshInfo.NbTetras (ORDER_LINEAR);
460 int nbQuadTetra = meshInfo.NbTetras (ORDER_QUADRATIC);
461 int nbLinPyra = meshInfo.NbPyramids(ORDER_LINEAR);
462 int nbQuadPyra = meshInfo.NbPyramids(ORDER_QUADRATIC);
463 int nbLinPrism = meshInfo.NbPrisms (ORDER_LINEAR);
464 int nbQuadPrism = meshInfo.NbPrisms (ORDER_QUADRATIC);
465 int nbLinHexa = meshInfo.NbHexas (ORDER_LINEAR);
466 int nbQuadHexa = meshInfo.NbHexas (ORDER_QUADRATIC);
467 int nbHexPrism = meshInfo.NbHexPrisms();
469 int nbLineGuess = int((4.0 / 3.0) * nbLinTetra + 2 * nbLinPrism + 2.5 * nbLinPyra + 3 * nbLinHexa);
470 int nbQuadEdgeGuess = int((4.0 / 3.0) * nbQuadTetra + 2 * nbQuadPrism + 2.5 * nbQuadPyra + 3 * nbQuadHexa);
471 int nbLinTriaGuess = 2 * nbLinTetra + nbLinPrism + 2 * nbLinPyra;
472 int nbQuadTriaGuess = 2 * nbQuadTetra + nbQuadPrism + 2 * nbQuadPyra;
473 int nbLinQuadGuess = int((2.0 / 3.0) * nbLinPrism + (1.0 / 2.0) * nbLinPyra + 3 * nbLinHexa);
474 int nbQuadQuadGuess = int((2.0 / 3.0) * nbQuadPrism + (1.0 / 2.0) * nbQuadPyra + 3 * nbQuadHexa);
476 int GuessSize[VTK_MAXTYPE];
477 GuessSize[VTK_LINE] = nbLineGuess;
478 GuessSize[VTK_QUADRATIC_EDGE] = nbQuadEdgeGuess;
479 GuessSize[VTK_TRIANGLE] = nbLinTriaGuess;
480 GuessSize[VTK_QUADRATIC_TRIANGLE] = nbQuadTriaGuess;
481 GuessSize[VTK_BIQUADRATIC_TRIANGLE] = nbQuadTriaGuess;
482 GuessSize[VTK_QUAD] = nbLinQuadGuess;
483 GuessSize[VTK_QUADRATIC_QUAD] = nbQuadQuadGuess;
484 GuessSize[VTK_BIQUADRATIC_QUAD] = nbQuadQuadGuess;
485 GuessSize[VTK_TETRA] = nbLinTetra;
486 GuessSize[VTK_QUADRATIC_TETRA] = nbQuadTetra;
487 GuessSize[VTK_PYRAMID] = nbLinPyra;
488 GuessSize[VTK_QUADRATIC_PYRAMID] = nbQuadPyra;
489 GuessSize[VTK_WEDGE] = nbLinPrism;
490 GuessSize[VTK_QUADRATIC_WEDGE] = nbQuadPrism;
491 GuessSize[VTK_HEXAHEDRON] = nbLinHexa;
492 GuessSize[VTK_QUADRATIC_HEXAHEDRON] = nbQuadHexa;
493 GuessSize[VTK_TRIQUADRATIC_HEXAHEDRON] = nbQuadHexa;
494 GuessSize[VTK_HEXAGONAL_PRISM] = nbHexPrism;
496 _downArray[VTK_LINE] ->allocate(nbLineGuess);
497 _downArray[VTK_QUADRATIC_EDGE] ->allocate(nbQuadEdgeGuess);
498 _downArray[VTK_TRIANGLE] ->allocate(nbLinTriaGuess);
499 _downArray[VTK_QUADRATIC_TRIANGLE] ->allocate(nbQuadTriaGuess);
500 _downArray[VTK_BIQUADRATIC_TRIANGLE] ->allocate(nbQuadTriaGuess);
501 _downArray[VTK_QUAD] ->allocate(nbLinQuadGuess);
502 _downArray[VTK_QUADRATIC_QUAD] ->allocate(nbQuadQuadGuess);
503 _downArray[VTK_BIQUADRATIC_QUAD] ->allocate(nbQuadQuadGuess);
504 _downArray[VTK_TETRA] ->allocate(nbLinTetra);
505 _downArray[VTK_QUADRATIC_TETRA] ->allocate(nbQuadTetra);
506 _downArray[VTK_PYRAMID] ->allocate(nbLinPyra);
507 _downArray[VTK_QUADRATIC_PYRAMID] ->allocate(nbQuadPyra);
508 _downArray[VTK_WEDGE] ->allocate(nbLinPrism);
509 _downArray[VTK_QUADRATIC_WEDGE] ->allocate(nbQuadPrism);
510 _downArray[VTK_HEXAHEDRON] ->allocate(nbLinHexa);
511 _downArray[VTK_QUADRATIC_HEXAHEDRON] ->allocate(nbQuadHexa);
512 _downArray[VTK_TRIQUADRATIC_HEXAHEDRON]->allocate(nbQuadHexa);
513 _downArray[VTK_HEXAGONAL_PRISM] ->allocate(nbHexPrism);
515 // --- iteration on vtkUnstructuredGrid cells, only faces
516 // for each vtk face:
517 // create a downward face entry with its downward id.
518 // compute vtk volumes, create downward volumes entry.
519 // mark face in downward volumes
520 // mark volumes in downward face
522 MESSAGE("--- iteration on vtkUnstructuredGrid cells, only faces");CHRONO(20);
523 int cellSize = this->Types->GetNumberOfTuples();
524 _cellIdToDownId.resize(cellSize, -1);
526 for (int i = 0; i < cellSize; i++)
528 int vtkFaceType = this->GetCellType(i);
529 if (SMDS_Downward::getCellDimension(vtkFaceType) == 2)
532 //ASSERT(_downArray[vtkFaceType]);
533 int connFaceId = _downArray[vtkFaceType]->addCell(vtkFaceId);
534 SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
535 downFace->setTempNodes(connFaceId, vtkFaceId);
536 int vols[2] = { -1, -1 };
537 int nbVolumes = downFace->computeVolumeIds(vtkFaceId, vols);
538 //MESSAGE("nbVolumes="<< nbVolumes);
539 for (int ivol = 0; ivol < nbVolumes; ivol++)
541 int vtkVolId = vols[ivol];
542 int vtkVolType = this->GetCellType(vtkVolId);
543 //ASSERT(_downArray[vtkVolType]);
544 int connVolId = _downArray[vtkVolType]->addCell(vtkVolId);
545 _downArray[vtkVolType]->addDownCell(connVolId, connFaceId, vtkFaceType);
546 _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId, vtkVolType);
547 // MESSAGE("Face " << vtkFaceId << " belongs to volume " << vtkVolId);
552 // --- iteration on vtkUnstructuredGrid cells, only volumes
553 // for each vtk volume:
554 // create downward volumes entry if not already done
555 // build a temporary list of faces described with their nodes
557 // compute the vtk volumes containing this face
558 // check if the face is already listed in the volumes (comparison of ordered list of nodes)
559 // if not, create a downward face entry (resizing of structure required)
560 // (the downward faces store a temporary list of nodes to ease the comparison)
561 // create downward volumes entry if not already done
562 // mark volumes in downward face
563 // mark face in downward volumes
566 MESSAGE("--- iteration on vtkUnstructuredGrid cells, only volumes");CHRONO(21);
568 for (int i = 0; i < cellSize; i++)
570 int vtkType = this->GetCellType(i);
571 if (SMDS_Downward::getCellDimension(vtkType) == 3)
575 // MESSAGE("vtk volume " << vtkVolId);
576 //ASSERT(_downArray[vtkType]);
577 /*int connVolId = */_downArray[vtkType]->addCell(vtkVolId);
579 // --- find all the faces of the volume, describe the faces by their nodes
581 SMDS_Down3D* downVol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
582 ListElemByNodesType facesWithNodes;
583 downVol->computeFacesWithNodes(vtkVolId, facesWithNodes);
584 // MESSAGE("vtk volume " << vtkVolId << " contains " << facesWithNodes.nbElems << " faces");
586 for (int iface = 0; iface < facesWithNodes.nbElems; iface++)
588 // --- find the volumes containing the face
591 int vtkFaceType = facesWithNodes.elems[iface].vtkType;
592 SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
593 int vols[2] = { -1, -1 };
594 int *nodes = &facesWithNodes.elems[iface].nodeIds[0];
595 int lg = facesWithNodes.elems[iface].nbNodes;
596 int nbVolumes = downFace->computeVolumeIdsFromNodesFace(nodes, lg, vols);
597 // MESSAGE("vtk volume " << vtkVolId << " face " << iface << " belongs to " << nbVolumes << " volumes");
599 // --- check if face is registered in the volumes
604 for (int ivol = 0; ivol < nbVolumes; ivol++)
606 int vtkVolId2 = vols[ivol];
607 int vtkVolType = this->GetCellType(vtkVolId2);
608 //ASSERT(_downArray[vtkVolType]);
609 int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
610 SMDS_Down3D* downVol2 = static_cast<SMDS_Down3D*> (_downArray[vtkVolType]);
611 connFaceId = downVol2->FindFaceByNodes(connVolId2, facesWithNodes.elems[iface]);
613 break; // --- face already created
616 // --- if face is not registered in the volumes, create face
621 connFaceId = _downArray[vtkFaceType]->addCell();
622 downFace->setTempNodes(connFaceId, facesWithNodes.elems[iface]);
625 // --- mark volumes in downward face and mark face in downward volumes
628 for (int ivol = 0; ivol < nbVolumes; ivol++)
630 int vtkVolId2 = vols[ivol];
631 int vtkVolType = this->GetCellType(vtkVolId2);
632 //ASSERT(_downArray[vtkVolType]);
633 int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
634 _downArray[vtkVolType]->addDownCell(connVolId2, connFaceId, vtkFaceType);
635 _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId2, vtkVolType);
636 // MESSAGE(" From volume " << vtkVolId << " face " << connFaceId << " belongs to volume " << vtkVolId2);
642 // --- iteration on vtkUnstructuredGrid cells, only edges
643 // for each vtk edge:
644 // create downward edge entry
645 // store the nodes id's in downward edge (redundant with vtkUnstructuredGrid)
646 // find downward faces
647 // (from vtk faces or volumes, get downward faces, they have a temporary list of nodes)
648 // mark edge in downward faces
649 // mark faces in downward edge
652 MESSAGE("--- iteration on vtkUnstructuredGrid cells, only edges");CHRONO(22);
654 for (int i = 0; i < cellSize; i++)
656 int vtkEdgeType = this->GetCellType(i);
657 if (SMDS_Downward::getCellDimension(vtkEdgeType) == 1)
660 //ASSERT(_downArray[vtkEdgeType]);
661 int connEdgeId = _downArray[vtkEdgeType]->addCell(vtkEdgeId);
662 SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
663 downEdge->setNodes(connEdgeId, vtkEdgeId);
665 int nbVtkCells = downEdge->computeVtkCells(connEdgeId, vtkIds);
667 unsigned char downTypes[1000];
668 int nbDownFaces = downEdge->computeFaces(connEdgeId, &vtkIds[0], nbVtkCells, downFaces, downTypes);
669 for (int n = 0; n < nbDownFaces; n++)
671 _downArray[downTypes[n]]->addDownCell(downFaces[n], connEdgeId, vtkEdgeType);
672 _downArray[vtkEdgeType]->addUpCell(connEdgeId, downFaces[n], downTypes[n]);
677 // --- iteration on downward faces (they are all listed now)
678 // for each downward face:
679 // build a temporary list of edges with their ordered list of nodes
681 // find all the vtk cells containing this edge
682 // then identify all the downward faces containing the edge, from the vtk cells
683 // check if the edge is already listed in the faces (comparison of ordered list of nodes)
684 // if not, create a downward edge entry with the node id's
685 // mark edge in downward faces
686 // mark downward faces in edge (size of list unknown, to be allocated)
688 CHRONOSTOP(22);CHRONO(23);
690 for (int vtkFaceType = 0; vtkFaceType < VTK_QUADRATIC_PYRAMID; vtkFaceType++)
692 if (SMDS_Downward::getCellDimension(vtkFaceType) != 2)
695 // --- find all the edges of the face, describe the edges by their nodes
697 SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
698 int maxId = downFace->getMaxId();
699 for (int faceId = 0; faceId < maxId; faceId++)
702 ListElemByNodesType edgesWithNodes;
703 downFace->computeEdgesWithNodes(faceId, edgesWithNodes);
704 // MESSAGE("downward face type " << vtkFaceType << " num " << faceId << " contains " << edgesWithNodes.nbElems << " edges");
707 for (int iedge = 0; iedge < edgesWithNodes.nbElems; iedge++)
710 // --- check if the edge is already registered by exploration of the faces
714 unsigned char vtkEdgeType = edgesWithNodes.elems[iedge].vtkType;
715 int *pts = &edgesWithNodes.elems[iedge].nodeIds[0];
716 SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
717 int nbVtkCells = downEdge->computeVtkCells(pts, vtkIds);
718 //CHRONOSTOP(41);CHRONO(42);
720 unsigned char downTypes[1000];
721 int nbDownFaces = downEdge->computeFaces(pts, &vtkIds[0], nbVtkCells, downFaces, downTypes);
726 for (int idf = 0; idf < nbDownFaces; idf++)
728 int faceId2 = downFaces[idf];
729 int faceType = downTypes[idf];
730 //ASSERT(_downArray[faceType]);
731 SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
732 connEdgeId = downFace2->FindEdgeByNodes(faceId2, edgesWithNodes.elems[iedge]);
734 break; // --- edge already created
737 // --- if edge is not registered in the faces, create edge
742 connEdgeId = _downArray[vtkEdgeType]->addCell();
743 downEdge->setNodes(connEdgeId, edgesWithNodes.elems[iedge].nodeIds);
747 // --- mark faces in downward edge and mark edge in downward faces
750 for (int idf = 0; idf < nbDownFaces; idf++)
752 int faceId2 = downFaces[idf];
753 int faceType = downTypes[idf];
754 //ASSERT(_downArray[faceType]);
755 //SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
756 _downArray[vtkEdgeType]->addUpCell(connEdgeId, faceId2, faceType);
757 _downArray[faceType]->addDownCell(faceId2, connEdgeId, vtkEdgeType);
758 // MESSAGE(" From face t:" << vtkFaceType << " " << faceId <<
759 // " edge " << connEdgeId << " belongs to face t:" << faceType << " " << faceId2);
765 CHRONOSTOP(23);CHRONO(24);
767 // compact downward connectivity structure: adjust downward arrays size, replace vector<vector int>> by a single vector<int>
768 // 3D first then 2D and last 1D to release memory before edge upCells reorganization, (temporary memory use)
770 for (int vtkType = VTK_QUADRATIC_PYRAMID; vtkType >= 0; vtkType--)
772 if (SMDS_Downward *down = _downArray[vtkType])
774 down->compactStorage();
780 for (int vtkType = 0; vtkType <= VTK_QUADRATIC_PYRAMID; vtkType++)
782 if (SMDS_Downward *down = _downArray[vtkType])
784 if (down->getMaxId())
786 MESSAGE("Cells of Type " << vtkType << " : number of entities, est: "
787 << GuessSize[vtkType] << " real: " << down->getMaxId());
790 }CHRONOSTOP(24);CHRONOSTOP(2);
794 /*! Get the neighbors of a cell.
795 * Only the neighbors having the dimension of the cell are taken into account
796 * (neighbors of a volume are the volumes sharing a face with this volume,
797 * neighbors of a face are the faces sharing an edge with this face...).
798 * @param neighborsVtkIds vector of neighbors vtk id's to fill (reserve enough space).
799 * @param downIds downward id's of cells of dimension n-1, to fill (reserve enough space).
800 * @param downTypes vtk types of cells of dimension n-1, to fill (reserve enough space).
801 * @param vtkId the vtk id of the cell
802 * @return number of neighbors
804 int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId, bool getSkin)
806 int vtkType = this->GetCellType(vtkId);
807 int cellDim = SMDS_Downward::getCellDimension(vtkType);
809 return 0; // TODO voisins des edges = edges connectees
810 int cellId = this->_cellIdToDownId[vtkId];
812 int nbCells = _downArray[vtkType]->getNumberOfDownCells(cellId);
813 const int *downCells = _downArray[vtkType]->getDownCells(cellId);
814 const unsigned char* downTyp = _downArray[vtkType]->getDownTypes(cellId);
816 // --- iteration on faces of the 3D cell (or edges on the 2D cell).
819 for (int i = 0; i < nbCells; i++)
821 int downId = downCells[i];
822 int cellType = downTyp[i];
823 int nbUp = _downArray[cellType]->getNumberOfUpCells(downId);
824 const int *upCells = _downArray[cellType]->getUpCells(downId);
825 const unsigned char* upTypes = _downArray[cellType]->getUpTypes(downId);
827 // ---for a volume, max 2 upCells, one is this cell, the other is a neighbor
828 // for a face, number of neighbors (connected faces) not known
830 for (int j = 0; j < nbUp; j++)
832 if ((upCells[j] == cellId) && (upTypes[j] == vtkType))
834 int vtkNeighbor = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
835 neighborsVtkIds[nb] = vtkNeighbor;
836 downIds[nb] = downId;
837 downTypes[nb] = cellType;
839 if (nb >= NBMAXNEIGHBORS)
841 INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
847 if (cellDim == 3 && nbUp == 1) // this face is on the skin of the volume
849 neighborsVtkIds[nb] = _downArray[cellType]->getVtkCellId(downId); // OK if skin present
850 downIds[nb] = downId;
851 downTypes[nb] = cellType;
853 if (nb >= NBMAXNEIGHBORS)
855 INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
864 /*! get the volumes containing a face or an edge of the grid
865 * The edge or face belongs to the vtkUnstructuredGrid
866 * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
867 * @param vtkId vtk id of the face or edge
869 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int vtkId)
871 int vtkType = this->GetCellType(vtkId);
872 int dim = SMDS_Downward::getCellDimension(vtkType);
874 unsigned char cellTypes[1000];
875 int downCellId[1000];
878 int downId = this->CellIdToDownId(vtkId);
881 MESSAGE("Downward structure not up to date: new edge not taken into account");
884 nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
885 const int *upCells = _downArray[vtkType]->getUpCells(downId);
886 const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
887 for (int i=0; i< nbFaces; i++)
889 cellTypes[i] = upTypes[i];
890 downCellId[i] = upCells[i];
896 cellTypes[0] = this->GetCellType(vtkId);
897 int downId = this->CellIdToDownId(vtkId);
900 MESSAGE("Downward structure not up to date: new face not taken into account");
903 downCellId[0] = downId;
907 for (int i=0; i<nbFaces; i++)
909 int vtkTypeFace = cellTypes[i];
910 int downId = downCellId[i];
911 int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
912 const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
913 const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
914 for (int j=0; j<nv; j++)
916 int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
918 volVtkIds[nbvol++] = vtkVolId;
924 /*! get the volumes containing a face or an edge of the downward structure
925 * The edge or face does not necessary belong to the vtkUnstructuredGrid
926 * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
927 * @param downId id in the downward structure
928 * @param downType type of cell
930 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int downId, unsigned char downType)
932 int vtkType = downType;
933 int dim = SMDS_Downward::getCellDimension(vtkType);
935 unsigned char cellTypes[1000];
936 int downCellId[1000];
939 nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
940 const int *upCells = _downArray[vtkType]->getUpCells(downId);
941 const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
942 for (int i=0; i< nbFaces; i++)
944 cellTypes[i] = upTypes[i];
945 downCellId[i] = upCells[i];
951 cellTypes[0] = vtkType;
952 downCellId[0] = downId;
956 for (int i=0; i<nbFaces; i++)
958 int vtkTypeFace = cellTypes[i];
959 int downId = downCellId[i];
960 int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
961 const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
962 const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
963 for (int j=0; j<nv; j++)
965 int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
967 volVtkIds[nbvol++] = vtkVolId;
973 /*! get the node id's of a cell.
974 * The cell is defined by it's downward connectivity id and type.
975 * @param nodeSet set of of vtk node id's to fill.
976 * @param downId downward connectivity id of the cell.
977 * @param downType type of cell.
979 void SMDS_UnstructuredGrid::GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType)
981 _downArray[downType]->getNodeIds(downId, nodeSet);
984 /*! change some nodes in cell without modifying type or internal connectivity.
985 * Nodes inverse connectivity is maintained up to date.
986 * @param vtkVolId vtk id of the cell
987 * @param localClonedNodeIds map old node id to new node id.
989 void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds)
992 vtkIdType *pts; // will refer to the point id's of the face
993 this->GetCellPoints(vtkVolId, npts, pts);
994 for (int i = 0; i < npts; i++)
996 if (localClonedNodeIds.count(pts[i]))
998 vtkIdType oldpt = pts[i];
999 pts[i] = localClonedNodeIds[oldpt];
1000 //MESSAGE(oldpt << " --> " << pts[i]);
1001 //this->RemoveReferenceToCell(oldpt, vtkVolId);
1002 //this->AddReferenceToCell(pts[i], vtkVolId);
1007 /*! reorder the nodes of a face
1008 * @param vtkVolId vtk id of a volume containing the face, to get an orientation for the face.
1009 * @param orderedNodes list of nodes to reorder (in out)
1010 * @return size of the list
1012 int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes)
1014 int vtkType = this->GetCellType(vtkVolId);
1015 dim = SMDS_Downward::getCellDimension(vtkType);
1018 SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
1019 int downVolId = this->_cellIdToDownId[vtkVolId];
1020 downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
1022 // else nothing to do;
1023 return orderedNodes.size();
1026 void SMDS_UnstructuredGrid::BuildLinks()
1028 // Remove the old links if they are already built
1031 this->Links->UnRegister(this);
1034 SMDS_CellLinks* links;
1035 this->Links = links = SMDS_CellLinks::New();
1036 this->Links->Allocate(this->GetNumberOfPoints());
1037 this->Links->Register(this);
1038 links->BuildLinks(this, this->Connectivity,this->GetCellTypesArray() );
1039 this->Links->Delete();
1042 void SMDS_UnstructuredGrid::DeleteLinks()
1044 // Remove the old links if they are already built
1047 this->Links->UnRegister(this);
1051 SMDS_CellLinks* SMDS_UnstructuredGrid::GetLinks()
1055 return static_cast< SMDS_CellLinks* >( this->Links );
1058 /*! Create a volume (prism or hexahedron) by duplication of a face.
1059 * Designed for use in creation of flat elements separating volume domains.
1060 * A face separating two domains is shared by two volume cells.
1061 * All the nodes are already created (for the two faces).
1062 * Each original Node is associated to corresponding nodes in the domains.
1063 * Some nodes may be duplicated for more than two domains, when domain separations intersect.
1064 * In that case, even some of the nodes to use for the original face may be changed.
1065 * @param vtkVolId: vtk id of a volume containing the face, to get an orientation for the face.
1066 * @param domain1: domain of the original face
1067 * @param domain2: domain of the duplicated face
1068 * @param originalNodes: the vtk node ids of the original face
1069 * @param nodeDomains: map(original id --> map(domain --> duplicated node id))
1070 * @return ok if success.
1073 SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
1076 std::set<int>& originalNodes,
1077 std::map<int, std::map<int, int> >& nodeDomains,
1078 std::map<int, std::map<long, int> >& nodeQuadDomains)
1080 //MESSAGE("extrudeVolumeFromFace " << vtkVolId);
1081 vector<vtkIdType> orderedOriginals;
1082 orderedOriginals.clear();
1083 set<int>::const_iterator it = originalNodes.begin();
1084 for (; it != originalNodes.end(); ++it)
1085 orderedOriginals.push_back(*it);
1088 int nbNodes = this->getOrderedNodesOfFace(vtkVolId, dim, orderedOriginals);
1089 vector<vtkIdType> orderedNodes;
1091 bool isQuadratic = false;
1092 switch (orderedOriginals.size())
1103 isQuadratic = false;
1109 long dom1 = domain1;
1110 long dom2 = domain2;
1111 long dom1_2; // for nodeQuadDomains
1112 if (domain1 < domain2)
1113 dom1_2 = dom1 + INT_MAX * dom2;
1115 dom1_2 = dom2 + INT_MAX * dom1;
1116 //cerr << "dom1=" << dom1 << " dom2=" << dom2 << " dom1_2=" << dom1_2 << endl;
1117 int ima = orderedOriginals.size();
1118 int mid = orderedOriginals.size() / 2;
1119 //cerr << "ima=" << ima << " mid=" << mid << endl;
1120 for (int i = 0; i < mid; i++)
1121 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1122 for (int i = 0; i < mid; i++)
1123 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1124 for (int i = mid; i < ima; i++)
1125 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1126 for (int i = mid; i < ima; i++)
1127 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1128 for (int i = 0; i < mid; i++)
1130 int oldId = orderedOriginals[i];
1132 if (nodeQuadDomains.count(oldId) && nodeQuadDomains[oldId].count(dom1_2))
1133 newId = nodeQuadDomains[oldId][dom1_2];
1136 double *coords = this->GetPoint(oldId);
1137 SMDS_MeshNode *newNode = _mesh->AddNode(coords[0], coords[1], coords[2]);
1138 newId = newNode->getVtkId();
1139 if (! nodeQuadDomains.count(oldId))
1141 std::map<long, int> emptyMap;
1142 nodeQuadDomains[oldId] = emptyMap;
1144 nodeQuadDomains[oldId][dom1_2] = newId;
1146 orderedNodes.push_back(newId);
1151 for (int i = 0; i < nbNodes; i++)
1152 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1154 for (int i = 0; i < nbNodes; i++)
1155 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1157 for (int i = nbNodes-1; i >= 0; i--)
1158 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1164 SMDS_MeshVolume *vol = _mesh->AddVolumeFromVtkIds(orderedNodes);
1169 SMDS_MeshFace *face = _mesh->AddFaceFromVtkIds(orderedNodes);
1173 // TODO update sub-shape list of elements and nodes
1177 //================================================================================
1179 * \brief Allocates data array for ball diameters
1180 * \param MaxVtkID - max ID of a ball element
1182 //================================================================================
1184 void SMDS_UnstructuredGrid::AllocateDiameters( vtkIdType MaxVtkID )
1186 SetBallDiameter( MaxVtkID, 0 );
1189 //================================================================================
1191 * \brief Sets diameter of a ball element
1192 * \param vtkID - vtk id of the ball element
1193 * \param diameter - diameter of the ball element
1195 //================================================================================
1197 void SMDS_UnstructuredGrid::SetBallDiameter( vtkIdType vtkID, double diameter )
1199 vtkDoubleArray* array = vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() );
1202 array = vtkDoubleArray::New();
1203 array->SetNumberOfComponents(1);
1204 vtkDataSet::CellData->SetScalars( array );
1206 array->InsertValue( vtkID, diameter );
1209 //================================================================================
1211 * \brief Returns diameter of a ball element
1212 * \param vtkID - vtk id of the ball element
1214 //================================================================================
1216 double SMDS_UnstructuredGrid::GetBallDiameter( vtkIdType vtkID ) const
1218 if ( vtkDataSet::CellData )
1219 return vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )->GetValue( vtkID );