]> Creatis software - creaVtk.git/blob - bbtk_creaVtk_PKG/src/bbcreaVtkCreateMeshFromPoints.cxx
463e632ca411f127a9a7f0f9e0cd8473252007c8
[creaVtk.git] / bbtk_creaVtk_PKG / src / bbcreaVtkCreateMeshFromPoints.cxx
1 //===== 
2 // Before editing this file, make sure it's a file of your own (i.e.: it wasn't generated from xml description; if so : your modifications will be lost)
3 //===== 
4 #include "bbcreaVtkCreateMeshFromPoints.h"
5 #include "bbcreaVtkPackage.h"
6
7 #include "vtkTriangleStrip.h"
8 #include <vtkMath.h>
9
10 namespace bbcreaVtk
11 {
12
13 BBTK_ADD_BLACK_BOX_TO_PACKAGE(creaVtk,CreateMeshFromPoints)
14 BBTK_BLACK_BOX_IMPLEMENTATION(CreateMeshFromPoints,bbtk::AtomicBlackBox);
15 //===== 
16 // Before editing this file, make sure it's a file of your own (i.e.: it wasn't generated from xml description; if so : your modifications will be lost)
17 //===== 
18 void CreateMeshFromPoints::Process()
19 {
20
21 // THE MAIN PROCESSING METHOD BODY
22 //   Here we simply set the input 'In' value to the output 'Out'
23 //   And print out the output value
24 // INPUT/OUTPUT ACCESSORS ARE OF THE FORM :
25 //    void bbSet{Input|Output}NAME(const TYPE&)
26 //    const TYPE& bbGet{Input|Output}NAME() const 
27 //    Where :
28 //    * NAME is the name of the input/output
29 //      (the one provided in the attribute 'name' of the tag 'input')
30 //    * TYPE is the C++ type of the input/output
31 //      (the one provided in the attribute 'type' of the tag 'input')
32
33 //    bbSetOutputOut( bbGetInputIn() );
34 //    std::cout << "Output value = " <<bbGetOutputOut() << std::endl;
35
36                 std::vector<double> lstX                = bbGetInputLstX();
37                 std::vector<double> lstY                = bbGetInputLstY();
38                 std::vector<double> lstZ                = bbGetInputLstZ();
39                 std::vector<int> lstIndexs              = bbGetInputLstIndexs();
40                 if ( (lstIndexs.size()<=1) || (lstX.size()==0) || (lstX.size()!=lstY.size()) || (lstY.size()!=lstZ.size()) )
41                 {
42                         printf("Warning! CreateMeshFromPoints::Process: List of points X Y Z  and LstIndexes is not correct\n");
43                         bbSetOutputOut(NULL);
44                 } else  {
45                         int ii,sizeSegment1,sizeSegment2;
46                         int endSegment;
47 //                      vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
48                         if (points!=NULL) points->Delete();
49                         points = vtkPoints::New();
50                         int i,sizeLstX  =       lstX.size();
51                         for (i=0;i<sizeLstX;i++)
52                         {
53                                 points->InsertNextPoint(lstX[i],lstY[i],lstZ[i]);
54                         } // for i
55 //                      vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
56                         if (cells!=NULL) cells->Delete();
57                         cells = vtkCellArray::New();
58                         int maxElements;
59                         int maxSegment1,maxSegment2;
60                         int iSeg1,iSeg2;
61                         int iGeneral    =       0;
62                         int     sizeLstIdexes=lstIndexs.size();
63                         for (i=0; i<sizeLstIdexes-1; i++ )
64                         {
65                                 sizeSegment1 = lstIndexs[i];
66                                 sizeSegment2 = lstIndexs[i+1];
67                                 vtkSmartPointer<vtkTriangleStrip> triangleStrip = vtkSmartPointer<vtkTriangleStrip>::New();
68                                 triangleStrip->GetPointIds()->SetNumberOfIds(sizeSegment1+sizeSegment2);
69                                 maxElements=sizeSegment1;
70                                 if (maxElements<sizeSegment2) maxElements=sizeSegment2;
71                                 maxSegment1     = iGeneral+sizeSegment1;
72                                 maxSegment2     = iGeneral+sizeSegment1+sizeSegment2;
73                                 iSeg1           = iGeneral;
74                                 iSeg2           = iGeneral+sizeSegment1;
75                                 for (ii=0; ii<maxElements; ii++)
76                                 {
77                                         triangleStrip->GetPointIds()->SetId(ii*2  ,iSeg1);
78                                         triangleStrip->GetPointIds()->SetId(ii*2+1,iSeg2);
79                                         iSeg1++;
80                                         iSeg2++;
81                     if (iSeg1>=maxSegment1) { iSeg1=maxSegment1-1; }
82                     if (iSeg2>=maxSegment2) { iSeg2=maxSegment2-1; }
83                                 } // for ii 
84                                 iGeneral=iGeneral+sizeSegment1;
85                                 cells->InsertNextCell(triangleStrip);
86                         } //for  LstIndexs
87
88                         int lastId1 = lstIndexs[0]-1;
89                         int lastId2 = sizeLstX - 1;
90                         int firstId2 = sizeLstX - lstIndexs[sizeLstIdexes - 1];
91                         bool face1open = lstX[0] != lstX[lastId1] && lstY[0] != lstY[lastId1] && lstZ[0] != lstZ[lastId1];
92                         bool face2open = lstX[firstId2] != lstX[lastId2] && lstY[firstId2] != lstY[lastId2] && lstZ[firstId2] != lstZ[lastId2];
93                         
94                         if(bbGetInputCloseSurface())
95                         {                       
96                                 //false = Open Contour
97                                 //true = Closed Contour
98                                 if(!face1open && !face2open)
99                                 {               
100                                         CloseContourSides(lstIndexs, true);
101                                 }else{
102                                         CloseOpenContourSurface(lstIndexs);
103                                 }
104                         }
105                     
106 //                      vtkPolyData *polydata = vtkPolyData::New();
107                         if (polydata!=NULL) polydata->Delete();
108                         polydata = vtkPolyData::New();
109                         polydata->SetPoints(points);
110                         polydata->SetStrips(cells);
111 //                      vtkCleanPolyData *clean=vtkCleanPolyData::New();
112                         if (clean!=NULL) clean->Delete();
113                         clean = vtkCleanPolyData::New();
114                         clean->SetInputData(polydata);
115                         clean->Update();
116 //                      vtkTriangleFilter *triangle = vtkTriangleFilter::New();
117                         if (triangle!=NULL) triangle->Delete();
118                         triangle = vtkTriangleFilter::New();
119                         triangle->SetInputData( clean->GetOutput() );
120                         triangle->Update();
121             bbSetOutputOut( triangle->GetOutput() );
122             //            bbSetOutputOut( clean->GetOutput() );
123                 }// if listXYZ size
124                 //printf("PG CreateMeshFromPoints::Process: End\n");
125 }
126
127 /**
128 * Closes the sides of the contour
129 * iterates in one way or the other, depending on the order of the points and calculated vectors.
130 * uPointOrder: Points are order in a U shape
131 * lstIndexs: number of points on each spline
132 */
133 void CreateMeshFromPoints::CloseContourSides(std::vector<int> lstIndexs, bool uPointOrder){
134         int     sizeLstIdexes = lstIndexs.size();
135         int sizePoints = bbGetInputLstX().size();
136
137         int firstIndex, end, centroidId, numPointsFace, contraryId;
138         int increment = uPointOrder?1:sizeLstIdexes;
139         double centroid[3];
140         
141         for(int facesIdx = 0; facesIdx < 2; facesIdx++){
142                 std::fill(std::begin(centroid), std::end(centroid), 0);
143                 if(facesIdx == 0)
144                 {
145                         firstIndex = 0;
146                         numPointsFace = uPointOrder?lstIndexs[0]: sizeLstIdexes;
147                         end = uPointOrder?firstIndex + numPointsFace:sizePoints - lstIndexs[sizeLstIdexes - 1] + 1;
148                         contraryId = sizePoints-1;
149                 }else{
150                         firstIndex = uPointOrder?sizePoints - lstIndexs[sizeLstIdexes-1]:lstIndexs[0]-1;
151                         numPointsFace = uPointOrder?lstIndexs[sizeLstIdexes-1]:sizeLstIdexes;
152                         end = uPointOrder?firstIndex + numPointsFace:sizePoints;
153                         contraryId = 0;
154                 }
155                 if(numPointsFace > 1)
156                 {
157                         bool validCentroid = CalcValidCentroid(centroid, firstIndex, end, increment, numPointsFace);
158                         if(validCentroid)
159                         {
160                                 bool normalOrder    = isPointingCorrectly(firstIndex, firstIndex+increment, centroid, contraryId);
161                                 centroidId          = points->InsertNextPoint(centroid[0], centroid[1], centroid[2]);
162                                 vtkSmartPointer<vtkTriangleStrip> triangleStrip1 = vtkSmartPointer<vtkTriangleStrip>::New();
163                                 triangleStrip1->GetPointIds()->SetNumberOfIds(numPointsFace*2+2);
164                                 if( normalOrder )
165                 {
166                                         int initial = firstIndex;
167                                         int triangleIndex = 0;
168                                         for(int index = initial; index < end; index+=increment){
169                                                 triangleStrip1->GetPointIds()->SetId(triangleIndex,index);
170                                                 triangleStrip1->GetPointIds()->SetId(triangleIndex+1,centroidId);
171                                                 if(index+increment >= end){
172                                                         triangleStrip1->GetPointIds()->SetId(triangleIndex+2,initial);
173                                                         triangleStrip1->GetPointIds()->SetId(triangleIndex+3,centroidId);
174                                                 }
175                                                 triangleIndex+=2;
176                                         }
177                                         cells->InsertNextCell(triangleStrip1);
178                                 } else {
179                                         int initial = firstIndex-1;
180                                         int triangleIndex = 0;
181                                         int triangleStripStart = end-1;
182                                         for(int index = triangleStripStart; index > initial; index-=increment){
183                                                 triangleStrip1->GetPointIds()->SetId(triangleIndex,index);
184                                                 triangleStrip1->GetPointIds()->SetId(triangleIndex+1,centroidId);
185                                                 if(index-increment <= initial){
186                                                         triangleStrip1->GetPointIds()->SetId(triangleIndex+2,triangleStripStart);
187                                                         triangleStrip1->GetPointIds()->SetId(triangleIndex+3,centroidId);
188                                                 }
189                                                 triangleIndex+=2;
190                                         }
191                                         cells->InsertNextCell(triangleStrip1);
192                                 }
193                         }//if validCentroid
194                 }//if numPointsFace
195         }//for facesIdx
196
197 }
198
199 /**
200 * Checks if the normal from firstPointId, secPointId and centroid points away 
201 * from the vector centroid to contrPointId.
202 * Used to check that the order used to create the new polygons is correct.
203 */
204 bool CreateMeshFromPoints::isPointingCorrectly( int firstPointId, int secPointId, double(&centroid)[3], int contrPointId) {
205
206         double firstPoint[3], secPoint[3], contrPoint[3];
207         points->GetPoint(firstPointId, firstPoint);
208         points->GetPoint(secPointId, secPoint);
209         
210         double firstVect[3], secVect[3], normal[3], contrVect[3];
211         
212         vtkMath::Subtract(firstPoint, centroid, firstVect);
213         vtkMath::Subtract(secPoint, centroid, secVect);
214         
215         points->GetPoint(contrPointId, contrPoint);
216         vtkMath::Subtract(contrPoint, centroid, contrVect);
217         
218         vtkMath::Cross(firstVect, secVect, normal);
219         double dotCalc;
220         dotCalc = vtkMath::Dot(normal, contrVect);
221         
222         return dotCalc>0;
223 }
224
225 /**
226 * Checks if points on each side of the shapes represent a curve.
227 */
228 bool CreateMeshFromPoints::CheckLinePointOrder(){
229         int sizePoints = bbGetInputLstX().size();
230         std::vector<int> lstIndexs = bbGetInputLstIndexs();
231         double point1[3], point2[3], point3[3];
232         double center[3];
233         double firstRadiusSum = 0;
234         double secondRadiusSum = 0;
235         for(int i = 0; i < lstIndexs[0] && lstIndexs[0] > 3; i+=3){
236                 if(i+3 <= lstIndexs[0]){
237                         points->GetPoint(i, point1);
238                         points->GetPoint(i+1, point2);
239                         points->GetPoint(i+2, point3);
240                         firstRadiusSum += vtkMath::Solve3PointCircle(point1, point2, point3, center);
241                 }
242         }
243         
244         for(int i = 0; i < sizePoints && lstIndexs.size() > 3; i+=(lstIndexs[0]*3)){
245                 if(i+(3*lstIndexs[0]) <= sizePoints){
246                         points->GetPoint(i, point1);
247                         points->GetPoint(i+lstIndexs[0], point2);
248                         points->GetPoint(i+(2*lstIndexs[0]), point3);
249                         secondRadiusSum += vtkMath::Solve3PointCircle(point1, point2, point3, center);
250                 }
251         }
252         
253         return firstRadiusSum > secondRadiusSum;
254 }
255
256 /**
257 * Closes an open contour
258 * lstIndexs: number of points on each spline
259 */
260 void CreateMeshFromPoints::CloseOpenContourSurface(std::vector<int> lstIndexs){
261         bool linePointOrder = CheckLinePointOrder();
262         CloseContourSides(lstIndexs, !linePointOrder);
263         CloseContourBottom(!linePointOrder);
264 }
265
266 /**
267 * Calculates centroid and checks if points are collinear.
268 * centroid: array to store calculation
269 * start: start index of points to use
270 * end: end index of points to use
271 * increment: increment to be used in point iteration
272 * numPoints: number of points used to calculate the centroid.
273 * Returns a bool indicating the validity of the centroid calculated.
274 * False = invalid centroid = all points are the same.
275 */
276 bool CreateMeshFromPoints::CalcValidCentroid(double(&centroid)[3], int start, int end, int increment, int numPoints){
277         double currPoint[3] = {}, prevPoint[3] = {}, middlePoint[3] = {}, firstPoint[3] = {};
278         double vector1[3], vector2[3];
279         bool samePoint = true;
280         int splineMidPoint = numPoints/2;
281         bool collinear = true;
282         
283         points->GetPoint(start, firstPoint);
284         points->GetPoint(splineMidPoint, middlePoint);
285         vtkMath::Subtract(middlePoint, firstPoint, vector1);
286         
287         for(int i = start; i < end; i+=increment){
288                 points->GetPoint(i, currPoint);
289                 if(samePoint && i > start && (currPoint[0] != prevPoint[0] || currPoint[1] != prevPoint[1] || currPoint[2] != prevPoint[2])){
290                         samePoint = false;
291                 }
292
293                 vtkMath::Subtract(currPoint, firstPoint, vector2);
294                 double angle = vtkMath::AngleBetweenVectors(vector1, vector2);
295                 if(angle > 0.0001 && collinear){
296                         collinear = false;
297                 }
298                 
299                 centroid[0] += currPoint[0];
300                 centroid[1] += currPoint[1];
301                 centroid[2] += currPoint[2];
302                 std::copy(std::begin(currPoint), std::end(currPoint), prevPoint);
303         }
304         
305         centroid[0] /= numPoints;
306         centroid[1] /= numPoints;
307         centroid[2] /= numPoints;
308         
309         return !samePoint && !collinear;
310 }
311
312 /**
313 * Closes the bottom of the given countour.
314 * Should only be used when its an open contour.
315 * uPointOrder: points are ordered in U shape
316 */
317 void CreateMeshFromPoints::CloseContourBottom(bool uPointOrder){
318         std::vector<int> lstIndexs = bbGetInputLstIndexs();
319         int sizeLstIdexes = lstIndexs.size();
320         int sizeLstX = bbGetInputLstX().size();
321         
322         vtkSmartPointer<vtkTriangleStrip> triangleStripBottom = vtkSmartPointer<vtkTriangleStrip>::New();
323         triangleStripBottom->GetPointIds()->SetNumberOfIds(sizeLstIdexes*2);
324         int triangleIndex = 0, currentId = 0, nextId = 0;
325         for(int splineIndex = 0; splineIndex < sizeLstIdexes;splineIndex++){
326                 triangleStripBottom->GetPointIds()->SetId(triangleIndex, currentId);
327                 nextId = uPointOrder?currentId + lstIndexs[splineIndex] - 1:sizeLstX - sizeLstIdexes + splineIndex;
328                 triangleStripBottom->GetPointIds()->SetId(triangleIndex+1, nextId);
329                 triangleIndex+=2;
330                 currentId = uPointOrder?nextId + 1: splineIndex+1;
331         }
332         cells->InsertNextCell(triangleStripBottom);
333 }
334
335 //===== 
336 // Before editing this file, make sure it's a file of your own (i.e.: it wasn't generated from xml description; if so : your modifications will be lost)
337 //===== 
338 void CreateMeshFromPoints::bbUserSetDefaultValues()
339 {
340
341 //  SET HERE THE DEFAULT INPUT/OUTPUT VALUES OF YOUR BOX 
342 //    Here we initialize the input 'In' to 0
343 //   bbSetInputIn(0);
344         bbSetInputCloseSurface(false);
345         points          = NULL;
346         cells           = NULL;
347         polydata        = NULL;
348         clean           = NULL;
349         triangle        = NULL;
350 }
351 //===== 
352 // Before editing this file, make sure it's a file of your own (i.e.: it wasn't generated from xml description; if so : your modifications will be lost)
353 //===== 
354 void CreateMeshFromPoints::bbUserInitializeProcessing()
355 {
356
357 //  THE INITIALIZATION METHOD BODY :
358 //    Here does nothing 
359 //    but this is where you should allocate the internal/output pointers 
360 //    if any 
361
362   
363 }
364 //===== 
365 // Before editing this file, make sure it's a file of your own (i.e.: it wasn't generated from xml description; if so : your modifications will be lost)
366 //===== 
367 void CreateMeshFromPoints::bbUserFinalizeProcessing()
368 {
369
370 //  THE FINALIZATION METHOD BODY :
371 //    Here does nothing 
372 //    but this is where you should desallocate the internal/output pointers 
373 //    if any
374   
375 }
376 }
377 // EO namespace bbcreaVtk
378
379