Salome HOME
PR: debug compact vtkUnstructuredGrid in progress
[modules/smesh.git] / src / SMDS / SMDS_UnstructuredGrid.cxx
1
2
3 #include "SMDS_UnstructuredGrid.hxx"
4 #include "utilities.h"
5
6 #include <vtkCellArray.h>
7 #include <vtkIdTypeArray.h>
8 #include <vtkUnsignedCharArray.h>
9
10 #include <list>
11
12 using namespace std;
13
14 SMDS_UnstructuredGrid* SMDS_UnstructuredGrid::New()
15 {
16         MESSAGE("SMDS_UnstructuredGrid::New");
17         return new SMDS_UnstructuredGrid();
18 }
19
20 SMDS_UnstructuredGrid::SMDS_UnstructuredGrid() : vtkUnstructuredGrid()
21 {
22 }
23
24 SMDS_UnstructuredGrid::~SMDS_UnstructuredGrid()
25 {
26 }
27
28
29 unsigned long SMDS_UnstructuredGrid::GetMTime()
30 {
31         unsigned long mtime = vtkUnstructuredGrid::GetMTime();
32         MESSAGE("vtkUnstructuredGrid::GetMTime: " << mtime);
33         return mtime;
34 }
35
36 void SMDS_UnstructuredGrid::UpdateInformation()
37 {
38         MESSAGE("SMDS_UnstructuredGrid::UpdateInformation");
39         return vtkUnstructuredGrid::UpdateInformation();
40 }
41
42 void SMDS_UnstructuredGrid::compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize,
43                                                                                 std::vector<int>& idCellsOldToNew, int newCellSize)
44 {
45         MESSAGE("------------------------- SMDS_UnstructuredGrid::compactGrid " << newNodeSize << " " << newCellSize);
46
47         int startHole = 0;
48         int endHole = 0;
49         int startBloc = 0;
50         int endBloc = 0;
51         int alreadyCopied = 0;
52         typedef enum {lookHoleStart, lookHoleEnd, lookBlocEnd} enumState;
53         enumState compactState = lookHoleStart;
54
55         // --- if newNodeSize, create a new compacted vtkPoints
56
57         vtkPoints *newPoints = 0;
58         if (newNodeSize)
59         {
60                 MESSAGE("-------------- compactGrid, newNodeSize");
61                 newPoints = vtkPoints::New();
62                 newPoints->Initialize();
63                 newPoints->Allocate(newNodeSize);
64                 newPoints->SetNumberOfPoints(newNodeSize);
65                 int oldNodeSize = idNodesOldToNew.size();
66
67                 for (int i=0; i< oldNodeSize; i++)
68                 {
69                         switch(compactState)
70                         {
71                         case lookHoleStart:
72                                 if (idNodesOldToNew[i] < 0)
73                                 {
74                                         startHole = i;
75                                         compactState = lookHoleEnd;
76                                 }
77                                 break;
78                         case lookHoleEnd:
79                                 if (idNodesOldToNew[i] >= 0)
80                                 {
81                                         endHole = i;
82                                         startBloc = i;
83                                         compactState = lookBlocEnd;
84                                 }
85                                 break;
86                         case lookBlocEnd:
87                                 if (idNodesOldToNew[i] < 0) endBloc = i; // see nbPoints below
88                                 else if (i == (oldNodeSize-1)) endBloc = i+1;
89                                 if (endBloc)
90                                 {
91                                         MESSAGE("-------------- newNodeSize, endbloc");
92                                         void *target = newPoints->GetVoidPointer(3*alreadyCopied);
93                                         void *source = this->Points->GetVoidPointer(3*startBloc);
94                                         int nbPoints = endBloc - startBloc;
95                                         memcpy(target, source, 3*sizeof(float)*nbPoints);
96                                         for (int j=startBloc; j<endBloc; j++)
97                                                 idNodesOldToNew[j] = alreadyCopied++;
98                                         compactState = lookHoleStart;
99                                         startHole = i;
100                                         endHole = 0;
101                                         startBloc = 0;
102                                         endBloc = 0;
103                                 }
104                                 break;
105                         }
106                 }
107                 if (!alreadyCopied) // no hole, but shorter, no need to modify idNodesOldToNew
108                 {
109                         MESSAGE("------------- newNodeSize, shorter")
110                         void *target = newPoints->GetVoidPointer(0);
111                         void *source = this->Points->GetVoidPointer(0);
112                         int nbPoints = newNodeSize;
113                         memcpy(target, source, 3*sizeof(float)*nbPoints);
114                 }
115         }
116
117         // --- create new compacted Connectivity, Locations and Types
118
119     int oldCellSize = this->Types->GetNumberOfTuples();
120
121         vtkCellArray *newConnectivity = vtkCellArray::New();
122         newConnectivity->Initialize();
123         int oldCellDataSize = this->Connectivity->GetData()->GetSize();
124         newConnectivity->Allocate(oldCellDataSize);
125         MESSAGE("oldCellSize="<< oldCellSize << " oldCellDataSize=" << oldCellDataSize);
126
127         vtkUnsignedCharArray *newTypes = vtkUnsignedCharArray::New();
128         newTypes->Initialize();
129         //newTypes->Allocate(oldCellSize);
130         newTypes->SetNumberOfValues(newCellSize);
131
132         vtkIdTypeArray *newLocations = vtkIdTypeArray::New();
133         newLocations->Initialize();
134         //newLocations->Allocate(oldCellSize);
135         newLocations->SetNumberOfValues(newCellSize);
136
137         startHole = 0;
138         endHole = 0;
139         startBloc = 0;
140         endBloc = 0;
141         alreadyCopied = 0;
142         compactState = lookHoleStart;
143
144         vtkIdType tmpid[50];
145         vtkIdType *pointsCell =&tmpid[0]; // --- points id to fill a new cell
146
147     for (int i=0; i<oldCellSize; i++)
148     {
149                 switch(compactState)
150                 {
151                 case lookHoleStart:
152                         if (this->Types->GetValue(i) == VTK_EMPTY_CELL)
153                         {
154                                 startHole = i;
155                                 compactState = lookHoleEnd;
156                         }
157                         break;
158                 case lookHoleEnd:
159                         if (this->Types->GetValue(i) != VTK_EMPTY_CELL)
160                         {
161                                 endHole = i;
162                                 startBloc = i;
163                                 compactState = lookBlocEnd;
164                         }
165                         break;
166                 case lookBlocEnd:
167                         if (this->Types->GetValue(i) == VTK_EMPTY_CELL) endBloc =i;
168                         else if (i == (oldCellSize-1)) endBloc = i+1;
169                         if (endBloc)
170                         {
171                         MESSAGE(" -------- newCellSize, endBloc");
172                                 for (int j=startBloc; j<endBloc; j++)
173                                 {
174                                         newTypes->SetValue(alreadyCopied, this->Types->GetValue(j));
175                                         idCellsOldToNew[j] = alreadyCopied;
176                                         int oldLoc = this->Locations->GetValue(j);
177                                         int nbpts;
178                                         this->Connectivity->GetCell(oldLoc, nbpts, pointsCell);
179                                         for (int l=0; l<nbpts; l++)
180                                                 pointsCell[l] = idNodesOldToNew[pointsCell[l]];
181                                         int newcnt = newConnectivity->InsertNextCell(nbpts, pointsCell);
182                                         int newLoc = newConnectivity->GetInsertLocation(nbpts);
183                                         newLocations->SetValue(alreadyCopied, newLoc);
184                                         alreadyCopied++;
185                                 }
186                         }
187                         break;
188                 }
189     }
190     if (!alreadyCopied) // no hole, but shorter
191     {
192         MESSAGE(" -------- newCellSize, shorter");
193                 for (int j=0; j<oldCellSize; j++)
194                 {
195                         newTypes->SetValue(alreadyCopied, this->Types->GetValue(j));
196                         idCellsOldToNew[j] = alreadyCopied;
197                         int oldLoc = this->Locations->GetValue(j);
198                         int nbpts;
199                         this->Connectivity->GetCell(oldLoc, nbpts, pointsCell);
200                         //MESSAGE(j << " " << alreadyCopied << " " << (int)this->Types->GetValue(j) << " " << oldLoc << " " << nbpts );
201                         for (int l=0; l<nbpts; l++)
202                         {
203                                 int oldval = pointsCell[l];
204                                 pointsCell[l] = idNodesOldToNew[oldval];
205                                 //MESSAGE("   " << oldval << " " << pointsCell[l]);
206                         }
207                         int newcnt = newConnectivity->InsertNextCell(nbpts, pointsCell);
208                         int newLoc = newConnectivity->GetInsertLocation(nbpts);
209                         //MESSAGE(newcnt << " " << newLoc);
210                         newLocations->SetValue(alreadyCopied, newLoc);
211                         alreadyCopied++;
212                 }
213     }
214
215     newConnectivity->Squeeze();
216     //newTypes->Squeeze();
217     //newLocations->Squeeze();
218
219     if (newNodeSize)
220     {
221         MESSAGE("------- newNodeSize, setPoints");
222         this->SetPoints(newPoints);
223     }
224     this->SetCells(newTypes, newLocations, newConnectivity);
225     this->BuildLinks();
226 }