]> Creatis software - cpMesh.git/blob - appli/examples/example_QuadEdgeDecimation.cxx
Read/Write image plugin integration
[cpMesh.git] / appli / examples / example_QuadEdgeDecimation.cxx
1 #include <cmath>
2 #include <cstdlib>
3 #include <iostream>
4 #include <string>
5
6 #include <cpm/DataStructures/QuadEdgeMesh.h>
7 #include <cpm/IO/MeshReader.h>
8 #include <itkRegularSphereMeshSource.h>
9 #include <vnl/vnl_math.h>
10
11 #include <vtkActor.h>
12 #include <vtkCallbackCommand.h>
13 #include <vtkProperty.h>
14 #include <vtkRenderer.h>
15 #include <vtkRenderWindow.h>
16 #include <vtkRenderWindowInteractor.h>
17 #include <vtkSmartPointer.h>
18 #include <cpm/VTK/MeshMapper.h>
19
20 #include <cpm/Algorithms/QuadEdge/DecimationCriteria.h>
21 #include <cpm/Algorithms/QuadEdge/SquaredEdgeLengthDecimationFilter.h>
22
23
24 // -------------------------------------------------------------------------
25 const unsigned int Dimension = 3;
26 typedef float TScalar;
27 typedef cpm::DataStructures::QuadEdgeMesh< TScalar, Dimension > TMesh;
28 typedef cpm::VTK::MeshMapper< TMesh > TMeshMapper;
29
30 typedef cpm::Algorithms::QuadEdge::NumberOfPointsCriterion< TMesh > TCriterion;
31 /*
32   typedef cpm::Algorithms::QuadEdge::NumberOfFacesCriterion< TMesh > TCriterion;
33   typedef cpm::Algorithms::QuadEdge::MaxMeasureBoundCriterion< TMesh > TCriterion;
34   typedef cpm::Algorithms::QuadEdge::MinMeasureBoundCriterion< TMesh > TCriterion;
35 */
36 typedef cpm::Algorithms::QuadEdge::SquaredEdgeLengthDecimationFilter< TMesh, TMesh, TCriterion > TDecimator;
37
38 // -------------------------------------------------------------------------
39 int main( int argc, char* argv[] )
40 {
41   if( argc < 2 )
42   {
43     std::cerr
44       << "Usage: " << argv[ 0 ]
45       << " input_mesh"
46       << std::endl;
47     return( 1 );
48
49   } // fi
50
51   typedef cpm::IO::MeshReader< TMesh > TReader;
52   TReader::Pointer reader = TReader::New( );
53   reader->SetFileName( argv[ 1 ] );
54   reader->Update( );
55   TMesh::Pointer mesh = reader->GetOutput( );
56
57   TCriterion::Pointer criterion = TCriterion::New( );
58   criterion->TopologicalChangeOn( );
59   criterion->SetNumberOfElements( mesh->GetNumberOfPoints( ) / 10 );
60   criterion->SetMeasureBound( 1 );
61
62   TDecimator::Pointer decimator = TDecimator::New( );
63   decimator->SetInput( mesh );
64   decimator->SetCriterion( criterion );
65   decimator->Update( );
66
67   // Map mesh
68   vtkSmartPointer< TMeshMapper > mapper =
69     vtkSmartPointer< TMeshMapper >::New( );
70   mapper->SetInputData( decimator->GetOutput( ) );
71
72   // Create actor
73   vtkSmartPointer< vtkActor > actor =
74     vtkSmartPointer< vtkActor >::New( );
75   actor->SetMapper( mapper );
76   actor->GetProperty( )->SetColor( 1, 1, 0 );
77   actor->GetProperty( )->SetOpacity( 1 );
78
79   // Configure visualization objects
80   vtkSmartPointer< vtkRenderer > renderer =
81     vtkSmartPointer< vtkRenderer >::New( );
82   renderer->SetBackground( 0.1, 0.3, 0.5 );
83
84   vtkSmartPointer< vtkRenderWindow > window =
85     vtkSmartPointer< vtkRenderWindow >::New( );
86   window->AddRenderer( renderer );
87   window->SetSize( 800, 800 );
88
89   // Set up the interaction
90   vtkSmartPointer< vtkRenderWindowInteractor > interactor =
91     vtkSmartPointer< vtkRenderWindowInteractor >::New( );
92   window->SetInteractor( interactor );
93
94   // Associate actors
95   renderer->AddActor( actor );
96
97   // Begin interaction
98   window->Render( );
99   interactor->Start( );
100
101   return( 0 );
102 }
103
104 // eof - $RCSfile$