]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/ITSSDDQAMaker.C
Coverity corrections.
[u/mrichter/AliRoot.git] / ITS / ITSSDDQAMaker.C
1 void ITSSDDQAMaker(char *iFile, Int_t MaxEvts=1000000, Int_t FirstEvt=0) {
2
3   //To have Baseline Histos uncomment parts with " // BL!!! " comment
4
5 cout << "SDD Quality Assurance Prototype Macro" << endl; 
6
7 const Int_t nSDDmodules= 260;
8 const Int_t imodoffset = 240;
9 const Int_t modtotSDD  = nSDDmodules*2;
10 const Int_t anode = 256;
11
12 Float_t xi = -0.5;
13 Float_t xf = xi + nSDDmodules;
14 TH1F *hModulePattern = new TH1F("hModulePattern","Modules pattern",nSDDmodules,xi,xf);
15 xf = xi + modtotSDD;
16 TH1F *hModuleSidePattern = new TH1F("hModuleSidePattern","Modules/Side pattern",modtotSDD,xi,xf);
17
18 TH2F *hModuleChargeMap[modtotSDD];  //260 dx e 260 sx  with A, T, Q
19 TH2F *hModuleCountsMap[modtotSDD];  //260 dx e 260 sx  with A, T, Ncounts
20 TH2I *hModuleCarlos = new TH2I("hModuleCarlos","hModuleCarlos",modtotSDD,xi,xf,101,-0.5,100.5);
21 /*
22   TH1F *hModuleSideBL[modtotSDD][anode];                                // BL!!!
23 */ 
24
25
26 //-------histograms definition 
27 Char_t *hisnam = new Char_t[50];
28 Char_t *histit = new Char_t[50];  
29 Char_t *hisnam2 = new Char_t[50];
30 Char_t *histit2 = new Char_t[50];    
31 Char_t *hisnam3 = new Char_t[50];
32  for(Int_t imod=0; imod<nSDDmodules;imod++){
33    for(Int_t isid=0;isid<2;isid++){
34      Int_t index=2*imod+isid;       //260*2 position
35
36      sprintf(hisnam,"chargeMap%d",index);
37      sprintf(histit,"Total Charge, module number %d",index);
38      hModuleChargeMap[index]=new TH2F(hisnam,histit,256,-0.5,255.5,256,-0.5,255.5);  
39
40      sprintf(hisnam2,"countsMap%d",index);
41      sprintf(histit2,"Number of Counts, module number %d",index);
42      hModuleCountsMap[index] = new TH2F(hisnam2,histit2,256,-0.5,255.5,256,-0.5,255.5);
43      /*
44        for(Int_t ianode=0; ianode<anode; ianode++){                         // BL!!!
45        sprintf(hisnam3,"BL_module_%d_%d",index,ianode); 
46        //cout<<hisnam3 <<endl;
47        hModuleSideBL[index][ianode] = new TH1F(hisnam3,hisnam3,256,0.,1024.);
48      }
49      */
50    }
51  }
52
53   TString strFile = iFile;
54   strFile += "?EventType=7";
55   AliRawReader *rd = new AliRawReaderDate(strFile.Data(),FirstEvt);  // open run
56   Int_t evCounter = 0;
57   do{                       // start loop on events
58     if(++evCounter > MaxEvts) { cout << MaxEvts << " events read, stop" << endl; evCounter--; break; }  
59     cout << "Read Event: " << evCounter+FirstEvt-1 << endl;
60
61     rd->RequireHeader(kFALSE);             
62     rd->Reset();                           // reset the current position to the beginning of the event
63  
64     Int_t nSkip = 0;                     // number of skipped signals
65     AliITSRawStreamSDD s(rd);            //This class provides access to ITS SDD digits in raw data.
66     Int_t iddl;
67     Int_t isddmod;
68     Int_t moduleSDD;
69     gStyle->SetPalette(1);
70     while(s.Next()){                     //read the next raw digit; returns kFALSE if there is no digit left
71       if(s.IsCompletedModule()) continue;
72       if(s.IsCompletedDDL()) continue;
73       iddl=rd->GetDDLID()-2; // -2 is temporary for test raw data
74
75         isddmod=s.GetModuleNumber(iddl,s.GetCarlosId());        //this is the FEE Carlos
76         //cout<<"DDLID= "<<iddl <<"; Module number= " <<isddmod  <<endl;
77         if(isddmod >= imodoffset) { 
78           hModulePattern->Fill(isddmod-imodoffset);             // 0 to 259    so  240 to 499
79           moduleSDD=2*(isddmod-imodoffset)+s.GetChannel();
80           hModuleSidePattern->Fill(moduleSDD);                  // 0 to 519
81           hModuleCarlos->Fill(isddmod-imodoffset,s.GetCarlosId());
82           //cout << "anode " << s.GetCoord1() << ", time bin: " << s.GetCoord2() << ", charge: " << s.GetSignal() << endl;        
83           Int_t coord1 = s.GetCoord1();
84           Int_t coord2 = s.GetCoord2();
85           Int_t signal = s.GetSignal();
86           hModuleChargeMap[moduleSDD]->Fill(coord2, coord1,signal);
87           hModuleCountsMap[moduleSDD]->Fill(coord2, coord1 );   
88           //hModuleSideBL[moduleSDD][coord1]->Fill(signal);             // BL  !!!
89         } else {
90           nSkip++;
91         }
92     }    
93     cout << "End of Event " << evCounter+FirstEvt-1 << ", " << nSkip << " wrong module numbers" << endl;
94   } while(rd->NextEvent()); // end loop on events
95   delete rd;
96
97   cout << "end after " << evCounter << " events" << endl;
98   /*
99   TNtuple *Baseline = new TNtuple("Baseline","Baseline","HalfModule:Anode:Mean:RMS");      // BL!!!
100   Float_t meanBL;
101   Float_t rmsBL;  
102   */
103   for(Int_t i=0; i<modtotSDD; i++){   
104     if(hModuleSidePattern->GetBinContent(i+1)){              //check if they're not empty
105       hModuleChargeMap[i]->GetXaxis()->SetTitle("Time Bin");
106       hModuleChargeMap[i]->GetYaxis()->SetTitle("Anode"); 
107       hModuleCountsMap[i]->GetXaxis()->SetTitle("Time Bin");
108       hModuleCountsMap[i]->GetYaxis()->SetTitle("Anode");  
109       /*
110       for(Int_t ianode=0; ianode<anode; ianode++ ){                                      // BL!!!
111         hModuleSideBL[i][ianode]->GetXaxis()->SetTitle("ADC counts");
112         hModuleSideBL[i][ianode]->GetYaxis()->SetTitle("#"); 
113         meanBL = hModuleSideBL[i][ianode]->GetMean();
114         rmsBL = hModuleSideBL[i][ianode]->GetRMS();
115         gaussfitBL = hModuleSideBL[i][ianode]->Fit("gaus");     
116         Baseline->Fill(i,ianode,meanBL,rmsBL);
117       }
118       */
119     }
120   }  
121   
122   hModuleSidePattern->GetXaxis()->SetTitle("2*(Module Number-1)+Side");
123   hModuleSidePattern->GetYaxis()->SetTitle("Counts");  
124   hModulePattern->GetXaxis()->SetTitle("Module Number");  
125   hModulePattern->GetYaxis()->SetTitle("Counts");  
126
127
128   //-------store Histograms
129   cout << "Store Histograms" << endl;
130   TString oFileName(iFile);
131   oFileName.Append(".root");
132   TFile *oFile = TFile::Open(oFileName,"recreate");
133   hModulePattern->Write();
134   hModuleSidePattern->Write();
135   hModuleCarlos->Write();
136   for(Int_t i=0; i<modtotSDD; i++){ 
137       if(hModuleSidePattern->GetBinContent(i+1)){            //check if they're not empty
138         hModuleChargeMap[i]->Write();
139         hModuleCountsMap[i]->Write();     
140         /* 
141         for(Int_t ianode=0; ianode<anode; ianode++ ){                           // BL!!!
142           hModuleSideBL[i][ianode]->Write();
143           Baseline->Write();
144         }
145         */
146       }
147   }
148   
149   oFile->Close();
150   cout << "Clear memory" << endl;
151   for(Int_t imod=0; imod<nSDDmodules;imod++){
152     for(Int_t isid=0;isid<2;isid++){
153       Int_t index=2*imod+isid;       //260*2 position
154       delete hModuleChargeMap[index];
155       delete hModuleCountsMap[index];       
156       /*
157       for(Int_t ianode=0; ianode<anode; ianode++ ){                              // BL!!!
158         delete hModuleSideBL[index][ianode]; 
159         delete Baseline;
160       }
161       */
162     }
163   }
164   delete hModulePattern;
165   delete hModuleSidePattern;
166   delete hModuleCarlos;
167
168 }