On 18/09/12 11:43, Luca Antiga wrote:
> Dear Roman,
>   thanks for sharing your code so we can look at it.
>
> The approach you use is to loop over all points in the dataset and check 
> which ones are used by just one cell.
> This is true for boundary points, but it's also true for points along a poly 
> line. In fact one possible strategy, as you
> already envisioned, is to convert polylines to lines, so that all 
> non-boundary points are used by two or more cells.

I figured out some code  that converts vtkPolyLines to vtkLines, see below. It 
seems to work fine in this special context.

> However, I think the cleanest way to go is to loop over cells (not points) 
> and for every cell consider the two endpoints
> and check whether one of the endpoints is used by only one cell (the same 
> cell you're on). This will get you the
> boundary you're after.

I agree, I just found the conversion from vtkPolyLines to vtkLines more 
versatile. Hope it results in the same, otherwise I'll try to implement code to 
extract the endpoints of the graph made up of vtkPolyLines

Best wishes
Roman


     ////vtkPolyLine2vtkLine(mesh);
     vtkCell *cell;
     vtkCell *edge;
     vtkPolyData *new_mesh= vtkPolyData::New();
     vtkPoints *points= vtkPoints::New();
     vtkCellArray *lines= vtkCellArray::New();

     points= mesh->GetPoints(); //GetVerts()

     for (vtkIdType i= 0; i < N; i++){
     cell= mesh->GetCell(i);
     if (cell->GetCellType() == VTK_POLY_LINE) {
         std::cerr << "Cell " << i << " is of type " << 
vtkCellTypes::GetClassNameFromTypeId(mesh->GetCell(i)->GetCellType()) << 
std::endl;
         vtkIdType noe= cell->GetNumberOfEdges(); //vtkPolyLine does not 
consist 
of edges!!!
         //vtkIdType noe= cell->GetNumberOfLines(); vtkCell has no member named 
GetNumberOfLines
         vtkIdType nop= cell->GetNumberOfPoints();
         std::cerr << "Cell contains " << noe << " edges and " << nop << " 
points." << std::endl;

         for (vtkIdType j= 1; j < nop; j++){ //start at 1 since a polyline has 
nop-1 lines!
         vtkIdType pi= cell->GetPointId(j);
         vtkLine *line= vtkLine::New();
         line->GetPointIds()->SetId(0,cell->GetPointId(j - 1));
         line->GetPointIds()->SetId(1,cell->GetPointId(j));
         //line->GetPointIds()->InsertNextId(edge->GetPointIds()->GetId(0));
         //line->GetPointIds()->InsertNextId(edge->GetPointIds()->GetId(1));
         lines->InsertNextCell(line);
         }
         }
     }

     new_mesh->SetPoints(points);
     new_mesh->SetLines(lines);


> Let me know if you need further hints.
>
> Best,
>
>
> Luca
>
>
> On Sep 13, 2012, at 9:58 PM, Dr. Roman Grothausmann wrote:
>
>> Dear Luca,
>>
>>
>> Thanks for Your reply. Below is my code with which I try to extract the end 
>> points, but somehow it finds that all points (from attached network from 
>> vmtknetworkextraction) are end points of vtkPolyLines. As You advice, I call 
>> BuildLinks and GetPointCells. What am I missing here? Do I need to convert 
>> vtkPolyLine to vtkLines? If so, how would I have to do that?
>>
>> Many thanks for Your help
>> Roman
>>
>>
>> #include<vtkXMLPolyDataReader.h>//for vtp-files (cannot contain 3D cells?)
>> //#include<vtkPolyDataReader.h>//for vtk-files (cannot contain 3D cells?)
>> #include<vtkPolyData.h>
>> #include<vtkCell.h>
>> #include<vtkCellType.h>
>> #include<vtkIdList.h>
>>
>> int main(int argc, char* argv[]){
>>
>>     if( argc<= 1 )
>>         {
>>         std::cerr<<  "Usage: "<<  argv[0];
>>         std::cerr<<  " inputGraph";
>>         //std::cerr<<  " outputPoints";
>>         std::cerr<<  std::endl;
>>         return EXIT_FAILURE;
>>         }
>>
>>     if(!(strcasestr(argv[1],".vtp"))) {
>>      std::cout<<  "The input should end with .vtp"<<  std::endl;
>>         return -1;
>>      }
>>
>>     //vtkPolyDataReader *reader = vtkPolyDataReader::New(); //*.vtk
>>     vtkXMLPolyDataReader *reader = vtkXMLPolyDataReader::New(); //*.vtp
>>
>>     reader->SetFileName(argv[1]);
>>     reader->Update();
>>
>>     vtkPolyData *mesh= reader->GetOutput();
>>     mesh->BuildLinks(); //for inverse topological queries
>>
>>     vtkIdType N= mesh->GetNumberOfCells();
>>     std::cerr<<  "Mesh contains "<<  N<<  " cells and "<<  
>> mesh->GetNumberOfPoints()<<  " points."<<  std::endl;
>>     double x[3];
>>
>>
>>     vtkIdType n= mesh->GetNumberOfPoints();
>>     //vtkIdType n= mesh->GetNumberOfVerts();
>>
>>     for (vtkIdType j= 0; j<  n; j++){
>>      vtkIdList *cellIdList= vtkIdList::New();
>>      mesh->GetPointCells(j, cellIdList); //obtain cells using a particular 
>> point. Make sure that routine BuildLinks() has been called.
>>      vtkIdType m= cellIdList->GetNumberOfIds();
>>      //for (vtkIdType k= 0; k<  m; k++){
>>      if (m<  1)
>>          std::cerr<<  "Point "<<  j<<  " is not part of any cell!"<<  
>> std::endl;
>>      else if (m == 1){
>>          //mesh->GetPoint(cell->GetPointId(j), x);
>>          mesh->GetPoint(j, x);
>>          //std::cout<<  x[0]<<  x[1]<<  x[2]<<  std::endl;
>>          printf("%d\t%f\t%f\t%f\n", j, x[0], x[1], x[2]);
>>          std::cerr<<  "Point "<<  j<<  " is just part of cell "<<  
>> cellIdList->GetId(0)<<  " which is of type "<<  
>> vtkCellTypes::GetClassNameFromTypeId(mesh->GetCell(cellIdList->GetId(0))->GetCellType())<<
>>   std::endl;    //  VTK_LINE || VTK_POLY_LINE) //from 
>> VTK/Filtering/vtkCellType.h, VTK_POLY_LINE: set of 1D lines
>>          }
>>      else
>>          std::cerr<<  "Point "<<  j<<  " is part of more than one cell!"<<  
>> std::endl;
>>      }
>>
>>     return EXIT_SUCCESS;
>>     }
>>
>>
>>
>> On 05/09/12 15:33, Luca Antiga wrote:
>>> Hi Roman,
>>>   in fact it does it. If you are acquainted with VTK then it is easy: just
>>> take the network output, cycle through all the cells and for each cell
>>> extract the coordinate of the endpoint if that point is not used by any
>>> other cell (you have to call GetPointCells, defined in the vtkPolyData
>>> class; don't forget to call BuildLinks beforehand so that you can
>>> make this kind of inverse topological queries).
>>> Best,
>>>
>>> Luca
>>>
>>>
>>> On Sep 5, 2012, at 10:12 AM, Dr. Roman Grothausmann wrote:
>>>
>>>> Dear Luca,
>>>>
>>>>
>>>> Thanks for Your reply. I'll give it  a try. From my understanding
>>>> vtkFeatureEdges does already this kind of thing, just for edges not for
>>>> points.
>>>>
>>>> Thanks for Your help
>>>> Roman
>>>>
>>>> On 03/09/12 11:39, Luca Antiga wrote:
>>>>> Dear Roman,
>>>>>   thanks for your appreciation of vmtk, it's good to have you onboard.
>>>>>
>>>>> It is indeed possible to use manually specified endpoints for 
>>>>> vmtkcenterlines,
>>>>> but there's some glue code that needs to be written in-between.
>>>>>
>>>>> vmtknetworkextraction outputs a polydata dataset composed of 
>>>>> interconnected
>>>>> lines, whose endpoints are the coordinates you're interested in. In order 
>>>>> to
>>>>> extract endpoints from the Network you should look for cells that have 
>>>>> either
>>>>> one of the Topology array values equal to -1 and extract the corresponding
>>>>> coordinates.
>>>>>
>>>>> Once you have the coordinates, you can specify them even at the command
>>>>> line using vmtkcenterlines using -seedselector pointlist -seeds x y z 
>>>>> -targets x0 y0 z0 x1 y1 z1 ...
>>>>>
>>>>> Let me know if you need support/advice for the extraction of the 
>>>>> endpoints. In
>>>>> fact this could be a nice addition to vmtk.
>>>>>
>>>>> Best,
>>>>>
>>>>>
>>>>> Luca
>>>>>
>>>>>
>>>>> On Sep 3, 2012, at 11:26 AM, Dr. Roman Grothausmann wrote:
>>>>>
>>>>>> Dear mailing list members,
>>>>>>
>>>>>>
>>>>>> vmtk is really cool. I'm just about getting to know its features
>>>>>> especially those for centre line extraction.
>>>>>> I found this email from Luca:
>>>>>> http://sourceforge.net/mailarchive/message.php?msg_id=29119553 to use
>>>>>> vmtknetworkextraction for this purpose.
>>>>>>
>>>>>> Now I wonder, if it is possible to use the endpoints of the result of
>>>>>> vmtknetworkextraction as the endpoints for vmtkcenterlines?
>>>>>> Or is there a way to save the manually selected start and end points of
>>>>>> vmtkcenterlines to use these when experimenting with other options of
>>>>>> vmtkcenterlines?
>>>>>>
>>>>>> Many thanks for vmtk and for any help or hints.
>>>>>> Roman
>>>>>>
>>>>>> --
>>>>>> Dr. Roman Grothausmann
>>>>>>
>>>>>> Tomographie und Digitale Bildverarbeitung
>>>>>> Tomography and Digital Image Analysis
>>>>>>
>>>>>> Institut für Funktionelle und Angewandte Anatomie, OE 4120
>>>>>> Medizinische Hochschule Hannover
>>>>>> Carl-Neuberg-Str. 1
>>>>>> 30625 Hannover
>>>>>>
>>>>>> Tel. +49 511 532-9574
>>>>>>
>>>>>>
>>>>>> ------------------------------------------------------------------------------
>>>>>> Live Security Virtual Conference
>>>>>> Exclusive live event will cover all the ways today's security and
>>>>>> threat landscape has changed and how IT managers can respond. Discussions
>>>>>> will include endpoint security, mobile security and the latest in malware
>>>>>> threats. http://www.accelacomm.com/jaw/sfrnl04242012/114/50122263/
>>>>>> _______________________________________________
>>>>>> vmtk-users mailing list
>>>>>> vmtk-users@lists.sourceforge.net
>>>>>> https://lists.sourceforge.net/lists/listinfo/vmtk-users
>>>>>
>>>>
>>>> --
>>>> Dr. Roman Grothausmann
>>>>
>>>> Tomographie und Digitale Bildverarbeitung
>>>> Tomography and Digital Image Analysis
>>>>
>>>> Institut für Funktionelle und Angewandte Anatomie, OE 4120
>>>> Medizinische Hochschule Hannover
>>>> Carl-Neuberg-Str. 1
>>>> 30625 Hannover
>>>>
>>>> Tel. +49 511 532-9574
>>>>
>>>> ------------------------------------------------------------------------------
>>>> Live Security Virtual Conference
>>>> Exclusive live event will cover all the ways today's security and
>>>> threat landscape has changed and how IT managers can respond. Discussions
>>>> will include endpoint security, mobile security and the latest in malware
>>>> threats. http://www.accelacomm.com/jaw/sfrnl04242012/114/50122263/
>>>> _______________________________________________
>>>> vmtk-users mailing list
>>>> vmtk-users@lists.sourceforge.net
>>>> https://lists.sourceforge.net/lists/listinfo/vmtk-users
>>>
>>
>> --
>> Dr. Roman Grothausmann
>>
>> Tomographie und Digitale Bildverarbeitung
>> Tomography and Digital Image Analysis
>>
>> Institut für Funktionelle und Angewandte Anatomie, OE 4120
>> Medizinische Hochschule Hannover
>> Carl-Neuberg-Str. 1
>> 30625 Hannover
>>
>> Tel. +49 511 532-9574
>> <Bronc01_bi4_p1_c-open_nw.vtp>
>

-- 
Dr. Roman Grothausmann

Tomographie und Digitale Bildverarbeitung
Tomography and Digital Image Analysis

Institut für Funktionelle und Angewandte Anatomie, OE 4120
Medizinische Hochschule Hannover
Carl-Neuberg-Str. 1
30625 Hannover

Tel. +49 511 532-9574

------------------------------------------------------------------------------
Don't let slow site performance ruin your business. Deploy New Relic APM
Deploy New Relic app performance management and know exactly
what is happening inside your Ruby, Python, PHP, Java, and .NET app
Try New Relic at no cost today and get our sweet Data Nerd shirt too!
http://p.sf.net/sfu/newrelic-dev2dev
_______________________________________________
vmtk-users mailing list
vmtk-users@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/vmtk-users

Reply via email to