1 // Copyright (C) 2010-2015 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 MESSAGE("SMDS_CellLinks::New");
44 return new SMDS_CellLinks();
47 void SMDS_CellLinks::ResizeForPoint(vtkIdType vtkID)
49 if ( vtkID > this->MaxId )
52 if ( vtkID >= this->Size )
53 vtkCellLinks::Resize( vtkID+SMDS_Mesh::chunkSize );
57 SMDS_CellLinks::SMDS_CellLinks() :
62 SMDS_CellLinks::~SMDS_CellLinks()
66 SMDS_UnstructuredGrid* SMDS_UnstructuredGrid::New()
68 MESSAGE("SMDS_UnstructuredGrid::New");
69 return new SMDS_UnstructuredGrid();
72 SMDS_UnstructuredGrid::SMDS_UnstructuredGrid() :
75 _cellIdToDownId.clear();
81 SMDS_UnstructuredGrid::~SMDS_UnstructuredGrid()
85 unsigned long SMDS_UnstructuredGrid::GetMTime()
87 unsigned long mtime = vtkUnstructuredGrid::GetMTime();
88 MESSAGE("vtkUnstructuredGrid::GetMTime: " << mtime);
91 // OUV_PORTING_VTK6: seems to be useless
93 void SMDS_UnstructuredGrid::Update()
95 MESSAGE("SMDS_UnstructuredGrid::Update");
96 return vtkUnstructuredGrid::Update();
99 void SMDS_UnstructuredGrid::UpdateInformation()
101 MESSAGE("SMDS_UnstructuredGrid::UpdateInformation");
102 return vtkUnstructuredGrid::UpdateInformation();
105 vtkPoints* SMDS_UnstructuredGrid::GetPoints()
107 // TODO erreur incomprehensible de la macro vtk GetPoints apparue avec la version paraview de fin aout 2010
108 //MESSAGE("*********************** SMDS_UnstructuredGrid::GetPoints " << this->Points << " " << vtkUnstructuredGrid::GetPoints());
112 //#ifdef VTK_HAVE_POLYHEDRON
113 int SMDS_UnstructuredGrid::InsertNextLinkedCell(int type, int npts, vtkIdType *pts)
115 if (type != VTK_POLYHEDRON)
116 return vtkUnstructuredGrid::InsertNextLinkedCell(type, npts, pts);
118 // --- type = VTK_POLYHEDRON
119 //MESSAGE("InsertNextLinkedCell VTK_POLYHEDRON");
120 int cellid = this->InsertNextCell(type, npts, pts);
122 set<vtkIdType> setOfNodes;
126 for (int nf = 0; nf < nbfaces; nf++)
128 int nbnodes = pts[i];
130 for (int k = 0; k < nbnodes; k++)
132 //MESSAGE(" cell " << cellid << " face " << nf << " node " << pts[i]);
133 setOfNodes.insert(pts[i]);
138 set<vtkIdType>::iterator it = setOfNodes.begin();
139 for (; it != setOfNodes.end(); ++it)
141 //MESSAGE("reverse link for node " << *it << " cell " << cellid);
142 this->Links->ResizeCellList(*it, 1);
143 this->Links->AddCellReference(cellid, *it);
150 void SMDS_UnstructuredGrid::setSMDS_mesh(SMDS_Mesh *mesh)
155 void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize,
156 std::vector<int>& idCellsOldToNew, int newCellSize)
158 MESSAGE("------------------------- SMDS_UnstructuredGrid::compactGrid " << newNodeSize << " " << newCellSize);//CHRONO(1);
159 int alreadyCopied = 0;
161 // --- if newNodeSize, create a new compacted vtkPoints
163 vtkPoints *newPoints = vtkPoints::New();
164 newPoints->SetDataType(VTK_DOUBLE);
165 newPoints->SetNumberOfPoints(newNodeSize);
168 MESSAGE("-------------- compactGrid, newNodeSize " << newNodeSize);
169 // rnv: to fix bug "21125: EDF 1233 SMESH: Degradation of precision in a test case for quadratic conversion"
170 // using double type for storing coordinates of nodes instead float.
171 int oldNodeSize = idNodesOldToNew.size();
174 while ( i < oldNodeSize )
176 // skip a hole if any
177 while ( i < oldNodeSize && idNodesOldToNew[i] < 0 )
180 // look for a block end
181 while ( i < oldNodeSize && idNodesOldToNew[i] >= 0 )
184 copyNodes(newPoints, idNodesOldToNew, alreadyCopied, startBloc, endBloc);
186 newPoints->Squeeze();
189 // --- create new compacted Connectivity, Locations and Types
191 int oldCellSize = this->Types->GetNumberOfTuples();
193 vtkCellArray *newConnectivity = vtkCellArray::New();
194 newConnectivity->Initialize();
195 int oldCellDataSize = this->Connectivity->GetData()->GetSize();
196 newConnectivity->Allocate(oldCellDataSize);
197 MESSAGE("oldCellSize="<< oldCellSize << " oldCellDataSize=" << oldCellDataSize);
199 vtkUnsignedCharArray *newTypes = vtkUnsignedCharArray::New();
200 newTypes->Initialize();
201 newTypes->SetNumberOfValues(newCellSize);
203 vtkIdTypeArray *newLocations = vtkIdTypeArray::New();
204 newLocations->Initialize();
205 newLocations->SetNumberOfValues(newCellSize);
207 // TODO some polyhedron may be huge (only in some tests)
208 vtkIdType tmpid[NBMAXNODESINCELL];
209 vtkIdType *pointsCell = &tmpid[0]; // --- points id to fill a new cell
213 while ( i < oldCellSize )
215 // skip a hole if any
216 while ( i < oldCellSize && this->Types->GetValue(i) == VTK_EMPTY_CELL )
219 // look for a block end
220 while ( i < oldCellSize && this->Types->GetValue(i) != VTK_EMPTY_CELL )
223 if ( endBloc > startBloc )
225 idCellsOldToNew, idNodesOldToNew,
226 newConnectivity, newLocations,
227 pointsCell, alreadyCopied,
230 newConnectivity->Squeeze();
232 if (1/*newNodeSize*/)
234 MESSAGE("------- newNodeSize, setPoints");
235 this->SetPoints(newPoints);
236 MESSAGE("NumberOfPoints: " << this->GetNumberOfPoints());
239 if (vtkDoubleArray* diameters =
240 vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )) // Balls
242 for (int oldCellID = 0; oldCellID < oldCellSize; oldCellID++)
244 if (this->Types->GetValue(oldCellID) == VTK_EMPTY_CELL)
246 int newCellId = idCellsOldToNew[ oldCellID ];
247 if (newTypes->GetValue(newCellId) == VTK_POLY_VERTEX)
248 diameters->SetValue( newCellId, diameters->GetValue( oldCellID ));
252 if (this->FaceLocations)
254 vtkIdTypeArray *newFaceLocations = vtkIdTypeArray::New();
255 newFaceLocations->Initialize();
256 newFaceLocations->Allocate(newTypes->GetSize());
257 vtkIdTypeArray *newFaces = vtkIdTypeArray::New();
258 newFaces->Initialize();
259 newFaces->Allocate(this->Faces->GetSize());
260 for (int i = 0; i < oldCellSize; i++)
262 if (this->Types->GetValue(i) == VTK_EMPTY_CELL)
264 int newCellId = idCellsOldToNew[i];
265 if (newTypes->GetValue(newCellId) == VTK_POLYHEDRON)
267 newFaceLocations->InsertNextValue(newFaces->GetMaxId()+1);
268 int oldFaceLoc = this->FaceLocations->GetValue(i);
269 int nCellFaces = this->Faces->GetValue(oldFaceLoc++);
270 newFaces->InsertNextValue(nCellFaces);
271 for (int n=0; n<nCellFaces; n++)
273 int nptsInFace = this->Faces->GetValue(oldFaceLoc++);
274 newFaces->InsertNextValue(nptsInFace);
275 for (int k=0; k<nptsInFace; k++)
277 int oldpt = this->Faces->GetValue(oldFaceLoc++);
278 newFaces->InsertNextValue(idNodesOldToNew[oldpt]);
284 newFaceLocations->InsertNextValue(-1);
287 newFaceLocations->Squeeze();
289 this->SetCells(newTypes, newLocations, newConnectivity, newFaceLocations, newFaces);
290 newFaceLocations->Delete();
295 this->SetCells(newTypes, newLocations, newConnectivity, FaceLocations, Faces);
300 newLocations->Delete();
301 newConnectivity->Delete();
305 void SMDS_UnstructuredGrid::copyNodes(vtkPoints *newPoints, std::vector<int>& idNodesOldToNew, int& alreadyCopied,
308 MESSAGE("copyNodes " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
309 void *target = newPoints->GetVoidPointer(3 * alreadyCopied);
310 void *source = this->Points->GetVoidPointer(3 * start);
311 int nbPoints = end - start;
314 memcpy(target, source, 3 * sizeof(double) * nbPoints);
315 for (int j = start; j < end; j++)
316 idNodesOldToNew[j] = alreadyCopied++; // old vtkId --> new vtkId
320 void SMDS_UnstructuredGrid::copyBloc(vtkUnsignedCharArray *newTypes,
321 std::vector<int>& idCellsOldToNew,
322 std::vector<int>& idNodesOldToNew,
323 vtkCellArray* newConnectivity,
324 vtkIdTypeArray* newLocations,
325 vtkIdType* pointsCell,
330 //MESSAGE("copyBloc " << alreadyCopied << " " << start << " " << end << " size: " << end - start << " total: " << alreadyCopied + end - start);
331 for (int j = start; j < end; j++)
333 newTypes->SetValue(alreadyCopied, this->Types->GetValue(j));
334 idCellsOldToNew[j] = alreadyCopied; // old vtkId --> new vtkId
335 vtkIdType oldLoc = this->Locations->GetValue(j);
337 vtkIdType *oldPtsCell = 0;
338 this->Connectivity->GetCell(oldLoc, nbpts, oldPtsCell);
339 assert(nbpts < NBMAXNODESINCELL);
340 //MESSAGE(j << " " << alreadyCopied << " " << (int)this->Types->GetValue(j) << " " << oldLoc << " " << nbpts );
341 for (int l = 0; l < nbpts; l++)
343 int oldval = oldPtsCell[l];
344 pointsCell[l] = idNodesOldToNew[oldval];
345 //MESSAGE(" " << oldval << " " << pointsCell[l]);
347 /*int newcnt = */newConnectivity->InsertNextCell(nbpts, pointsCell);
348 int newLoc = newConnectivity->GetInsertLocation(nbpts);
349 //MESSAGE(newcnt << " " << newLoc);
350 newLocations->SetValue(alreadyCopied, newLoc);
355 int SMDS_UnstructuredGrid::CellIdToDownId(int vtkCellId)
357 if ((vtkCellId < 0) || (vtkCellId >= (int)_cellIdToDownId.size()))
359 //MESSAGE("SMDS_UnstructuredGrid::CellIdToDownId structure not up to date: vtkCellId="
360 // << vtkCellId << " max="<< _cellIdToDownId.size());
363 return _cellIdToDownId[vtkCellId];
366 void SMDS_UnstructuredGrid::setCellIdToDownId(int vtkCellId, int downId)
368 // ASSERT((vtkCellId >= 0) && (vtkCellId < _cellIdToDownId.size()));
369 _cellIdToDownId[vtkCellId] = downId;
372 void SMDS_UnstructuredGrid::CleanDownwardConnectivity()
374 for (size_t i = 0; i < _downArray.size(); i++)
377 delete _downArray[i];
380 _cellIdToDownId.clear();
383 /*! Build downward connectivity: to do only when needed because heavy memory load.
384 * Downward connectivity is no more valid if vtkUnstructuredGrid is modified.
387 void SMDS_UnstructuredGrid::BuildDownwardConnectivity(bool withEdges)
389 MESSAGE("SMDS_UnstructuredGrid::BuildDownwardConnectivity");CHRONO(2);
390 // TODO calcul partiel sans edges
392 // --- erase previous data if any
394 this->CleanDownwardConnectivity();
396 // --- create SMDS_Downward structures (in _downArray vector[vtkCellType])
398 _downArray.resize(VTK_MAXTYPE + 1, 0);
400 _downArray[VTK_LINE] = new SMDS_DownEdge(this);
401 _downArray[VTK_QUADRATIC_EDGE] = new SMDS_DownQuadEdge(this);
402 _downArray[VTK_TRIANGLE] = new SMDS_DownTriangle(this);
403 _downArray[VTK_QUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
404 _downArray[VTK_BIQUADRATIC_TRIANGLE] = new SMDS_DownQuadTriangle(this);
405 _downArray[VTK_QUAD] = new SMDS_DownQuadrangle(this);
406 _downArray[VTK_QUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
407 _downArray[VTK_BIQUADRATIC_QUAD] = new SMDS_DownQuadQuadrangle(this);
408 _downArray[VTK_TETRA] = new SMDS_DownTetra(this);
409 _downArray[VTK_QUADRATIC_TETRA] = new SMDS_DownQuadTetra(this);
410 _downArray[VTK_PYRAMID] = new SMDS_DownPyramid(this);
411 _downArray[VTK_QUADRATIC_PYRAMID] = new SMDS_DownQuadPyramid(this);
412 _downArray[VTK_WEDGE] = new SMDS_DownPenta(this);
413 _downArray[VTK_QUADRATIC_WEDGE] = new SMDS_DownQuadPenta(this);
414 _downArray[VTK_HEXAHEDRON] = new SMDS_DownHexa(this);
415 _downArray[VTK_QUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
416 _downArray[VTK_TRIQUADRATIC_HEXAHEDRON] = new SMDS_DownQuadHexa(this);
417 _downArray[VTK_HEXAGONAL_PRISM] = new SMDS_DownPenta(this);
419 // --- get detailed info of number of cells of each type, allocate SMDS_downward structures
421 const SMDS_MeshInfo &meshInfo = _mesh->GetMeshInfo();
423 int nbLinTetra = meshInfo.NbTetras (ORDER_LINEAR);
424 int nbQuadTetra = meshInfo.NbTetras (ORDER_QUADRATIC);
425 int nbLinPyra = meshInfo.NbPyramids(ORDER_LINEAR);
426 int nbQuadPyra = meshInfo.NbPyramids(ORDER_QUADRATIC);
427 int nbLinPrism = meshInfo.NbPrisms (ORDER_LINEAR);
428 int nbQuadPrism = meshInfo.NbPrisms (ORDER_QUADRATIC);
429 int nbLinHexa = meshInfo.NbHexas (ORDER_LINEAR);
430 int nbQuadHexa = meshInfo.NbHexas (ORDER_QUADRATIC);
431 int nbHexPrism = meshInfo.NbHexPrisms();
433 int nbLineGuess = int((4.0 / 3.0) * nbLinTetra + 2 * nbLinPrism + 2.5 * nbLinPyra + 3 * nbLinHexa);
434 int nbQuadEdgeGuess = int((4.0 / 3.0) * nbQuadTetra + 2 * nbQuadPrism + 2.5 * nbQuadPyra + 3 * nbQuadHexa);
435 int nbLinTriaGuess = 2 * nbLinTetra + nbLinPrism + 2 * nbLinPyra;
436 int nbQuadTriaGuess = 2 * nbQuadTetra + nbQuadPrism + 2 * nbQuadPyra;
437 int nbLinQuadGuess = int((2.0 / 3.0) * nbLinPrism + (1.0 / 2.0) * nbLinPyra + 3 * nbLinHexa);
438 int nbQuadQuadGuess = int((2.0 / 3.0) * nbQuadPrism + (1.0 / 2.0) * nbQuadPyra + 3 * nbQuadHexa);
440 int GuessSize[VTK_MAXTYPE];
441 GuessSize[VTK_LINE] = nbLineGuess;
442 GuessSize[VTK_QUADRATIC_EDGE] = nbQuadEdgeGuess;
443 GuessSize[VTK_TRIANGLE] = nbLinTriaGuess;
444 GuessSize[VTK_QUADRATIC_TRIANGLE] = nbQuadTriaGuess;
445 GuessSize[VTK_BIQUADRATIC_TRIANGLE] = nbQuadTriaGuess;
446 GuessSize[VTK_QUAD] = nbLinQuadGuess;
447 GuessSize[VTK_QUADRATIC_QUAD] = nbQuadQuadGuess;
448 GuessSize[VTK_BIQUADRATIC_QUAD] = nbQuadQuadGuess;
449 GuessSize[VTK_TETRA] = nbLinTetra;
450 GuessSize[VTK_QUADRATIC_TETRA] = nbQuadTetra;
451 GuessSize[VTK_PYRAMID] = nbLinPyra;
452 GuessSize[VTK_QUADRATIC_PYRAMID] = nbQuadPyra;
453 GuessSize[VTK_WEDGE] = nbLinPrism;
454 GuessSize[VTK_QUADRATIC_WEDGE] = nbQuadPrism;
455 GuessSize[VTK_HEXAHEDRON] = nbLinHexa;
456 GuessSize[VTK_QUADRATIC_HEXAHEDRON] = nbQuadHexa;
457 GuessSize[VTK_TRIQUADRATIC_HEXAHEDRON] = nbQuadHexa;
458 GuessSize[VTK_HEXAGONAL_PRISM] = nbHexPrism;
460 _downArray[VTK_LINE] ->allocate(nbLineGuess);
461 _downArray[VTK_QUADRATIC_EDGE] ->allocate(nbQuadEdgeGuess);
462 _downArray[VTK_TRIANGLE] ->allocate(nbLinTriaGuess);
463 _downArray[VTK_QUADRATIC_TRIANGLE] ->allocate(nbQuadTriaGuess);
464 _downArray[VTK_BIQUADRATIC_TRIANGLE] ->allocate(nbQuadTriaGuess);
465 _downArray[VTK_QUAD] ->allocate(nbLinQuadGuess);
466 _downArray[VTK_QUADRATIC_QUAD] ->allocate(nbQuadQuadGuess);
467 _downArray[VTK_BIQUADRATIC_QUAD] ->allocate(nbQuadQuadGuess);
468 _downArray[VTK_TETRA] ->allocate(nbLinTetra);
469 _downArray[VTK_QUADRATIC_TETRA] ->allocate(nbQuadTetra);
470 _downArray[VTK_PYRAMID] ->allocate(nbLinPyra);
471 _downArray[VTK_QUADRATIC_PYRAMID] ->allocate(nbQuadPyra);
472 _downArray[VTK_WEDGE] ->allocate(nbLinPrism);
473 _downArray[VTK_QUADRATIC_WEDGE] ->allocate(nbQuadPrism);
474 _downArray[VTK_HEXAHEDRON] ->allocate(nbLinHexa);
475 _downArray[VTK_QUADRATIC_HEXAHEDRON] ->allocate(nbQuadHexa);
476 _downArray[VTK_TRIQUADRATIC_HEXAHEDRON]->allocate(nbQuadHexa);
477 _downArray[VTK_HEXAGONAL_PRISM] ->allocate(nbHexPrism);
479 // --- iteration on vtkUnstructuredGrid cells, only faces
480 // for each vtk face:
481 // create a downward face entry with its downward id.
482 // compute vtk volumes, create downward volumes entry.
483 // mark face in downward volumes
484 // mark volumes in downward face
486 MESSAGE("--- iteration on vtkUnstructuredGrid cells, only faces");CHRONO(20);
487 int cellSize = this->Types->GetNumberOfTuples();
488 _cellIdToDownId.resize(cellSize, -1);
490 for (int i = 0; i < cellSize; i++)
492 int vtkFaceType = this->GetCellType(i);
493 if (SMDS_Downward::getCellDimension(vtkFaceType) == 2)
496 //ASSERT(_downArray[vtkFaceType]);
497 int connFaceId = _downArray[vtkFaceType]->addCell(vtkFaceId);
498 SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
499 downFace->setTempNodes(connFaceId, vtkFaceId);
500 int vols[2] = { -1, -1 };
501 int nbVolumes = downFace->computeVolumeIds(vtkFaceId, vols);
502 //MESSAGE("nbVolumes="<< nbVolumes);
503 for (int ivol = 0; ivol < nbVolumes; ivol++)
505 int vtkVolId = vols[ivol];
506 int vtkVolType = this->GetCellType(vtkVolId);
507 //ASSERT(_downArray[vtkVolType]);
508 int connVolId = _downArray[vtkVolType]->addCell(vtkVolId);
509 _downArray[vtkVolType]->addDownCell(connVolId, connFaceId, vtkFaceType);
510 _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId, vtkVolType);
511 // MESSAGE("Face " << vtkFaceId << " belongs to volume " << vtkVolId);
516 // --- iteration on vtkUnstructuredGrid cells, only volumes
517 // for each vtk volume:
518 // create downward volumes entry if not already done
519 // build a temporary list of faces described with their nodes
521 // compute the vtk volumes containing this face
522 // check if the face is already listed in the volumes (comparison of ordered list of nodes)
523 // if not, create a downward face entry (resizing of structure required)
524 // (the downward faces store a temporary list of nodes to ease the comparison)
525 // create downward volumes entry if not already done
526 // mark volumes in downward face
527 // mark face in downward volumes
530 MESSAGE("--- iteration on vtkUnstructuredGrid cells, only volumes");CHRONO(21);
532 for (int i = 0; i < cellSize; i++)
534 int vtkType = this->GetCellType(i);
535 if (SMDS_Downward::getCellDimension(vtkType) == 3)
539 // MESSAGE("vtk volume " << vtkVolId);
540 //ASSERT(_downArray[vtkType]);
541 /*int connVolId = */_downArray[vtkType]->addCell(vtkVolId);
543 // --- find all the faces of the volume, describe the faces by their nodes
545 SMDS_Down3D* downVol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
546 ListElemByNodesType facesWithNodes;
547 downVol->computeFacesWithNodes(vtkVolId, facesWithNodes);
548 // MESSAGE("vtk volume " << vtkVolId << " contains " << facesWithNodes.nbElems << " faces");
550 for (int iface = 0; iface < facesWithNodes.nbElems; iface++)
552 // --- find the volumes containing the face
555 int vtkFaceType = facesWithNodes.elems[iface].vtkType;
556 SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
557 int vols[2] = { -1, -1 };
558 int *nodes = &facesWithNodes.elems[iface].nodeIds[0];
559 int lg = facesWithNodes.elems[iface].nbNodes;
560 int nbVolumes = downFace->computeVolumeIdsFromNodesFace(nodes, lg, vols);
561 // MESSAGE("vtk volume " << vtkVolId << " face " << iface << " belongs to " << nbVolumes << " volumes");
563 // --- check if face is registered in the volumes
568 for (int ivol = 0; ivol < nbVolumes; ivol++)
570 int vtkVolId2 = vols[ivol];
571 int vtkVolType = this->GetCellType(vtkVolId2);
572 //ASSERT(_downArray[vtkVolType]);
573 int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
574 SMDS_Down3D* downVol2 = static_cast<SMDS_Down3D*> (_downArray[vtkVolType]);
575 connFaceId = downVol2->FindFaceByNodes(connVolId2, facesWithNodes.elems[iface]);
577 break; // --- face already created
580 // --- if face is not registered in the volumes, create face
585 connFaceId = _downArray[vtkFaceType]->addCell();
586 downFace->setTempNodes(connFaceId, facesWithNodes.elems[iface]);
589 // --- mark volumes in downward face and mark face in downward volumes
592 for (int ivol = 0; ivol < nbVolumes; ivol++)
594 int vtkVolId2 = vols[ivol];
595 int vtkVolType = this->GetCellType(vtkVolId2);
596 //ASSERT(_downArray[vtkVolType]);
597 int connVolId2 = _downArray[vtkVolType]->addCell(vtkVolId2);
598 _downArray[vtkVolType]->addDownCell(connVolId2, connFaceId, vtkFaceType);
599 _downArray[vtkFaceType]->addUpCell(connFaceId, connVolId2, vtkVolType);
600 // MESSAGE(" From volume " << vtkVolId << " face " << connFaceId << " belongs to volume " << vtkVolId2);
606 // --- iteration on vtkUnstructuredGrid cells, only edges
607 // for each vtk edge:
608 // create downward edge entry
609 // store the nodes id's in downward edge (redundant with vtkUnstructuredGrid)
610 // find downward faces
611 // (from vtk faces or volumes, get downward faces, they have a temporary list of nodes)
612 // mark edge in downward faces
613 // mark faces in downward edge
616 MESSAGE("--- iteration on vtkUnstructuredGrid cells, only edges");CHRONO(22);
618 for (int i = 0; i < cellSize; i++)
620 int vtkEdgeType = this->GetCellType(i);
621 if (SMDS_Downward::getCellDimension(vtkEdgeType) == 1)
624 //ASSERT(_downArray[vtkEdgeType]);
625 int connEdgeId = _downArray[vtkEdgeType]->addCell(vtkEdgeId);
626 SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
627 downEdge->setNodes(connEdgeId, vtkEdgeId);
629 int nbVtkCells = downEdge->computeVtkCells(connEdgeId, vtkIds);
631 unsigned char downTypes[1000];
632 int nbDownFaces = downEdge->computeFaces(connEdgeId, &vtkIds[0], nbVtkCells, downFaces, downTypes);
633 for (int n = 0; n < nbDownFaces; n++)
635 _downArray[downTypes[n]]->addDownCell(downFaces[n], connEdgeId, vtkEdgeType);
636 _downArray[vtkEdgeType]->addUpCell(connEdgeId, downFaces[n], downTypes[n]);
641 // --- iteration on downward faces (they are all listed now)
642 // for each downward face:
643 // build a temporary list of edges with their ordered list of nodes
645 // find all the vtk cells containing this edge
646 // then identify all the downward faces containing the edge, from the vtk cells
647 // check if the edge is already listed in the faces (comparison of ordered list of nodes)
648 // if not, create a downward edge entry with the node id's
649 // mark edge in downward faces
650 // mark downward faces in edge (size of list unknown, to be allocated)
652 CHRONOSTOP(22);CHRONO(23);
654 for (int vtkFaceType = 0; vtkFaceType < VTK_QUADRATIC_PYRAMID; vtkFaceType++)
656 if (SMDS_Downward::getCellDimension(vtkFaceType) != 2)
659 // --- find all the edges of the face, describe the edges by their nodes
661 SMDS_Down2D* downFace = static_cast<SMDS_Down2D*> (_downArray[vtkFaceType]);
662 int maxId = downFace->getMaxId();
663 for (int faceId = 0; faceId < maxId; faceId++)
666 ListElemByNodesType edgesWithNodes;
667 downFace->computeEdgesWithNodes(faceId, edgesWithNodes);
668 // MESSAGE("downward face type " << vtkFaceType << " num " << faceId << " contains " << edgesWithNodes.nbElems << " edges");
671 for (int iedge = 0; iedge < edgesWithNodes.nbElems; iedge++)
674 // --- check if the edge is already registered by exploration of the faces
678 unsigned char vtkEdgeType = edgesWithNodes.elems[iedge].vtkType;
679 int *pts = &edgesWithNodes.elems[iedge].nodeIds[0];
680 SMDS_Down1D* downEdge = static_cast<SMDS_Down1D*> (_downArray[vtkEdgeType]);
681 int nbVtkCells = downEdge->computeVtkCells(pts, vtkIds);
682 //CHRONOSTOP(41);CHRONO(42);
684 unsigned char downTypes[1000];
685 int nbDownFaces = downEdge->computeFaces(pts, &vtkIds[0], nbVtkCells, downFaces, downTypes);
690 for (int idf = 0; idf < nbDownFaces; idf++)
692 int faceId2 = downFaces[idf];
693 int faceType = downTypes[idf];
694 //ASSERT(_downArray[faceType]);
695 SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
696 connEdgeId = downFace2->FindEdgeByNodes(faceId2, edgesWithNodes.elems[iedge]);
698 break; // --- edge already created
701 // --- if edge is not registered in the faces, create edge
706 connEdgeId = _downArray[vtkEdgeType]->addCell();
707 downEdge->setNodes(connEdgeId, edgesWithNodes.elems[iedge].nodeIds);
711 // --- mark faces in downward edge and mark edge in downward faces
714 for (int idf = 0; idf < nbDownFaces; idf++)
716 int faceId2 = downFaces[idf];
717 int faceType = downTypes[idf];
718 //ASSERT(_downArray[faceType]);
719 //SMDS_Down2D* downFace2 = static_cast<SMDS_Down2D*> (_downArray[faceType]);
720 _downArray[vtkEdgeType]->addUpCell(connEdgeId, faceId2, faceType);
721 _downArray[faceType]->addDownCell(faceId2, connEdgeId, vtkEdgeType);
722 // MESSAGE(" From face t:" << vtkFaceType << " " << faceId <<
723 // " edge " << connEdgeId << " belongs to face t:" << faceType << " " << faceId2);
729 CHRONOSTOP(23);CHRONO(24);
731 // compact downward connectivity structure: adjust downward arrays size, replace vector<vector int>> by a single vector<int>
732 // 3D first then 2D and last 1D to release memory before edge upCells reorganization, (temporary memory use)
734 for (int vtkType = VTK_QUADRATIC_PYRAMID; vtkType >= 0; vtkType--)
736 if (SMDS_Downward *down = _downArray[vtkType])
738 down->compactStorage();
744 for (int vtkType = 0; vtkType <= VTK_QUADRATIC_PYRAMID; vtkType++)
746 if (SMDS_Downward *down = _downArray[vtkType])
748 if (down->getMaxId())
750 MESSAGE("Cells of Type " << vtkType << " : number of entities, est: "
751 << GuessSize[vtkType] << " real: " << down->getMaxId());
754 }CHRONOSTOP(24);CHRONOSTOP(2);
758 /*! Get the neighbors of a cell.
759 * Only the neighbors having the dimension of the cell are taken into account
760 * (neighbors of a volume are the volumes sharing a face with this volume,
761 * neighbors of a face are the faces sharing an edge with this face...).
762 * @param neighborsVtkIds vector of neighbors vtk id's to fill (reserve enough space).
763 * @param downIds downward id's of cells of dimension n-1, to fill (reserve enough space).
764 * @param downTypes vtk types of cells of dimension n-1, to fill (reserve enough space).
765 * @param vtkId the vtk id of the cell
766 * @return number of neighbors
768 int SMDS_UnstructuredGrid::GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId, bool getSkin)
770 int vtkType = this->GetCellType(vtkId);
771 int cellDim = SMDS_Downward::getCellDimension(vtkType);
773 return 0; // TODO voisins des edges = edges connectees
774 int cellId = this->_cellIdToDownId[vtkId];
776 int nbCells = _downArray[vtkType]->getNumberOfDownCells(cellId);
777 const int *downCells = _downArray[vtkType]->getDownCells(cellId);
778 const unsigned char* downTyp = _downArray[vtkType]->getDownTypes(cellId);
780 // --- iteration on faces of the 3D cell (or edges on the 2D cell).
783 for (int i = 0; i < nbCells; i++)
785 int downId = downCells[i];
786 int cellType = downTyp[i];
787 int nbUp = _downArray[cellType]->getNumberOfUpCells(downId);
788 const int *upCells = _downArray[cellType]->getUpCells(downId);
789 const unsigned char* upTypes = _downArray[cellType]->getUpTypes(downId);
791 // ---for a volume, max 2 upCells, one is this cell, the other is a neighbor
792 // for a face, number of neighbors (connected faces) not known
794 for (int j = 0; j < nbUp; j++)
796 if ((upCells[j] == cellId) && (upTypes[j] == vtkType))
798 int vtkNeighbor = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
799 neighborsVtkIds[nb] = vtkNeighbor;
800 downIds[nb] = downId;
801 downTypes[nb] = cellType;
803 if (nb >= NBMAXNEIGHBORS)
805 INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
811 if (cellDim == 3 && nbUp == 1) // this face is on the skin of the volume
813 neighborsVtkIds[nb] = _downArray[cellType]->getVtkCellId(downId); // OK if skin present
814 downIds[nb] = downId;
815 downTypes[nb] = cellType;
817 if (nb >= NBMAXNEIGHBORS)
819 INFOS("SMDS_UnstructuredGrid::GetNeighbors problem: NBMAXNEIGHBORS=" <<NBMAXNEIGHBORS << " not enough");
828 /*! get the volumes containing a face or an edge of the grid
829 * The edge or face belongs to the vtkUnstructuredGrid
830 * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
831 * @param vtkId vtk id of the face or edge
833 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int vtkId)
835 int vtkType = this->GetCellType(vtkId);
836 int dim = SMDS_Downward::getCellDimension(vtkType);
838 unsigned char cellTypes[1000];
839 int downCellId[1000];
842 int downId = this->CellIdToDownId(vtkId);
845 MESSAGE("Downward structure not up to date: new edge not taken into account");
848 nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
849 const int *upCells = _downArray[vtkType]->getUpCells(downId);
850 const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
851 for (int i=0; i< nbFaces; i++)
853 cellTypes[i] = upTypes[i];
854 downCellId[i] = upCells[i];
860 cellTypes[0] = this->GetCellType(vtkId);
861 int downId = this->CellIdToDownId(vtkId);
864 MESSAGE("Downward structure not up to date: new face not taken into account");
867 downCellId[0] = downId;
871 for (int i=0; i<nbFaces; i++)
873 int vtkTypeFace = cellTypes[i];
874 int downId = downCellId[i];
875 int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
876 const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
877 const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
878 for (int j=0; j<nv; j++)
880 int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
882 volVtkIds[nbvol++] = vtkVolId;
888 /*! get the volumes containing a face or an edge of the downward structure
889 * The edge or face does not necessary belong to the vtkUnstructuredGrid
890 * @param volVtkIds vector of parent volume ids to fill (reserve enough space!)
891 * @param downId id in the downward structure
892 * @param downType type of cell
894 int SMDS_UnstructuredGrid::GetParentVolumes(int* volVtkIds, int downId, unsigned char downType)
896 int vtkType = downType;
897 int dim = SMDS_Downward::getCellDimension(vtkType);
899 unsigned char cellTypes[1000];
900 int downCellId[1000];
903 nbFaces = _downArray[vtkType]->getNumberOfUpCells(downId);
904 const int *upCells = _downArray[vtkType]->getUpCells(downId);
905 const unsigned char* upTypes = _downArray[vtkType]->getUpTypes(downId);
906 for (int i=0; i< nbFaces; i++)
908 cellTypes[i] = upTypes[i];
909 downCellId[i] = upCells[i];
915 cellTypes[0] = vtkType;
916 downCellId[0] = downId;
920 for (int i=0; i<nbFaces; i++)
922 int vtkTypeFace = cellTypes[i];
923 int downId = downCellId[i];
924 int nv = _downArray[vtkTypeFace]->getNumberOfUpCells(downId);
925 const int *upCells = _downArray[vtkTypeFace]->getUpCells(downId);
926 const unsigned char* upTypes = _downArray[vtkTypeFace]->getUpTypes(downId);
927 for (int j=0; j<nv; j++)
929 int vtkVolId = _downArray[upTypes[j]]->getVtkCellId(upCells[j]);
931 volVtkIds[nbvol++] = vtkVolId;
937 /*! get the node id's of a cell.
938 * The cell is defined by it's downward connectivity id and type.
939 * @param nodeSet set of of vtk node id's to fill.
940 * @param downId downward connectivity id of the cell.
941 * @param downType type of cell.
943 void SMDS_UnstructuredGrid::GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType)
945 _downArray[downType]->getNodeIds(downId, nodeSet);
948 /*! change some nodes in cell without modifying type or internal connectivity.
949 * Nodes inverse connectivity is maintained up to date.
950 * @param vtkVolId vtk id of the cell
951 * @param localClonedNodeIds map old node id to new node id.
953 void SMDS_UnstructuredGrid::ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds)
956 vtkIdType *pts; // will refer to the point id's of the face
957 this->GetCellPoints(vtkVolId, npts, pts);
958 for (int i = 0; i < npts; i++)
960 if (localClonedNodeIds.count(pts[i]))
962 vtkIdType oldpt = pts[i];
963 pts[i] = localClonedNodeIds[oldpt];
964 //MESSAGE(oldpt << " --> " << pts[i]);
965 //this->RemoveReferenceToCell(oldpt, vtkVolId);
966 //this->AddReferenceToCell(pts[i], vtkVolId);
971 /*! reorder the nodes of a face
972 * @param vtkVolId vtk id of a volume containing the face, to get an orientation for the face.
973 * @param orderedNodes list of nodes to reorder (in out)
974 * @return size of the list
976 int SMDS_UnstructuredGrid::getOrderedNodesOfFace(int vtkVolId, int& dim, std::vector<vtkIdType>& orderedNodes)
978 int vtkType = this->GetCellType(vtkVolId);
979 dim = SMDS_Downward::getCellDimension(vtkType);
982 SMDS_Down3D *downvol = static_cast<SMDS_Down3D*> (_downArray[vtkType]);
983 int downVolId = this->_cellIdToDownId[vtkVolId];
984 downvol->getOrderedNodesOfFace(downVolId, orderedNodes);
986 // else nothing to do;
987 return orderedNodes.size();
990 void SMDS_UnstructuredGrid::BuildLinks()
992 // Remove the old links if they are already built
995 this->Links->UnRegister(this);
998 this->Links = SMDS_CellLinks::New();
999 this->Links->Allocate(this->GetNumberOfPoints());
1000 this->Links->Register(this);
1001 this->Links->BuildLinks(this, this->Connectivity);
1002 this->Links->Delete();
1005 /*! Create a volume (prism or hexahedron) by duplication of a face.
1006 * Designed for use in creation of flat elements separating volume domains.
1007 * A face separating two domains is shared by two volume cells.
1008 * All the nodes are already created (for the two faces).
1009 * Each original Node is associated to corresponding nodes in the domains.
1010 * Some nodes may be duplicated for more than two domains, when domain separations intersect.
1011 * In that case, even some of the nodes to use for the original face may be changed.
1012 * @param vtkVolId: vtk id of a volume containing the face, to get an orientation for the face.
1013 * @param domain1: domain of the original face
1014 * @param domain2: domain of the duplicated face
1015 * @param originalNodes: the vtk node ids of the original face
1016 * @param nodeDomains: map(original id --> map(domain --> duplicated node id))
1017 * @return ok if success.
1019 SMDS_MeshCell* SMDS_UnstructuredGrid::extrudeVolumeFromFace(int vtkVolId,
1022 std::set<int>& originalNodes,
1023 std::map<int, std::map<int, int> >& nodeDomains,
1024 std::map<int, std::map<long, int> >& nodeQuadDomains)
1026 //MESSAGE("extrudeVolumeFromFace " << vtkVolId);
1027 vector<vtkIdType> orderedOriginals;
1028 orderedOriginals.clear();
1029 set<int>::const_iterator it = originalNodes.begin();
1030 for (; it != originalNodes.end(); ++it)
1031 orderedOriginals.push_back(*it);
1034 int nbNodes = this->getOrderedNodesOfFace(vtkVolId, dim, orderedOriginals);
1035 vector<vtkIdType> orderedNodes;
1037 bool isQuadratic = false;
1038 switch (orderedOriginals.size())
1049 isQuadratic = false;
1055 long dom1 = domain1;
1056 long dom2 = domain2;
1057 long dom1_2; // for nodeQuadDomains
1058 if (domain1 < domain2)
1059 dom1_2 = dom1 + INT_MAX * dom2;
1061 dom1_2 = dom2 + INT_MAX * dom1;
1062 //cerr << "dom1=" << dom1 << " dom2=" << dom2 << " dom1_2=" << dom1_2 << endl;
1063 int ima = orderedOriginals.size();
1064 int mid = orderedOriginals.size() / 2;
1065 //cerr << "ima=" << ima << " mid=" << mid << endl;
1066 for (int i = 0; i < mid; i++)
1067 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1068 for (int i = 0; i < mid; i++)
1069 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1070 for (int i = mid; i < ima; i++)
1071 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1072 for (int i = mid; i < ima; i++)
1073 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1074 for (int i = 0; i < mid; i++)
1076 int oldId = orderedOriginals[i];
1078 if (nodeQuadDomains.count(oldId) && nodeQuadDomains[oldId].count(dom1_2))
1079 newId = nodeQuadDomains[oldId][dom1_2];
1082 double *coords = this->GetPoint(oldId);
1083 SMDS_MeshNode *newNode = _mesh->AddNode(coords[0], coords[1], coords[2]);
1084 newId = newNode->getVtkId();
1085 if (! nodeQuadDomains.count(oldId))
1087 std::map<long, int> emptyMap;
1088 nodeQuadDomains[oldId] = emptyMap;
1090 nodeQuadDomains[oldId][dom1_2] = newId;
1092 orderedNodes.push_back(newId);
1097 for (int i = 0; i < nbNodes; i++)
1098 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain1]);
1100 for (int i = 0; i < nbNodes; i++)
1101 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1103 for (int i = nbNodes-1; i >= 0; i--)
1104 orderedNodes.push_back(nodeDomains[orderedOriginals[i]][domain2]);
1110 SMDS_MeshVolume *vol = _mesh->AddVolumeFromVtkIds(orderedNodes);
1115 SMDS_MeshFace *face = _mesh->AddFaceFromVtkIds(orderedNodes);
1119 // TODO update sub-shape list of elements and nodes
1123 //================================================================================
1125 * \brief Allocates data array for ball diameters
1126 * \param MaxVtkID - max ID of a ball element
1128 //================================================================================
1130 void SMDS_UnstructuredGrid::AllocateDiameters( vtkIdType MaxVtkID )
1132 SetBallDiameter( MaxVtkID, 0 );
1135 //================================================================================
1137 * \brief Sets diameter of a ball element
1138 * \param vtkID - vtk id of the ball element
1139 * \param diameter - diameter of the ball element
1141 //================================================================================
1143 void SMDS_UnstructuredGrid::SetBallDiameter( vtkIdType vtkID, double diameter )
1145 vtkDoubleArray* array = vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() );
1148 array = vtkDoubleArray::New();
1149 array->SetNumberOfComponents(1);
1150 vtkDataSet::CellData->SetScalars( array );
1152 array->InsertValue( vtkID, diameter );
1155 //================================================================================
1157 * \brief Returns diameter of a ball element
1158 * \param vtkID - vtk id of the ball element
1160 //================================================================================
1162 double SMDS_UnstructuredGrid::GetBallDiameter( vtkIdType vtkID ) const
1164 if ( vtkDataSet::CellData )
1165 return vtkDoubleArray::SafeDownCast( vtkDataSet::CellData->GetScalars() )->GetValue( vtkID );