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 SMDS_CellLinks::SMDS_CellLinks() :
61 SMDS_CellLinks::~SMDS_CellLinks()
65 SMDS_UnstructuredGrid* SMDS_UnstructuredGrid::New()
67 return new SMDS_UnstructuredGrid();
70 SMDS_UnstructuredGrid::SMDS_UnstructuredGrid() :
73 _cellIdToDownId.clear();
79 SMDS_UnstructuredGrid::~SMDS_UnstructuredGrid()
83 unsigned long SMDS_UnstructuredGrid::GetMTime()
85 unsigned long mtime = vtkUnstructuredGrid::GetMTime();
88 // OUV_PORTING_VTK6: seems to be useless
90 void SMDS_UnstructuredGrid::Update()
92 return vtkUnstructuredGrid::Update();
95 void SMDS_UnstructuredGrid::UpdateInformation()
97 return vtkUnstructuredGrid::UpdateInformation();
100 vtkPoints* SMDS_UnstructuredGrid::GetPoints()
102 // TODO erreur incomprehensible de la macro vtk GetPoints apparue avec la version paraview de fin aout 2010
106 //#ifdef VTK_HAVE_POLYHEDRON
107 int SMDS_UnstructuredGrid::InsertNextLinkedCell(int type, int npts, vtkIdType *pts)
109 if (type != VTK_POLYHEDRON)
110 return vtkUnstructuredGrid::InsertNextLinkedCell(type, npts, pts);
112 // --- type = VTK_POLYHEDRON
113 int cellid = this->InsertNextCell(type, npts, pts);
115 set<vtkIdType> setOfNodes;
119 for (int nf = 0; nf < nbfaces; nf++)
121 int nbnodes = pts[i];
123 for (int k = 0; k < nbnodes; k++)
125 setOfNodes.insert(pts[i]);
130 set<vtkIdType>::iterator it = setOfNodes.begin();
131 for (; it != setOfNodes.end(); ++it)
133 this->Links->ResizeCellList(*it, 1);
134 this->Links->AddCellReference(cellid, *it);
141 void SMDS_UnstructuredGrid::setSMDS_mesh(SMDS_Mesh *mesh)
146 void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize,
147 std::vector<int>& idCellsOldToNew, int newCellSize)
149 //MESSAGE("------------------------- compactGrid " << newNodeSize << " " << newCellSize);//CHRONO(1);
150 int alreadyCopied = 0;
152 // --- if newNodeSize, create a new compacted vtkPoints
154 vtkPoints *newPoints = vtkPoints::New();
155 newPoints->SetDataType(VTK_DOUBLE);
156 newPoints->SetNumberOfPoints(newNodeSize);
159 // rnv: to fix bug "21125: EDF 1233 SMESH: Degradation of precision in a test case for quadratic conversion"
160 // using double type for storing coordinates of nodes instead float.
161 int oldNodeSize = idNodesOldToNew.size();
164 while ( i < oldNodeSize )
166 // skip a hole if any
167 while ( i < oldNodeSize && idNodesOldToNew[i] < 0 )
170 // look for a block end
171 while ( i < oldNodeSize && idNodesOldToNew[i] >= 0 )
174 copyNodes(newPoints, idNodesOldToNew, alreadyCopied, startBloc, endBloc);
176 newPoints->Squeeze();
179 // --- create new compacted Connectivity, Locations and Types
181 int oldCellSize = this->Types->GetNumberOfTuples();
183 vtkCellArray *newConnectivity = vtkCellArray::New();
184 newConnectivity->Initialize();
185 int oldCellDataSize = this->Connectivity->GetData()->GetSize();
186 newConnectivity->Allocate(oldCellDataSize);
188 vtkUnsignedCharArray *newTypes = vtkUnsignedCharArray::New();
189 newTypes->Initialize();
190 newTypes->SetNumberOfValues(newCellSize);
192 vtkIdTypeArray *newLocations = vtkIdTypeArray::New();
193 newLocations->Initialize();
194 newLocations->SetNumberOfValues(newCellSize);
196 // TODO some polyhedron may be huge (only in some tests)
197 vtkIdType tmpid[NBMAXNODESINCELL];
198 vtkIdType *pointsCell = &tmpid[0]; // --- points id to fill a new cell
202 while ( i < oldCellSize )
204 // skip a hole if any
205 while ( i < oldCellSize && this->Types->GetValue(i) == VTK_EMPTY_CELL )
208 // look for a block end
209 while ( i < oldCellSize && this->Types->GetValue(i) != VTK_EMPTY_CELL )
212 if ( endBloc > startBloc )
214 idCellsOldToNew, idNodesOldToNew,
215 newConnectivity, newLocations,
216 pointsCell, alreadyCopied,
219 newConnectivity->Squeeze();
221 if (1/*newNodeSize*/)
223 this->SetPoints(newPoints);
226 if (vtkDoubleArray* diameters =
227 vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )) // Balls
229 for (int oldCellID = 0; oldCellID < oldCellSize; oldCellID++)
231 if (this->Types->GetValue(oldCellID) == VTK_EMPTY_CELL)
233 int newCellId = idCellsOldToNew[ oldCellID ];
234 if (newTypes->GetValue(newCellId) == VTK_POLY_VERTEX)
235 diameters->SetValue( newCellId, diameters->GetValue( oldCellID ));
239 if (this->FaceLocations)
241 vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New();
242 newFaceLocations->Initialize();
243 newFaceLocations->Allocate(newTypes->GetSize());
244 vtkIdTypeArray *newFaces = vtkIdTypeArray::New();
245 newFaces->Initialize();
246 newFaces->Allocate(this->Faces->GetSize());
247 for (int i = 0; i < oldCellSize; i++)
249 if (this->Types->GetValue(i) == VTK_EMPTY_CELL)
251 int newCellId = idCellsOldToNew[i];
252 if (newTypes->GetValue(newCellId) == VTK_POLYHEDRON)
254 newFaceLocations->InsertNextValue(newFaces->GetMaxId()+1);
255 int oldFaceLoc = this->FaceLocations->GetValue(i);
256 int nCellFaces = this->Faces->GetValue(oldFaceLoc++);
257 newFaces->InsertNextValue(nCellFaces);
258 for (int n=0; n<nCellFaces; n++)
260 int nptsInFace = this->Faces->GetValue(oldFaceLoc++);
261 newFaces->InsertNextValue(nptsInFace);
262 for (int k=0; k<nptsInFace; k++)
264 int oldpt = this->Faces->GetValue(oldFaceLoc++);
265 newFaces->InsertNextValue(idNodesOldToNew[oldpt]);
271 newFaceLocations->InsertNextValue(-1);
274 newFaceLocations->Squeeze();
276 this->SetCells(newTypes, newLocations, newConnectivity, newFaceLocations, newFaces);
277 newFaceLocations->Delete();
282 this->SetCells(newTypes, newLocations, newConnectivity, FaceLocations, Faces);
287 newLocations->Delete();
288 newConnectivity->Delete();
292 void SMDS_UnstructuredGrid::copyNodes(vtkPoints * newPoints,
293 std::vector<int>& idNodesOldToNew,
298 void *target = newPoints->GetVoidPointer(3 * alreadyCopied);
299 void *source = this->Points->GetVoidPointer(3 * start);
300 int nbPoints = end - start;
303 memcpy(target, source, 3 * sizeof(double) * nbPoints);
304 for (int j = start; j < end; j++)
305 idNodesOldToNew[j] = alreadyCopied++; // old vtkId --> new vtkId
309 void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes,
310 std::vector<int>& idCellsOldToNew,
311 std::vector<int>& idNodesOldToNew,
312 vtkCellArray* newConnectivity,
313 vtkIdTypeArray* newLocations,
314 vtkIdType* pointsCell,
319 for (int j = start; j < end; j++)
321 newTypes->SetValue(alreadyCopied, this->Types->GetValue(j));
322 idCellsOldToNew[j] = alreadyCopied; // old vtkId --> new vtkId
323 vtkIdType oldLoc = this->Locations->GetValue(j);
325 vtkIdType *oldPtsCell = 0;
326 this->Connectivity->GetCell(oldLoc, nbpts, oldPtsCell);
327 assert(nbpts < NBMAXNODESINCELL);
328 for (int l = 0; l < nbpts; l++)
330 int oldval = oldPtsCell[l];
331 pointsCell[l] = idNodesOldToNew[oldval];
333 /*int newcnt = */newConnectivity->InsertNextCell(nbpts, pointsCell);
334 int newLoc = newConnectivity->GetInsertLocation(nbpts);
335 newLocations->SetValue(alreadyCopied, newLoc);
340 int SMDS_UnstructuredGrid::CellIdToDownId(int vtkCellId)
342 if ((vtkCellId < 0) || (vtkCellId >= (int)_cellIdToDownId.size()))
346 return _cellIdToDownId[vtkCellId];
349 void SMDS_UnstructuredGrid::setCellIdToDownId(int vtkCellId, int downId)
351 // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
352 _cellIdToDownId[vtkCellId] = downId;
355 void SMDS_UnstructuredGrid::CleanDownwardConnectivity()
357 for (size_t i = 0; i < _downArray.size(); i++)
360 delete _downArray[i];
363 _cellIdToDownId.clear();
366 /*! Build downward connectivity: to do only when needed because heavy memory load.
367 * Downward connectivity is no more valid if vtkUnstructuredGrid is modified.
370 void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
372 MESSAGE("SMDS_UnstructuredGrid::BuildDownwardConnectivity");CHRONO(2);
373 // TODO calcul partiel sans edges
375 // --- erase previous data if any
377 this->CleanDownwardConnectivity();
379 // --- create SMDS_Downward structures (in _downArray vector[vtkCellType])
381 _downArray.resize(VTK_MAXTYPE + 1, 0);
383 _downArray[VTK_LINE] = new SMDS_DownEdge(this);
384 _downArray[VTK_QUADRATIC_EDGE] = new SMDS_DownQuadEdge(this);
385 _downArray[VTK_TRIANGLE] = new SMDS_DownTriangle(this);
386 _downArray[VTK_QUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
387 _downArray[VTK_BIQUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
388 _downArray[VTK_QUAD] = new SMDS_DownQuadrangle(this);
389 _downArray[VTK_QUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
390 _downArray[VTK_BIQUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
391 _downArray[VTK_TETRA] = new SMDS_DownTetra(this);
392 _downArray[VTK_QUADRATIC_TETRA] = new SMDS_DownQuadTetra(this);
393 _downArray[VTK_PYRAMID] = new SMDS_DownPyramid(this);
394 _downArray[VTK_QUADRATIC_PYRAMID] = new SMDS_DownQuadPyramid(this);
395 _downArray[VTK_WEDGE] = new SMDS_DownPenta(this);
396 _downArray[VTK_QUADRATIC_WEDGE] = new SMDS_DownQuadPenta(this);
397 _downArray[VTK_HEXAHEDRON] = new SMDS_DownHexa(this);
398 _downArray[VTK_QUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
399 _downArray[VTK_TRIQUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
400 _downArray[VTK_HEXAGONAL_PRISM] = new SMDS_DownPenta(this);
402 // --- get detailed info of number of cells of each type, allocate SMDS_downward structures
404 const SMDS_MeshInfo &meshInfo = _mesh->GetMeshInfo();
406 int nbLinTetra = meshInfo.NbTetras (ORDER_LINEAR);
407 int nbQuadTetra = meshInfo.NbTetras (ORDER_QUADRATIC);
408 int nbLinPyra = meshInfo.NbPyramids(ORDER_LINEAR);
409 int nbQuadPyra = meshInfo.NbPyramids(ORDER_QUADRATIC);
410 int nbLinPrism = meshInfo.NbPrisms (ORDER_LINEAR);
411 int nbQuadPrism = meshInfo.NbPrisms (ORDER_QUADRATIC);
412 int nbLinHexa = meshInfo.NbHexas (ORDER_LINEAR);
413 int nbQuadHexa = meshInfo.NbHexas (ORDER_QUADRATIC);
414 int nbHexPrism = meshInfo.NbHexPrisms();
416 int nbLineGuess = int((4.0 / 3.0) * nbLinTetra + 2 * nbLinPrism + 2.5 * nbLinPyra + 3 * nbLinHexa);
417 int nbQuadEdgeGuess = int((4.0 / 3.0) * nbQuadTetra + 2 * nbQuadPrism + 2.5 * nbQuadPyra + 3 * nbQuadHexa);
418 int nbLinTriaGuess = 2 * nbLinTetra + nbLinPrism + 2 * nbLinPyra;
419 int nbQuadTriaGuess = 2 * nbQuadTetra + nbQuadPrism + 2 * nbQuadPyra;
420 int nbLinQuadGuess = int((2.0 / 3.0) * nbLinPrism + (1.0 / 2.0) * nbLinPyra + 3 * nbLinHexa);
421 int nbQuadQuadGuess = int((2.0 / 3.0) * nbQuadPrism + (1.0 / 2.0) * nbQuadPyra + 3 * nbQuadHexa);
423 int GuessSize[VTK_MAXTYPE];
424 GuessSize[VTK_LINE] = nbLineGuess;
425 GuessSize[VTK_QUADRATIC_EDGE] = nbQuadEdgeGuess;
426 GuessSize[VTK_TRIANGLE] = nbLinTriaGuess;
427 GuessSize[VTK_QUADRATIC_TRIANGLE] = nbQuadTriaGuess;
428 GuessSize[VTK_BIQUADRATIC_TRIANGLE] = nbQuadTriaGuess;
429 GuessSize[VTK_QUAD] = nbLinQuadGuess;
430 GuessSize[VTK_QUADRATIC_QUAD] = nbQuadQuadGuess;
431 GuessSize[VTK_BIQUADRATIC_QUAD] = nbQuadQuadGuess;
432 GuessSize[VTK_TETRA] = nbLinTetra;
433 GuessSize[VTK_QUADRATIC_TETRA] = nbQuadTetra;
434 GuessSize[VTK_PYRAMID] = nbLinPyra;
435 GuessSize[VTK_QUADRATIC_PYRAMID] = nbQuadPyra;
436 GuessSize[VTK_WEDGE] = nbLinPrism;
437 GuessSize[VTK_QUADRATIC_WEDGE] = nbQuadPrism;
438 GuessSize[VTK_HEXAHEDRON] = nbLinHexa;
439 GuessSize[VTK_QUADRATIC_HEXAHEDRON] = nbQuadHexa;
440 GuessSize[VTK_TRIQUADRATIC_HEXAHEDRON] = nbQuadHexa;
441 GuessSize[VTK_HEXAGONAL_PRISM] = nbHexPrism;
443 _downArray[VTK_LINE] ->allocate(nbLineGuess);
444 _downArray[VTK_QUADRATIC_EDGE] ->allocate(nbQuadEdgeGuess);
445 _downArray[VTK_TRIANGLE] ->allocate(nbLinTriaGuess);
446 _downArray[VTK_QUADRATIC_TRIANGLE] ->allocate(nbQuadTriaGuess);
447 _downArray[VTK_BIQUADRATIC_TRIANGLE] ->allocate(nbQuadTriaGuess);
448 _downArray[VTK_QUAD] ->allocate(nbLinQuadGuess);
449 _downArray[VTK_QUADRATIC_QUAD] ->allocate(nbQuadQuadGuess);
450 _downArray[VTK_BIQUADRATIC_QUAD] ->allocate(nbQuadQuadGuess);
451 _downArray[VTK_TETRA] ->allocate(nbLinTetra);
452 _downArray[VTK_QUADRATIC_TETRA] ->allocate(nbQuadTetra);
453 _downArray[VTK_PYRAMID] ->allocate(nbLinPyra);
454 _downArray[VTK_QUADRATIC_PYRAMID] ->allocate(nbQuadPyra);
455 _downArray[VTK_WEDGE] ->allocate(nbLinPrism);
456 _downArray[VTK_QUADRATIC_WEDGE] ->allocate(nbQuadPrism);
457 _downArray[VTK_HEXAHEDRON] ->allocate(nbLinHexa);
458 _downArray[VTK_QUADRATIC_HEXAHEDRON] ->allocate(nbQuadHexa);
459 _downArray[VTK_TRIQUADRATIC_HEXAHEDRON]->allocate(nbQuadHexa);
460 _downArray[VTK_HEXAGONAL_PRISM] ->allocate(nbHexPrism);
462 // --- iteration on vtkUnstructuredGrid cells, only faces
463 // for each vtk face:
464 // create a downward face entry with its downward id.
465 // compute vtk volumes, create downward volumes entry.
466 // mark face in downward volumes
467 // mark volumes in downward face
469 MESSAGE("--- iteration on vtkUnstructuredGrid cells, only faces");CHRONO(20);
470 int cellSize = this->Types->GetNumberOfTuples();
471 _cellIdToDownId.resize(cellSize, -1);
473 for (int i = 0; i < cellSize; i++)
475 int vtkFaceType = this->GetCellType(i);
476 if (SMDS_Downward::getCellDimension(vtkFaceType) == 2)
479 //ASSERT(_downArray[vtkFaceType]);
480 int connFaceId = _downArray[vtkFaceType]->addCell(vtkFaceId);
481 SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
482 downFace->setTempNodes(connFaceId, vtkFaceId);
483 int vols[2] = { -1, -1 };
484 int nbVolumes = downFace->computeVolumeIds(vtkFaceId, vols);
485 //MESSAGE("nbVolumes="<< nbVolumes);
486 for (int ivol = 0; ivol < nbVolumes; ivol++)
488 int vtkVolId = vols[ivol];
489 int vtkVolType = this->GetCellType(vtkVolId);
490 //ASSERT(_downArray[vtkVolType]);
491 int connVolId = _downArray[vtkVolType]->addCell(vtkVolId);
492 _downArray[vtkVolType]->addDownCell(connVolId, connFaceId, vtkFaceType);
493 _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId, vtkVolType);
494 // MESSAGE("Face " << vtkFaceId << " belongs to volume " << vtkVolId);
499 // --- iteration on vtkUnstructuredGrid cells, only volumes
500 // for each vtk volume:
501 // create downward volumes entry if not already done
502 // build a temporary list of faces described with their nodes
504 // compute the vtk volumes containing this face
505 // check if the face is already listed in the volumes (comparison of ordered list of nodes)
506 // if not, create a downward face entry (resizing of structure required)
507 // (the downward faces store a temporary list of nodes to ease the comparison)
508 // create downward volumes entry if not already done
509 // mark volumes in downward face
510 // mark face in downward volumes
513 MESSAGE("--- iteration on vtkUnstructuredGrid cells, only volumes");CHRONO(21);
515 for (int i = 0; i < cellSize; i++)
517 int vtkType = this->GetCellType(i);
518 if (SMDS_Downward::getCellDimension(vtkType) == 3)
522 // MESSAGE("vtk volume " << vtkVolId);
523 //ASSERT(_downArray[vtkType]);
524 /*int connVolId = */_downArray[vtkType]->addCell(vtkVolId);
526 // --- find all the faces of the volume, describe the faces by their nodes
528 SMDS_Down3D* downVol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
529 ListElemByNodesType facesWithNodes;
530 downVol->computeFacesWithNodes(vtkVolId, facesWithNodes);
531 // MESSAGE("vtk volume " << vtkVolId << " contains " << facesWithNodes.nbElems << " faces");
533 for (int iface = 0; iface < facesWithNodes.nbElems; iface++)
535 // --- find the volumes containing the face
538 int vtkFaceType = facesWithNodes.elems[iface].vtkType;
539 SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
540 int vols[2] = { -1, -1 };
541 int *nodes = &facesWithNodes.elems[iface].nodeIds[0];
542 int lg = facesWithNodes.elems[iface].nbNodes;
543 int nbVolumes = downFace->computeVolumeIdsFromNodesFace(nodes, lg, vols);
544 // MESSAGE("vtk volume " << vtkVolId << " face " << iface << " belongs to " << nbVolumes << " volumes");
546 // --- check if face is registered in the volumes
551 for (int ivol = 0; ivol < nbVolumes; ivol++)
553 int vtkVolId2 = vols[ivol];
554 int vtkVolType = this->GetCellType(vtkVolId2);
555 //ASSERT(_downArray[vtkVolType]);
556 int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
557 SMDS_Down3D* downVol2 = static_cast<SMDS_Down3D*> (_downArray[vtkVolType]);
558 connFaceId = downVol2->FindFaceByNodes(connVolId2, facesWithNodes.elems[iface]);
560 break; // --- face already created
563 // --- if face is not registered in the volumes, create face
568 connFaceId = _downArray[vtkFaceType]->addCell();
569 downFace->setTempNodes(connFaceId, facesWithNodes.elems[iface]);
572 // --- mark volumes in downward face and mark face in downward volumes
575 for (int ivol = 0; ivol < nbVolumes; ivol++)
577 int vtkVolId2 = vols[ivol];
578 int vtkVolType = this->GetCellType(vtkVolId2);
579 //ASSERT(_downArray[vtkVolType]);
580 int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
581 _downArray[vtkVolType]->addDownCell(connVolId2, connFaceId, vtkFaceType);
582 _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId2, vtkVolType);
583 // MESSAGE(" From volume " << vtkVolId << " face " << connFaceId << " belongs to volume " << vtkVolId2);
589 // --- iteration on vtkUnstructuredGrid cells, only edges
590 // for each vtk edge:
591 // create downward edge entry
592 // store the nodes id's in downward edge (redundant with vtkUnstructuredGrid)
593 // find downward faces
594 // (from vtk faces or volumes, get downward faces, they have a temporary list of nodes)
595 // mark edge in downward faces
596 // mark faces in downward edge
599 MESSAGE("--- iteration on vtkUnstructuredGrid cells, only edges");CHRONO(22);
601 for (int i = 0; i < cellSize; i++)
603 int vtkEdgeType = this->GetCellType(i);
604 if (SMDS_Downward::getCellDimension(vtkEdgeType) == 1)
607 //ASSERT(_downArray[vtkEdgeType]);
608 int connEdgeId = _downArray[vtkEdgeType]->addCell(vtkEdgeId);
609 SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
610 downEdge->setNodes(connEdgeId, vtkEdgeId);
612 int nbVtkCells = downEdge->computeVtkCells(connEdgeId, vtkIds);
614 unsigned char downTypes[1000];
615 int nbDownFaces = downEdge->computeFaces(connEdgeId, &vtkIds[0], nbVtkCells, downFaces, downTypes);
616 for (int n = 0; n < nbDownFaces; n++)
618 _downArray[downTypes[n]]->addDownCell(downFaces[n], connEdgeId, vtkEdgeType);
619 _downArray[vtkEdgeType]->addUpCell(connEdgeId, downFaces[n], downTypes[n]);
624 // --- iteration on downward faces (they are all listed now)
625 // for each downward face:
626 // build a temporary list of edges with their ordered list of nodes
628 // find all the vtk cells containing this edge
629 // then identify all the downward faces containing the edge, from the vtk cells
630 // check if the edge is already listed in the faces (comparison of ordered list of nodes)
631 // if not, create a downward edge entry with the node id's
632 // mark edge in downward faces
633 // mark downward faces in edge (size of list unknown, to be allocated)
635 CHRONOSTOP(22);CHRONO(23);
637 for (int vtkFaceType = 0; vtkFaceType < VTK_QUADRATIC_PYRAMID; vtkFaceType++)
639 if (SMDS_Downward::getCellDimension(vtkFaceType) != 2)
642 // --- find all the edges of the face, describe the edges by their nodes
644 SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
645 int maxId = downFace->getMaxId();
646 for (int faceId = 0; faceId < maxId; faceId++)
649 ListElemByNodesType edgesWithNodes;
650 downFace->computeEdgesWithNodes(faceId, edgesWithNodes);
651 // MESSAGE("downward face type " << vtkFaceType << " num " << faceId << " contains " << edgesWithNodes.nbElems << " edges");
654 for (int iedge = 0; iedge < edgesWithNodes.nbElems; iedge++)
657 // --- check if the edge is already registered by exploration of the faces
661 unsigned char vtkEdgeType = edgesWithNodes.elems[iedge].vtkType;
662 int *pts = &edgesWithNodes.elems[iedge].nodeIds[0];
663 SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
664 int nbVtkCells = downEdge->computeVtkCells(pts, vtkIds);
665 //CHRONOSTOP(41);CHRONO(42);
667 unsigned char downTypes[1000];
668 int nbDownFaces = downEdge->computeFaces(pts, &vtkIds[0], nbVtkCells, downFaces, downTypes);
673 for (int idf = 0; idf < nbDownFaces; idf++)
675 int faceId2 = downFaces[idf];
676 int faceType = downTypes[idf];
677 //ASSERT(_downArray[faceType]);
678 SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
679 connEdgeId = downFace2->FindEdgeByNodes(faceId2, edgesWithNodes.elems[iedge]);
681 break; // --- edge already created
684 // --- if edge is not registered in the faces, create edge
689 connEdgeId = _downArray[vtkEdgeType]->addCell();
690 downEdge->setNodes(connEdgeId, edgesWithNodes.elems[iedge].nodeIds);
694 // --- mark faces in downward edge and mark edge in downward faces
697 for (int idf = 0; idf < nbDownFaces; idf++)
699 int faceId2 = downFaces[idf];
700 int faceType = downTypes[idf];
701 //ASSERT(_downArray[faceType]);
702 //SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
703 _downArray[vtkEdgeType]->addUpCell(connEdgeId, faceId2, faceType);
704 _downArray[faceType]->addDownCell(faceId2, connEdgeId, vtkEdgeType);
705 // MESSAGE(" From face t:" << vtkFaceType << " " << faceId <<
706 // " edge " << connEdgeId << " belongs to face t:" << faceType << " " << faceId2);
712 CHRONOSTOP(23);CHRONO(24);
714 // compact downward connectivity structure: adjust downward arrays size, replace vector<vector int>> by a single vector<int>
715 // 3D first then 2D and last 1D to release memory before edge upCells reorganization, (temporary memory use)
717 for (int vtkType = VTK_QUADRATIC_PYRAMID; vtkType >= 0; vtkType--)
719 if (SMDS_Downward *down = _downArray[vtkType])
721 down->compactStorage();
727 for (int vtkType = 0; vtkType <= VTK_QUADRATIC_PYRAMID; vtkType++)
729 if (SMDS_Downward *down = _downArray[vtkType])
731 if (down->getMaxId())
733 MESSAGE("Cells of Type " << vtkType << " : number of entities, est: "
734 << GuessSize[vtkType] << " real: " << down->getMaxId());
737 }CHRONOSTOP(24);CHRONOSTOP(2);
741 /*! Get the neighbors of a cell.
742 * Only the neighbors having the dimension of the cell are taken into account
743 * (neighbors of a volume are the volumes sharing a face with this volume,
744 * neighbors of a face are the faces sharing an edge with this face...).
745 * @param neighborsVtkIds vector of neighbors vtk id's to fill (reserve enough space).
746 * @param downIds downward id's of cells of dimension n-1, to fill (reserve enough space).
747 * @param downTypes vtk types of cells of dimension n-1, to fill (reserve enough space).
748 * @param vtkId the vtk id of the cell
749 * @return number of neighbors
751 int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId, bool getSkin)
753 int vtkType = this->GetCellType(vtkId);
754 int cellDim = SMDS_Downward::getCellDimension(vtkType);
756 return 0; // TODO voisins des edges = edges connectees
757 int cellId = this->_cellIdToDownId[vtkId];
759 int nbCells = _downArray[vtkType]->getNumberOfDownCells(cellId);
760 const int *downCells = _downArray[vtkType]->getDownCells(cellId);
761 const unsigned char* downTyp = _downArray[vtkType]->getDownTypes(cellId);
763 // --- iteration on faces of the 3D cell (or edges on the 2D cell).
766 for (int i = 0; i < nbCells; i++)
768 int downId = downCells[i];
769 int cellType = downTyp[i];
770 int nbUp = _downArray[cellType]->getNumberOfUpCells(downId);
771 const int *upCells = _downArray[cellType]->getUpCells(downId);
772 const unsigned char* upTypes = _downArray[cellType]->getUpTypes(downId);
774 // ---for a volume, max 2 upCells, one is this cell, the other is a neighbor
775 // for a face, number of neighbors (connected faces) not known
777 for (int j = 0; j < nbUp; j++)
779 if ((upCells[j] == cellId) && (upTypes[j] == vtkType))
781 int vtkNeighbor = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
782 neighborsVtkIds[nb] = vtkNeighbor;
783 downIds[nb] = downId;
784 downTypes[nb] = cellType;
786 if (nb >= NBMAXNEIGHBORS)
788 INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
794 if (cellDim == 3 && nbUp == 1) // this face is on the skin of the volume
796 neighborsVtkIds[nb] = _downArray[cellType]->getVtkCellId(downId); // OK if skin present
797 downIds[nb] = downId;
798 downTypes[nb] = cellType;
800 if (nb >= NBMAXNEIGHBORS)
802 INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
811 /*! get the volumes containing a face or an edge of the grid
812 * The edge or face belongs to the vtkUnstructuredGrid
813 * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
814 * @param vtkId vtk id of the face or edge
816 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int vtkId)
818 int vtkType = this->GetCellType(vtkId);
819 int dim = SMDS_Downward::getCellDimension(vtkType);
821 unsigned char cellTypes[1000];
822 int downCellId[1000];
825 int downId = this->CellIdToDownId(vtkId);
828 MESSAGE("Downward structure not up to date: new edge not taken into account");
831 nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
832 const int *upCells = _downArray[vtkType]->getUpCells(downId);
833 const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
834 for (int i=0; i< nbFaces; i++)
836 cellTypes[i] = upTypes[i];
837 downCellId[i] = upCells[i];
843 cellTypes[0] = this->GetCellType(vtkId);
844 int downId = this->CellIdToDownId(vtkId);
847 MESSAGE("Downward structure not up to date: new face not taken into account");
850 downCellId[0] = downId;
854 for (int i=0; i<nbFaces; i++)
856 int vtkTypeFace = cellTypes[i];
857 int downId = downCellId[i];
858 int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
859 const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
860 const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
861 for (int j=0; j<nv; j++)
863 int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
865 volVtkIds[nbvol++] = vtkVolId;
871 /*! get the volumes containing a face or an edge of the downward structure
872 * The edge or face does not necessary belong to the vtkUnstructuredGrid
873 * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
874 * @param downId id in the downward structure
875 * @param downType type of cell
877 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int downId, unsigned char downType)
879 int vtkType = downType;
880 int dim = SMDS_Downward::getCellDimension(vtkType);
882 unsigned char cellTypes[1000];
883 int downCellId[1000];
886 nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
887 const int *upCells = _downArray[vtkType]->getUpCells(downId);
888 const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
889 for (int i=0; i< nbFaces; i++)
891 cellTypes[i] = upTypes[i];
892 downCellId[i] = upCells[i];
898 cellTypes[0] = vtkType;
899 downCellId[0] = downId;
903 for (int i=0; i<nbFaces; i++)
905 int vtkTypeFace = cellTypes[i];
906 int downId = downCellId[i];
907 int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
908 const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
909 const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
910 for (int j=0; j<nv; j++)
912 int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
914 volVtkIds[nbvol++] = vtkVolId;
920 /*! get the node id's of a cell.
921 * The cell is defined by it's downward connectivity id and type.
922 * @param nodeSet set of of vtk node id's to fill.
923 * @param downId downward connectivity id of the cell.
924 * @param downType type of cell.
926 void SMDS_UnstructuredGrid::GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType)
928 _downArray[downType]->getNodeIds(downId, nodeSet);
931 /*! change some nodes in cell without modifying type or internal connectivity.
932 * Nodes inverse connectivity is maintained up to date.
933 * @param vtkVolId vtk id of the cell
934 * @param localClonedNodeIds map old node id to new node id.
936 void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds)
939 vtkIdType *pts; // will refer to the point id's of the face
940 this->GetCellPoints(vtkVolId, npts, pts);
941 for (int i = 0; i < npts; i++)
943 if (localClonedNodeIds.count(pts[i]))
945 vtkIdType oldpt = pts[i];
946 pts[i] = localClonedNodeIds[oldpt];
947 //MESSAGE(oldpt << " --> " << pts[i]);
948 //this->RemoveReferenceToCell(oldpt, vtkVolId);
949 //this->AddReferenceToCell(pts[i], vtkVolId);
954 /*! reorder the nodes of a face
955 * @param vtkVolId vtk id of a volume containing the face, to get an orientation for the face.
956 * @param orderedNodes list of nodes to reorder (in out)
957 * @return size of the list
959 int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes)
961 int vtkType = this->GetCellType(vtkVolId);
962 dim = SMDS_Downward::getCellDimension(vtkType);
965 SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
966 int downVolId = this->_cellIdToDownId[vtkVolId];
967 downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
969 // else nothing to do;
970 return orderedNodes.size();
973 void SMDS_UnstructuredGrid::BuildLinks()
975 // Remove the old links if they are already built
978 this->Links->UnRegister(this);
981 this->Links = SMDS_CellLinks::New();
982 this->Links->Allocate(this->GetNumberOfPoints());
983 this->Links->Register(this);
984 this->Links->BuildLinks(this, this->Connectivity);
985 this->Links->Delete();
988 /*! Create a volume (prism or hexahedron) by duplication of a face.
989 * Designed for use in creation of flat elements separating volume domains.
990 * A face separating two domains is shared by two volume cells.
991 * All the nodes are already created (for the two faces).
992 * Each original Node is associated to corresponding nodes in the domains.
993 * Some nodes may be duplicated for more than two domains, when domain separations intersect.
994 * In that case, even some of the nodes to use for the original face may be changed.
995 * @param vtkVolId: vtk id of a volume containing the face, to get an orientation for the face.
996 * @param domain1: domain of the original face
997 * @param domain2: domain of the duplicated face
998 * @param originalNodes: the vtk node ids of the original face
999 * @param nodeDomains: map(original id --> map(domain --> duplicated node id))
1000 * @return ok if success.
1003 SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
1006 std::set<int>& originalNodes,
1007 std::map<int, std::map<int, int> >& nodeDomains,
1008 std::map<int, std::map<long, int> >& nodeQuadDomains)
1010 //MESSAGE("extrudeVolumeFromFace " << vtkVolId);
1011 vector<vtkIdType> orderedOriginals;
1012 orderedOriginals.clear();
1013 set<int>::const_iterator it = originalNodes.begin();
1014 for (; it != originalNodes.end(); ++it)
1015 orderedOriginals.push_back(*it);
1018 int nbNodes = this->getOrderedNodesOfFace(vtkVolId, dim, orderedOriginals);
1019 vector<vtkIdType> orderedNodes;
1021 bool isQuadratic = false;
1022 switch (orderedOriginals.size())
1033 isQuadratic = false;
1039 long dom1 = domain1;
1040 long dom2 = domain2;
1041 long dom1_2; // for nodeQuadDomains
1042 if (domain1 < domain2)
1043 dom1_2 = dom1 + INT_MAX * dom2;
1045 dom1_2 = dom2 + INT_MAX * dom1;
1046 //cerr << "dom1=" << dom1 << " dom2=" << dom2 << " dom1_2=" << dom1_2 << endl;
1047 int ima = orderedOriginals.size();
1048 int mid = orderedOriginals.size() / 2;
1049 //cerr << "ima=" << ima << " mid=" << mid << endl;
1050 for (int i = 0; i < mid; i++)
1051 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1052 for (int i = 0; i < mid; i++)
1053 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1054 for (int i = mid; i < ima; i++)
1055 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1056 for (int i = mid; i < ima; i++)
1057 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1058 for (int i = 0; i < mid; i++)
1060 int oldId = orderedOriginals[i];
1062 if (nodeQuadDomains.count(oldId) && nodeQuadDomains[oldId].count(dom1_2))
1063 newId = nodeQuadDomains[oldId][dom1_2];
1066 double *coords = this->GetPoint(oldId);
1067 SMDS_MeshNode *newNode = _mesh->AddNode(coords[0], coords[1], coords[2]);
1068 newId = newNode->getVtkId();
1069 if (! nodeQuadDomains.count(oldId))
1071 std::map<long, int> emptyMap;
1072 nodeQuadDomains[oldId] = emptyMap;
1074 nodeQuadDomains[oldId][dom1_2] = newId;
1076 orderedNodes.push_back(newId);
1081 for (int i = 0; i < nbNodes; i++)
1082 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1084 for (int i = 0; i < nbNodes; i++)
1085 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1087 for (int i = nbNodes-1; i >= 0; i--)
1088 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1094 SMDS_MeshVolume *vol = _mesh->AddVolumeFromVtkIds(orderedNodes);
1099 SMDS_MeshFace *face = _mesh->AddFaceFromVtkIds(orderedNodes);
1103 // TODO update sub-shape list of elements and nodes
1107 //================================================================================
1109 * \brief Allocates data array for ball diameters
1110 * \param MaxVtkID - max ID of a ball element
1112 //================================================================================
1114 void SMDS_UnstructuredGrid::AllocateDiameters( vtkIdType MaxVtkID )
1116 SetBallDiameter( MaxVtkID, 0 );
1119 //================================================================================
1121 * \brief Sets diameter of a ball element
1122 * \param vtkID - vtk id of the ball element
1123 * \param diameter - diameter of the ball element
1125 //================================================================================
1127 void SMDS_UnstructuredGrid::SetBallDiameter( vtkIdType vtkID, double diameter )
1129 vtkDoubleArray* array = vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() );
1132 array = vtkDoubleArray::New();
1133 array->SetNumberOfComponents(1);
1134 vtkDataSet::CellData->SetScalars( array );
1136 array->InsertValue( vtkID, diameter );
1139 //================================================================================
1141 * \brief Returns diameter of a ball element
1142 * \param vtkID - vtk id of the ball element
1144 //================================================================================
1146 double SMDS_UnstructuredGrid::GetBallDiameter( vtkIdType vtkID ) const
1148 if ( vtkDataSet::CellData )
1149 return vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )->GetValue( vtkID );