]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/mapping/macros/testMotifTypeIterators.C
Corrected list of motif types for station2 (Ivana)
[u/mrichter/AliRoot.git] / MUON / mapping / macros / testMotifTypeIterators.C
1 // $Id$
2 //
3 // Test macro for reading motif type data and iterate over them.
4
5 void testMotifTypeIterators(AliMpStationType station = kStation1,
6                             AliMpPlaneType plane = kBendingPlane)
7 {
8   TString names;
9   TString names2;
10   Int_t nv =0;
11   if ( station == kStation1 )
12     if ( plane == kBendingPlane ) 
13       names ="ABCDEFGHI";
14     else
15       names = "ABCDEFGHIJKLMN";
16   else if ( station == kStation2 ) 
17     if ( plane == kBendingPlane ) {
18       names ="ABCDEFGHIJKLMNOPQRSTUVWXY";
19       names2 ="abcdefghimnptuvvvvv";
20       nv = 5;
21     }  
22     else {
23       names = "ABCEFGHIJKLMN";
24       names2 ="abcdefgijklmnopqrstuwvvvvv";
25       nv = 5;
26     }  
27   Int_t nofMotifs = names.Length() + names2.Length(); 
28   // cout << " nofMotifs: " << nofMotifs << endl;   
29     
30   TH2C* histos[] = new TH2C* [nofMotifs];
31   TCanvas* canv[] = new TCanvas* [1+(nofMotifs-1)/4];
32   Int_t i;
33   for (i=0;i<1+(nofMotifs-1)/4;++i){
34     // canv[i] = new TCanvas(Form("canv%d",i),"Iterator viewing...");
35                // CINT limitation on DEC
36                
37     TString cname("canv"); cname += i;
38     canv[i] = new TCanvas(cname.Data(),"Iterator viewing...");
39     
40     canv[i]->Divide(2,2); 
41   }
42     
43   AliMpReader r(station, plane);
44   //r.SetVerboseLevel(2);
45
46   for (i=0;i<nofMotifs;++i){
47
48     // Get motif name
49     TString mname;
50     if (i<names.Length())
51       mname = names(i, 1);
52     else {
53       mname = names2(i-names.Length(), 1); 
54       if (mname == "v")
55         mname += i - names.Length() - (names2.Length()-nv-1);
56       else         
57         mname += "1";
58     }   
59     //if (i==36) continue;  
60         // break for these motifs (St2, BP) - to be investigated
61    
62     AliMpMotifType *mt = r.BuildMotifType(mname);
63
64     canv[i/4]->cd(1+ (i%4));
65     //histos[i] = new TH2C(Form("h%d",i),Form("Motif type %c",names[i]),
66     //                     mt->GetNofPadsX(),-0.5,mt->GetNofPadsX()-0.5,
67     //                     mt->GetNofPadsY(),-0.5,mt->GetNofPadsY()-0.5);
68                // CINT limitation on DEC
69
70     TString hname("h"); hname += i;
71
72     histos[i] = new TH2C(hname.Data(), mname.Data(),
73                          mt->GetNofPadsX(),-0.5,mt->GetNofPadsX()-0.5,
74                          mt->GetNofPadsY(),-0.5,mt->GetNofPadsY()-0.5);
75
76     cout<<"Motif Type "<<mt->GetID()<<endl;
77     cout<<"--------------------------------"<<endl;
78     Int_t num=0;
79
80     AliMpMotifTypePadIterator it = AliMpMotifTypePadIterator(mt);
81
82     for (it.First(); ! it.IsDone(); it.Next()) {
83       cout << "Iterator " << num << ' '<< it.CurrentItem().GetIndices() << endl;
84       ++num;
85       histos[i]->Fill(it.CurrentItem().GetIndices().GetFirst(),
86                       it.CurrentItem().GetIndices().GetSecond(),num);
87     }
88
89     //delete mt;
90     histos[i]->Draw("text");
91     canv[i/4]->Update();
92   }
93 }