]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STAT/Macros/kDTreeTest.C
1. Removed obsolete test functions
[u/mrichter/AliRoot.git] / STAT / Macros / kDTreeTest.C
1 /*
2   Test macro for TKDTree
3   
4   //Initialization:
5
6   gSystem->AddIncludePath("-I$ALICE_ROOT/STAT")
7   gSystem->Load("$ALICE_ROOT/lib/tgt_linux/libSTAT.so");
8   .L $ALICE_ROOT/STAT/Macros/kDTreeTest.C+ 
9
10   
11   TestBuild();       // test build function of kdTree for memory leaks
12   TestSpeed();       // test the CPU consumption to build kdTree
13   TestkdtreeIF();    // test functionality of the kdTree
14   TestSizeIF();      // test the size of kdtree - search application - Alice TPC tracker situation
15   //
16 */
17
18 #include <malloc.h>
19 #include "TSystem.h"
20 #include "TMatrixD.h"
21 #include "TRandom.h"
22 #include "TGraph.h"
23 #include "TStopwatch.h"
24 #include "TKDTree.h"
25
26
27
28
29 void TestBuild(const Int_t npoints = 1000000, const Int_t bsize = 100);
30 void TestSpeed(const Int_t npower2 = 20, const Int_t bsize = 10);
31
32 void TestkdtreeIF(Int_t npoints=1000, Int_t bsize=9, Int_t nloop=1000, Int_t mode = 2);
33 void TestSizeIF(Int_t nsec=36, Int_t nrows=159, Int_t npoints=1000,  Int_t bsize=10, Int_t mode=1);
34
35
36
37 Float_t Mem()
38 {
39   // get mem info
40   ProcInfo_t procInfo;
41   gSystem->GetProcInfo(&procInfo);
42   return procInfo.fMemVirtual;
43 }
44
45 //______________________________________________________________________
46 void kDTreeTest()
47 {
48   //
49   //
50   //
51   printf("\n\tTesting kDTree memory usage ...\n");
52   TestBuild();  
53   printf("\n\tTesting kDTree speed ...\n");
54   TestSpeed();
55 }
56
57 //______________________________________________________________________
58 void TestBuild(const Int_t npoints, const Int_t bsize){  
59   //
60   // Test kdTree for memory leaks
61   //
62   Float_t *data0 =  new Float_t[npoints*2];
63   Float_t *data[2];
64   data[0] = &data0[0];
65   data[1] = &data0[npoints];
66   for (Int_t i=0;i<npoints;i++) {
67     data[1][i]= gRandom->Rndm();
68     data[0][i]= gRandom->Rndm();
69   }
70   Float_t before =Mem();
71   TKDTreeIF *kdtree = new TKDTreeIF(npoints, 2, bsize, data);
72         Float_t after = Mem();
73   printf("Memory usage %f KB\n",after-before);
74         delete kdtree;
75   Float_t end = Mem();
76   printf("Memory leak %f KB\n", end-before);
77         return; 
78 }
79
80 //______________________________________________________________________
81 void TestSpeed(const Int_t npower2, const Int_t bsize)
82 {
83   //
84   // Test of building time of kdTree
85   //
86   if(npower2 < 10){
87     printf("Please specify a power of 2 greater than 10\n");
88     return;
89   }
90   
91   Int_t npoints = Int_t(pow(2., npower2))*bsize;
92   Float_t *data0 =  new Float_t[npoints*2];
93   Float_t *data[2];
94   data[0] = &data0[0];
95   data[1] = &data0[npoints];
96   for (Int_t i=0;i<npoints;i++) {
97     data[1][i]= gRandom->Rndm();
98     data[0][i]= gRandom->Rndm();
99   }
100   
101   TGraph *g = new TGraph(npower2-10);
102   g->SetMarkerStyle(7);
103   TStopwatch timer;
104   Int_t tpoints;
105   TKDTreeIF *kdtree = 0x0;
106   for(int i=10; i<npower2; i++){
107     tpoints = Int_t(pow(2., i))*bsize;
108     timer.Start(kTRUE);
109     kdtree = new TKDTreeIF(tpoints, 2, bsize, data);
110     timer.Stop();
111     g->SetPoint(i-10, i, timer.CpuTime());
112     printf("npoints [%d] nodes [%d] cpu time %f [s]\n", tpoints, kdtree->GetNNodes(), timer.CpuTime());
113     //timer.Print("u");
114     delete kdtree;
115   }
116   g->Draw("apl");
117   return;
118 }
119
120 //______________________________________________________________________
121 void TestSizeIF(Int_t nsec, Int_t nrows, Int_t npoints,  Int_t bsize, Int_t mode)
122 {
123   //
124   // Test size to build kdtree
125   //
126   Float_t before =Mem();
127   for (Int_t isec=0; isec<nsec;isec++)
128     for (Int_t irow=0;irow<nrows;irow++){
129       TestkdtreeIF(npoints,1,mode,bsize);
130     }
131   Float_t after = Mem();
132   printf("Memory usage %f\n",after-before);
133 }
134
135
136
137
138 void  TestkdtreeIF(Int_t npoints, Int_t bsize, Int_t nloop, Int_t mode)
139 {
140 //
141 // Test speed and functionality of 2D kdtree.
142 // Input parametrs:
143 // npoints - number of data points
144 // bsize   - bucket size
145 // nloop   - number of loops
146 // mode    - tasks to be performed by the kdTree
147 //         - 0  : time building the tree
148 //
149
150  
151   Float_t rangey  = 100;
152   Float_t rangez  = 100;
153   Float_t drangey = 0.1;
154   Float_t drangez = 0.1;
155
156   //
157   Float_t *data0 =  new Float_t[npoints*2];
158   Float_t *data[2];
159   data[0] = &data0[0];
160   data[1] = &data0[npoints];
161   Int_t i;   
162   for (i=0; i<npoints; i++){
163     data[0][i]          = gRandom->Uniform(-rangey, rangey);
164     data[1][i]          = gRandom->Uniform(-rangez, rangez);
165   }
166   TStopwatch timer;
167   
168   // check time build
169   printf("building kdTree ...\n");
170   timer.Start(kTRUE);
171   TKDTreeIF *kdtree = new TKDTreeIF(npoints, 2, bsize, data);
172   timer.Stop();
173   timer.Print();
174   if(mode == 0) return;
175   
176   Float_t countern=0;
177   Float_t counteriter  = 0;
178   Float_t counterfound = 0;
179   
180   if (mode ==2){
181     if (nloop) timer.Start(kTRUE);
182     Int_t res[npoints];
183     Int_t nfound = 0;
184     for (Int_t kloop = 0;kloop<nloop;kloop++){
185       if (kloop==0){
186         counteriter = 0;
187         counterfound= 0;
188         countern    = 0;
189       }
190       for (Int_t i=0;i<npoints;i++){
191         Float_t point[2]={data[0][i],data[1][i]};
192         Float_t delta[2]={drangey,drangez};
193         Int_t iter  =0;
194         nfound =0;
195         Int_t bnode =0;
196         //kdtree->FindBNode(point,delta, bnode);
197         //continue;
198         kdtree->FindInRangeA(point,delta,res,nfound,iter,bnode);
199         if (kloop==0){
200           //Bool_t isOK = kTRUE;
201           Bool_t isOK = kFALSE;
202           for (Int_t ipoint=0;ipoint<nfound;ipoint++)
203             if (res[ipoint]==i) isOK =kTRUE;
204           counteriter+=iter;
205           counterfound+=nfound;
206           if (isOK) {
207             countern++;
208           }else{
209             printf("Bug\n");
210           }
211         }
212       }
213     }
214     
215     if (nloop){
216       timer.Stop();
217       timer.Print();
218     }
219   }
220   delete [] data0;
221   
222   counteriter/=npoints;
223   counterfound/=npoints;
224   if (nloop) printf("Find nearest point:\t%f\t%f\t%f\n",countern, counteriter, counterfound);
225 }