]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDQADataMakerRec.cxx
HV sector finder optimized + QA histograms added
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDQADataMakerRec.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16
17 /* $Id$ */
18
19 // --- ROOT system ---
20 #include <TClonesArray.h>
21 #include <TFile.h> 
22 #include <TH1F.h> 
23 #include <TH2F.h>
24 #include <Riostream.h>
25 // --- Standard library ---
26
27 // --- AliRoot header files ---
28 #include "AliESDCaloCluster.h"
29 #include "AliESDEvent.h"
30 #include "AliQAChecker.h"
31 #include "AliLog.h"
32 #include "AliHMPIDDigit.h"
33 #include "AliHMPIDHit.h"
34 #include "AliHMPIDCluster.h"
35 #include "AliHMPIDQADataMakerRec.h"
36 #include "AliHMPIDQAChecker.h"
37 #include "AliHMPIDParam.h"
38 #include "AliHMPIDRawStream.h"
39 #include "AliLog.h"
40
41 //.
42 // HMPID AliHMPIDQADataMakerRec base class
43 // for QA of reconstruction
44 // here also errors are calculated
45 //.
46
47 ClassImp(AliHMPIDQADataMakerRec)
48            
49 //____________________________________________________________________________ 
50   AliHMPIDQADataMakerRec::AliHMPIDQADataMakerRec() : 
51   AliQADataMakerRec(AliQA::GetDetName(AliQA::kHMPID), "HMPID Quality Assurance Data Maker")
52 {
53   // ctor
54 }
55
56 //____________________________________________________________________________ 
57 AliHMPIDQADataMakerRec::AliHMPIDQADataMakerRec(const AliHMPIDQADataMakerRec& qadm) :
58   AliQADataMakerRec() 
59 {
60   //copy ctor 
61   SetName((const char*)qadm.GetName()) ; 
62   SetTitle((const char*)qadm.GetTitle()); 
63 }
64
65 //__________________________________________________________________
66 AliHMPIDQADataMakerRec& AliHMPIDQADataMakerRec::operator = (const AliHMPIDQADataMakerRec& qadm )
67 {
68   // Equal operator.
69   this->~AliHMPIDQADataMakerRec();
70   new(this) AliHMPIDQADataMakerRec(qadm);
71   return *this;
72 }
73  
74 //____________________________________________________________________________ 
75
76 void AliHMPIDQADataMakerRec::InitRecPoints()
77 {
78   // create cluster histograms in RecPoint subdir
79
80   TH1F *hCluEvt=new TH1F("CluPerEvt","Cluster multiplicity"   ,100,0,100);
81   TH1F *hCluChi2  =new TH1F("CluChi2"  ,"Chi2 "               ,1000,0,100);
82   TH1F *hCluFlg   =new TH1F("CluFlg"   ,"Cluster flag"        ,14,-1.5,12.5);
83   TH1F *hCluSize     =new TH1F("CluSize"  ,"Cluster size        ",100,0,100);
84   TH1F *hCluSizeMip  =new TH1F("CluSizeMip"  ,"Cluster size  (cluster Q > 100 ADC)   ",100,0,100);
85   TH1F *hCluQ     =new TH1F("CluQ"     ,"Cluster charge (ADC)",5000,0,5000);
86   TH1F *hCluQSect[42];
87   for(Int_t iCh =0; iCh <7; iCh++){
88    for(Int_t iSect =0; iSect < 6; iSect++){
89       hCluQSect[iCh*6+iSect] = new TH1F(Form("QClusCh%iSect%i",iCh,iSect) ,Form("Cluster charge (ADC) in Chamber %i and sector %i",iCh,iSect),5000,0,5000);
90      }  
91   }
92
93   Add2RecPointsList(hCluEvt , 0);
94   Add2RecPointsList(hCluChi2, 1);
95   Add2RecPointsList(hCluFlg , 2);
96   Add2RecPointsList(hCluSize, 3);
97   Add2RecPointsList(hCluSizeMip, 4);
98   Add2RecPointsList(hCluQ   , 5);
99   for(Int_t i=0; i< 42; i++)  Add2RecPointsList(hCluQSect[i],i+6);
100 }
101 //____________________________________________________________________________
102
103 void AliHMPIDQADataMakerRec::InitRaws()
104 {
105 //
106 // Booking QA histo for Raw data
107 //
108   const Int_t kNerr = (Int_t)AliHMPIDRawStream::kSumErr+1;
109   TH1F *hqPad[14], *hSumErr[14];
110
111   for(Int_t iddl =0; iddl<AliHMPIDRawStream::kNDDL; iddl++) {
112     hqPad[iddl] = new TH1F(Form("hqPadDDL%i",iddl), Form("Pad Q Entries at DDL %i",iddl), 500,0,5000);
113     Add2RawsList(hqPad[iddl],iddl);
114     hSumErr[iddl] = new TH1F(Form("SumErrDDL%i",iddl), Form("Error summary for ddl %i",iddl), 2*kNerr,0,2*kNerr);
115     hSumErr[iddl]->SetYTitle("%");
116     
117     for(Int_t ilabel=0; ilabel< kNerr; ilabel++) {
118       hSumErr[iddl]->GetXaxis()->CenterLabels(kTRUE);
119       //hSumErr[iddl]->GetXaxis()->SetBinLabel((2*ilabel+1),Form("%i  %s",ilabel+1,hnames[ilabel]));
120       hSumErr[iddl]->GetXaxis()->SetBinLabel((2*ilabel+1),Form("%i  %s",ilabel+1,AliHMPIDRawStream::GetErrName(ilabel)));
121       }
122       
123     Add2RawsList(hSumErr[iddl],iddl+14);
124  }
125   TH1F *hNevRaws = new TH1F("NevRaws","Events per DDL",15,0,15);
126   Add2RawsList(hNevRaws,28);
127 }
128 //____________________________________________________________________________
129 void AliHMPIDQADataMakerRec::InitESDs()
130 {
131   //
132   //Booking ESDs histograms
133    TH2F*  hCkovP  = new TH2F("CkovP" , "#theta_{c}, [rad];P, [GeV]"   , 150,      0,  7  ,100, 0, 1)   ;
134    TH2F*  hSigP   = new TH2F("SigP"  ,"#sigma_{#theta_c} [mrad];[GeV]", 150,      0,  7  ,100, 0, 1)   ;
135    TH2F*  hDifXY  = new TH2F("DifXY" ,"diff"                          , 200,    -10, 10  ,200,-10,10)  ;
136    TH2F*  hMvsP = new TH2F("MvsP","Reconstructed Mass vs P",60,0,6,1000,0,1) ;
137    TH1F*  hPid[5];
138    hPid[0] = new TH1F("PidE" ,"electron response"              , 101, -0.005,1.005)             ;
139    hPid[1] = new TH1F("PidMu","#mu response"                   , 101, -0.005,1.005)             ;
140    hPid[2] = new TH1F("PidPi","#pi response"                   , 101, -0.005,1.005)             ;
141    hPid[3] = new TH1F("PidK" ,"K response"                     , 101, -0.005,1.005)             ;
142    hPid[4] = new TH1F("PidP" ,"p response"                         ,101, -0.005,1.005)             ;
143
144 Add2ESDsList(hCkovP,0);
145 Add2ESDsList(hSigP ,1);
146 Add2ESDsList(hDifXY,2);
147 Add2ESDsList(hMvsP,3);
148 for(Int_t i=0; i< 5; i++) Add2ESDsList(hPid[i],i+4);
149
150
151 }
152 //____________________________________________________________________________
153 void AliHMPIDQADataMakerRec::MakeRaws(AliRawReader *rawReader)
154 {
155 //
156 // Filling Raws QA histos
157 //
158     AliHMPIDRawStream stream(rawReader);
159
160     while(stream.Next())
161      {
162        UInt_t ddl=stream.GetDDLNumber(); //returns 0,1,2 ... 13
163
164        for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) {
165        GetRawsData(ddl)->Fill(stream.GetChargeArray()[iPad]);}
166                                                            
167        GetRawsData(28)->Fill(ddl);
168
169        for(Int_t iErr =1; iErr<(Int_t)AliHMPIDRawStream::kSumErr; iErr++){
170  
171          Int_t numOfErr = stream.GetErrors(ddl,iErr);
172
173          GetRawsData(ddl+14)->Fill(iErr,numOfErr);
174        }
175      }
176    stream.Delete();
177 }
178 //___________________________________________________________________________
179 void AliHMPIDQADataMakerRec::MakeRecPoints(TTree * clustersTree)
180 {
181   //
182   //filling QA histos for clusters
183   //
184
185   AliHMPIDParam *pPar =AliHMPIDParam::Instance();
186   TClonesArray *clusters = new TClonesArray("AliHMPIDCluster");
187   for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){
188     TBranch *branch = clustersTree->GetBranch(Form("HMPID%d",i));
189     branch->SetAddress(&clusters);
190     branch->GetEntry(0);
191     GetRecPointsData(0)->Fill(i,clusters->GetEntries());
192     TIter next(clusters);
193     AliHMPIDCluster *clu;
194     while ( (clu = dynamic_cast<AliHMPIDCluster *>(next())) ) {
195       GetRecPointsData(1)->Fill(clu->Chi2());
196       GetRecPointsData(2)->Fill(clu->Status());
197       GetRecPointsData(3)->Fill(clu->Size());
198       if(clu->Q()>100) GetRecPointsData(3)->Fill(clu->Size());
199       GetRecPointsData(5)->Fill(clu->Q());
200       Int_t hv =  pPar->InHVSector(clu->Y()); 
201       GetRecPointsData(6+clu->Ch()*6+hv)->Fill(clu->Q());     
202     }
203   }
204
205   clusters->Delete();
206   delete clusters;
207 }
208
209 //____________________________________________________________________________
210 void AliHMPIDQADataMakerRec::MakeESDs(AliESDEvent * esd)
211 {
212   //
213   //fills QA histos for ESD
214   //
215
216   for(Int_t iTrk = 0 ; iTrk < esd->GetNumberOfTracks() ; iTrk++){
217     AliESDtrack *pTrk = esd->GetTrack(iTrk) ;
218     GetESDsData(0)->Fill(pTrk->GetP(),pTrk->GetHMPIDsignal());
219     GetESDsData(1)->Fill( pTrk->GetP(),TMath::Sqrt(pTrk->GetHMPIDchi2()));
220     Float_t xm,ym; Int_t q,np;  
221     pTrk->GetHMPIDmip(xm,ym,q,np);                       //mip info
222     Float_t xRad,yRad,th,ph;        
223     pTrk->GetHMPIDtrk(xRad,yRad,th,ph);              //track info at the middle of the radiator
224     Float_t xPc = xRad+9.25*TMath::Tan(th)*TMath::Cos(ph); // temporar: linear extrapol (B=0!)
225     Float_t yPc = yRad+9.25*TMath::Tan(th)*TMath::Sin(ph); // temporar:          "
226     GetESDsData(2)->Fill(xm-xPc,ym-yPc); //track info
227     if(pTrk->GetHMPIDsignal()>0) {
228      Double_t a = 1.292*1.292*TMath::Cos(pTrk->GetHMPIDsignal())*TMath::Cos(pTrk->GetHMPIDsignal())-1.;
229      if(a > 0) {
230     Double_t mass = pTrk->P()*TMath::Sqrt(1.292*1.292*TMath::Cos(pTrk->GetHMPIDsignal())*TMath::Cos(pTrk->GetHMPIDsignal())-1.);
231     GetESDsData(3)->Fill( pTrk->GetP(),mass);
232      }
233     }
234    Double_t pid[5] ;      pTrk->GetHMPIDpid(pid) ;
235     for(Int_t i = 0 ; i < 5 ; i++) GetESDsData(4+i)->Fill(pid[i]) ;
236   }
237 }
238 //____________________________________________________________________________
239 void AliHMPIDQADataMakerRec::StartOfDetectorCycle()
240 {
241   //Detector specific actions at start of cycle
242   
243 }
244
245 void AliHMPIDQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t , TObjArray *)
246 {
247   //Detector specific actions at end of cycle
248   // do the QA checking
249
250  //  AliQAChecker::Instance()->Run(AliQA::kHMPID, task, obj);
251 }
252