]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/Hqa.C
e00d08562fef2109a7082485e7466476ca0ac519
[u/mrichter/AliRoot.git] / HMPID / Hqa.C
1 #include <TSystem.h>
2 #include <TFile.h>
3 #include <TTree.h>
4 #include <TClonesArray.h>
5 #include <TObjArray.h>
6 #include <TH2F.h>
7 #include <TCanvas.h>
8
9 #include <AliHMPIDHit.h>
10 #include <AliHMPIDCluster.h>
11
12 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
13 TObjArray *CreateContainer(const char *classname,TTree *pTree)
14 {
15   TObjArray *pOA=new TObjArray(7); pOA->SetOwner();
16   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
17     TClonesArray *pCA=new TClonesArray(classname);
18     pOA->AddAt(pCA,iCh);    
19     pTree->SetBranchAddress(Form("HMPID%i",iCh),&pCA);
20   }
21   return pOA;
22 }
23 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
24 TH1F *hHitQdc; TH2F *hHitMap[7];  
25 void Hits(Int_t mode,TTree *pTree=0x0)
26 {
27   switch(mode){
28     case 1:
29       hHitQdc=new TH1F("HitQdc","Hit Qdc all chamber;QDC",500,0,4000);
30       for(Int_t iCh=0;iCh<7;iCh++) hHitMap[iCh]=new TH2F(Form("HitMap%i",iCh),Form("Ch%i;x_{Hit};y_{Hit}",iCh),160,0,160,160,0,160);      
31       return;
32     case 2:
33       if(pTree==0) return;
34       TClonesArray *pHits=new TClonesArray("AliHMPIDHit");  pTree->SetBranchAddress("HMPID",&pHits);  
35       for(Int_t iEnt=0;iEnt<pTree->GetEntriesFast();iEnt++){//entries loop
36         pTree->GetEntry(iEnt);
37         for(Int_t iHit=0;iHit<pHits->GetEntriesFast();iHit++){//hits loop
38           AliHMPIDHit *pHit = (AliHMPIDHit*)pHits->UncheckedAt(iHit);
39           hHitMap[pHit->Ch()]->Fill(pHit->LorsX(),pHit->LorsY());
40           hHitQdc->Fill(pHit->Q());
41         }//hits loop      
42       }//entries loop
43       delete pHits;
44       return;
45     case 3:
46       TCanvas *c1=new TCanvas("HitCan","Hits",1280,800); c1->Divide(3,3);
47   
48       for(Int_t iCh=0;iCh<7;iCh++){
49         if(iCh==6) c1->cd(1); if(iCh==5) c1->cd(2);
50         if(iCh==4) c1->cd(4); if(iCh==3) c1->cd(5); if(iCh==2) c1->cd(6);
51                               if(iCh==1) c1->cd(8); if(iCh==0) c1->cd(9);
52         hHitMap[iCh]->Draw();
53       }  
54       c1->cd(3); gPad->SetLogy(); hHitQdc->Draw();
55       return;
56   }
57 }//Hits()
58 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
59 TH1F *hCluEvt,*hCluChi2,*hCluFlg,*hCluSize;
60 void Clus(Int_t mode, TTree *pTree=0x0)
61 {
62   switch(mode){
63     case 1:
64       hCluEvt=new TH1F("CluPerEvt","# clusters per event",11,-0.5,10.5);
65       hCluChi2  =new TH1F("CluChi2"  ,"Chi2 "               ,1000,0,100);
66       hCluFlg   =new TH1F("CluFlg"   ,"Cluster flag"        ,14,-1.5,12.5);                       hCluFlg->SetFillColor(5);
67       hCluSize  =new TH1F("CluSize"  ,"Raw cluster size    ",100,0,100);
68       return;
69     case 2:      
70       if(pTree==0) return;
71       TObjArray *pLst=CreateContainer("AliHMPIDCluster",pTree); pTree->GetEntry(0);
72       for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){//chambers loop
73         TClonesArray *pClus=(TClonesArray *)pLst->UncheckedAt(iCh);
74         hCluEvt->Fill(pClus->GetEntriesFast());
75         for(Int_t iClu=0;iClu<pClus->GetEntriesFast();iClu++){//clusters loop
76           AliHMPIDCluster *pClu=(AliHMPIDCluster*)pClus->UncheckedAt(iClu);
77           hCluFlg->Fill(pClu->Status());  hCluChi2->Fill(pClu->Chi2());  hCluSize->Fill(pClu->Size());
78         }
79       }
80       delete pLst;           
81       return;  
82     case 3:
83       TCanvas *c1=new TCanvas("CluComCan","Clusters in common",1280,800); c1->Divide(3,3);
84       c1->cd(1); hCluEvt->Draw();       c1->cd(2); hCluChi2->Draw(); 
85       c1->cd(3); hCluFlg->Draw();       c1->cd(4); hCluSize->Draw(); 
86       return;
87   }//switch
88 }//Clus()
89 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
90 void Hqa()
91 {
92
93   TFile *fh=0; if(gSystem->IsFileInIncludePath("HMPID.Hits.root"))      fh=TFile::Open("HMPID.Hits.root"     ,"read");if(fh->IsZombie()) fh=0;
94   TFile *fd=0; if(gSystem->IsFileInIncludePath("HMPID.Digits.root"))    fd=TFile::Open("HMPID.Digits.root"   ,"read");if(fd->IsZombie()) fd=0;
95   TFile *fc=0; if(gSystem->IsFileInIncludePath("HMPID.RecPoints.root")) fc=TFile::Open("HMPID.RecPoints.root","read");if(fc->IsZombie()) fc=0;
96   if(fh==0 && fd==0 && fc==0){Printf("Nothing to do!"); return;}
97   if(fh) Hits(1); if(fc) Clus(1);  //book
98   Int_t iEvt=0;
99   while(1){
100     TTree *th=0; if(fh) th=(TTree*)fh->Get(Form("Event%i/TreeH",iEvt));
101     TTree *td=0; if(fd) td=(TTree*)fd->Get(Form("Event%i/TreeD",iEvt));
102     TTree *tc=0; if(fc) tc=(TTree*)fc->Get(Form("Event%i/TreeR",iEvt));
103     Hits(2,th);   Clus(2,tc); //fill
104     if(th==0 && td==0 && tc==0) break;
105     iEvt++;
106     Printf("Event %i processed",iEvt);
107   }
108   if(fh) Hits(3);  if(fc) Clus(3); //plot
109 }
110
111 /*
112
113 void BlahQA()
114 {
115   TFile *pFile = new TFile("Hqa.root","recreate");
116   AliHMPIDCluster::DoCorrSin(kFALSE);
117   AliHMPIDDigit::fSigmas=4;
118   Int_t nEvts=10000;
119   TLegend *lQ=new TLegend(0.5,0.5,0.9,0.9);
120
121   TH1F *hQ7eV  =new TH1F("hQ7eV"  ,"" ,300,-50,2000);     hQ7eV  ->SetLineColor(kRed);      lQ->AddEntry(hQ7eV  ,"Ckov 7 eV");  hQ7eV->SetStats(0);
122   TH1F *hQ200eV=new TH1F("hQ200eV","" ,300,-50,2000);     hQ200eV->SetLineColor(kBlack);    lQ->AddEntry(hQ200eV,"mip 200 eV");
123   TH1F *hQ500eV=new TH1F("hQ500eV","" ,300,-50,2000);     hQ500eV->SetLineColor(kCyan);     lQ->AddEntry(hQ500eV,"mip 500 eV");
124   TH1F *hQ900eV=new TH1F("hQ900eV","" ,300,-50,2000);     hQ900eV->SetLineColor(kGreen);    lQ->AddEntry(hQ900eV,"mip 900 eV");
125   
126   
127   gStyle->SetOptStat(10);
128   TH2F *pCluMapSi1  =new TH2F("cluMapSi1","Size 1 map"       ,1700,-10,160,1700,-10,160);
129   TH2F *pCluMapLo0  =new TH2F("cluMNoLo0","Loc Max 0 map"    ,1700,-10,160,1700,-10,160);     
130   TH2F *pCluMapLo1  =new TH2F("cluMapLo1","Loc Max 1 map"    ,1700,-10,160,1700,-10,160);        
131   TH2F *pCluMapUnf  =new TH2F("cluMapUnf","Unfolded map"     ,1700,-10,160,1700,-10,160);         
132   TH2F *pCluMapEdg  =new TH2F("cluMapEdg","On edge map"      ,1700,-10,160,1700,-10,160);         
133   TH2F *pCluMapCoG  =new TH2F("cluMapCoG","CoG  map"         ,1700,-10,160,1700,-10,160);            
134   TH2F *pCluMapEmp  =new TH2F("cluMapEmp","undefined-empty"  ,1700,-10,160,1700,-10,160);      
135   TH2F *pCluMapNoLoc=new TH2F("cluMapNoLoc","no loc maxima"  ,1700,-10,160,1700,-10,160);      
136   TH2F *pCluMapAbn  =new TH2F("cluMapAbn","abnormal fit"     ,1700,-10,160,1700,-10,160);      
137   TH2F *pCluMapNot  =new TH2F("cluMapNot","Raw Clusters"     ,1700,-10,160,1700,-10,160);      
138   TH2F *pCluMapMax  =new TH2F("cluMapMax","N. locs excceds"  ,1700,-10,160,1700,-10,160);      
139
140   TH1F *hHitCluDifX = new TH1F("hHitCluDifX" ,";entries;x_{Hit}-x_{Clu} [cm]"   ,1000,-1,1);          hHitCluDifX->Sumw2();    hHitCluDifX->SetFillColor(kYellow);
141 //  TH2F *hHitCluDifXv= new TH2F("hHitCluDifXv",";x_{Hit};x_{Hit}-x_{Clu} [cm]"   ,500,-0.5,0.5,1000,-0.2,0.2);hHitCluDifXv->Sumw2();
142   TProfile *hHitCluDifXv= new TProfile("hHitCluDifXv",";x_{Hit};x_{Hit}-x_{Clu} [cm]"   ,500,-0.5,0.5);
143   TH1F *hHitCluDifY = new TH1F("hHitCluDifY" ,";entries;y_{Hit}-y_{Clu} [cm]"   ,1000,-1,1);          hHitCluDifY->Sumw2();    hHitCluDifY->SetFillColor(kYellow);
144   TH2F *hHitCluDifXY= new TH2F("hHitCluDifXY",";x_{Hit}-x_{Clu};y_{Hit}-y_{Clu}",1000,-1,1,1000,-1,1);hHitCluDifXY->Sumw2();  
145   TH1F *hHitCluDifQ = new TH1F("hHitCluDifQ" ,";entries;(Q_{Clu}-Q_{Hit})/Q_{Hit}" ,200 ,-200,200);   hHitCluDifQ->Sumw2();    hHitCluDifQ->SetFillColor(kYellow);
146   
147  
148   Float_t e200=200e-9,e500=500e-9,e900=900e-9,e7=7e-9;//predefined  Eloss
149   
150   AliHMPIDHit hit(0,0,kProton,0,0,0);
151   for(int i=0;i<5000;i++){
152     hQ200eV->Fill(hit.QdcTot(e200));  hQ500eV->Fill(hit.QdcTot(e500));   hQ900eV->Fill(hit.QdcTot(e900));  hQ7eV->Fill(hit.QdcTot(e7));
153   }  
154   TClonesArray hits("AliHMPIDHit");  TClonesArray sdigs("AliHMPIDDigit");
155   TObjArray digs(7); for(Int_t i=0;i<7;i++) digs.AddAt(new TClonesArray("AliHMPIDDigit"),i);
156   TObjArray clus(7); for(Int_t i=0;i<7;i++) clus.AddAt(new TClonesArray("AliHMPIDCluster"),i);
157   
158     
159   for(Int_t iEvt=0;iEvt<nEvts;iEvt++){//events loop
160     if(iEvt%50==0)Printf("============> iEvt = %d ",iEvt);
161     
162     Int_t ch,pid; Float_t e,hitx,hity,hitq;
163 //    Int_t nHits=(type==999)?1:40;
164     Int_t nHits=1;
165     for(Int_t iHit=0;iHit<nHits;iHit++){//hits loop for the current event
166       switch(iHit){
167         case 0:  ch=0;pid=kProton;e=e200;
168                  hitx=AliHMPIDParam::SizePadX()*(6+gRandom->Rndm());
169                  hity=AliHMPIDParam::SizePadY()*(6+gRandom->Rndm());break; //mip ramdomly distributed in one pad
170         case 1:  ch=0;pid=kProton;e=e200;hitx=0.4;hity=0.42;break; //mip in left-hand bottom coner of chamber 0
171         case 2:  ch=0;pid=kProton;e=e200;hitx=0.4;hity=30  ;break; //mip on left edge of chamber 0
172         case 3:  ch=0;pid=kProton;e=e200;hitx=40; hity=0.42;break; //mip on bottom edge of chamber 0
173         default: ch=gRandom->Rndm()*6; pid=(gRandom->Rndm()>0.9)? kProton:50000050;
174                   if(pid==kProton) 
175                     e=gRandom->Rndm()*900e-9; 
176                   else
177                     e=5.5e-9+3e-9*gRandom->Rndm();
178                   hitx=gRandom->Rndm()*AliHMPIDParam::SizeAllX(); hity=gRandom->Rndm()*AliHMPIDParam::SizeAllY();break; //random hit
179       }
180       new(hits[iHit]) AliHMPIDHit(ch,e,pid,iHit,hitx,hity);                          
181       hitq=e;
182     }//hits loop
183     
184     AliHMPIDv1::Hit2Sdi(&hits,&sdigs);
185     AliHMPIDDigitizer::DoNoise(kFALSE);
186     AliHMPIDDigitizer::Sdi2Dig(&sdigs,&digs);
187     AliHMPIDReconstructor::Dig2Clu(&digs,&clus,kTRUE);
188
189 // From here normal procedure for QA
190
191     Int_t kMaxCh=(nHits==1)?0:AliHMPIDParam::kMaxCh;
192     for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=kMaxCh;iCh++){//chambers loop
193       TClonesArray *pDigs=(TClonesArray *)digs.UncheckedAt(iCh);
194       TClonesArray *pClus=(TClonesArray *)clus.UncheckedAt(iCh);
195       
196 //      if(pClus->GetEntriesFast()>nHits) {hits.Print();}
197 //      if(pClus->GetEntriesFast()>nHits) {pDigs->Print();Printf("----S D I G I T S-------------");}
198 //      if(pClus->GetEntriesFast()>nHits) {sdigs.Print();Printf("-------------------------------");}
199       hCluPerEvt->Fill(pClus->GetEntriesFast());
200       for(Int_t iClu=0;iClu<pClus->GetEntriesFast();iClu++){//clusters loop
201         AliHMPIDCluster *pClu=(AliHMPIDCluster*)pClus->UncheckedAt(iClu);
202         Float_t clux=pClu->X(); Float_t cluy=pClu->Y(); Float_t cluq=pClu->Q();
203         hCluFlg->Fill(pClu->Status());
204         hCluChi2->Fill(pClu->Chi2());
205         hCluRawSize->Fill(pClu->Size());
206          
207         switch(pClu->Status()){
208             case AliHMPIDCluster::kSi1:   pCluMapSi1->Fill(clux,cluy);  break;
209             case AliHMPIDCluster::kLo1:   pCluMapLo1->Fill(clux,cluy);  break;
210             case AliHMPIDCluster::kUnf:   pCluMapUnf->Fill(clux,cluy);  break; 
211             case AliHMPIDCluster::kMax:   pCluMapMax->Fill(clux,cluy);  break;
212             case AliHMPIDCluster::kEdg:   pCluMapEdg->Fill(clux,cluy);  break;
213             case AliHMPIDCluster::kCoG:   pCluMapCoG->Fill(clux,cluy);  break;
214             case AliHMPIDCluster::kNoLoc: pCluMapNoLoc->Fill(clux,cluy);break;
215             case AliHMPIDCluster::kAbn:   pCluMapAbn->Fill(clux,cluy);  break;
216             case AliHMPIDCluster::kNot:   pCluMapNot->Fill(clux,cluy);  break;
217             default:     pCluMapEmp->Fill(clux,cluy); break;            
218         }
219         
220         hHitCluDifX->Fill(hitx-clux); hHitCluDifY->Fill(hity-cluy); hHitCluDifXY->Fill(hitx-clux,hity-cluy); hHitCluDifQ->Fill(100*(cluq-hitq)/hitq);
221         // distorsion due to feedback photons
222         Int_t pc,px,py;
223         AliHMPIDParam::Lors2Pad(hitx,hity,pc,px,py);
224         Float_t padCenterX = AliHMPIDParam::LorsX(pc,px);
225         if(pClu->Size()>1)hHitCluDifXv->Fill(hitx-padCenterX,(hitx-clux));        
226         //
227       }//clusters loop
228     }//chambers loop
229       
230     hits.Delete();  sdigs.Delete();  for(int i = 0;i<7;i++){((TClonesArray*)digs.At(i))->Delete();((TClonesArray*)clus.At(i))->Delete();}
231   }//events loop      
232 }
233
234
235
236 */