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