]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDQADataMakerRec.cxx
input for DA from Michal
[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 <TProfile.h>
25 #include <Riostream.h>
26 // --- Standard library ---
27
28 // --- AliRoot header files ---
29 #include "AliESDCaloCluster.h"
30 #include "AliESDEvent.h"
31 #include "AliQAChecker.h"
32 #include "AliLog.h"
33 #include "AliHMPIDDigit.h"
34 #include "AliHMPIDHit.h"
35 #include "AliHMPIDCluster.h"
36 #include "AliHMPIDQADataMakerRec.h"
37 #include "AliHMPIDQAChecker.h"
38 #include "AliHMPIDParam.h"
39 #include "AliHMPIDRawStream.h"
40 #include "AliLog.h"
41
42 //.
43 // HMPID AliHMPIDQADataMakerRec base class
44 // for QA of reconstruction
45 // here also errors are calculated
46 //.
47
48 ClassImp(AliHMPIDQADataMakerRec)
49            
50 //____________________________________________________________________________ 
51   AliHMPIDQADataMakerRec::AliHMPIDQADataMakerRec() : 
52   AliQADataMakerRec(AliQA::GetDetName(AliQA::kHMPID), "HMPID Quality Assurance Data Maker")
53 {
54   fEvtRaw=0;
55   // ctor
56 }
57
58 //____________________________________________________________________________ 
59 AliHMPIDQADataMakerRec::AliHMPIDQADataMakerRec(const AliHMPIDQADataMakerRec& qadm) :
60   AliQADataMakerRec() 
61 {
62   //copy ctor 
63   SetName((const char*)qadm.GetName()) ; 
64   SetTitle((const char*)qadm.GetTitle()); 
65 }
66
67 //__________________________________________________________________
68 AliHMPIDQADataMakerRec& AliHMPIDQADataMakerRec::operator = (const AliHMPIDQADataMakerRec& qadm )
69 {
70   // Equal operator.
71   this->~AliHMPIDQADataMakerRec();
72   new(this) AliHMPIDQADataMakerRec(qadm);
73   return *this;
74 }
75  
76 //____________________________________________________________________________ 
77
78 void AliHMPIDQADataMakerRec::InitRecPoints()
79 {
80   // create cluster histograms in RecPoint subdir
81
82   TProfile *hCluMult          = new TProfile("CluMult"   ,"Cluster multiplicity per chamber"    , 16, -1 , 7  , 0, 500);
83   Add2RecPointsList(hCluMult    , 0);
84   
85   TH2F *hCluFlg          = new TH2F("CluFlg"      ,"Cluster flag  "                              ,  56  ,-1.5, 12.5, 70, -0.5, 6.5);
86   Add2RecPointsList(hCluFlg    , 1);
87   
88   TH1F *hCluSizeMip[7], *hCluSizePho[7];
89   
90   TH1F *hCluQSect[42], *hCluQSectZoom[42];
91   
92   for(Int_t iCh =0; iCh <7; iCh++){
93    hCluSizeMip[iCh] = new TH1F(Form("CluSizeMipCh%i",iCh),Form("Cluster size  MIP  (cluster Q > 100 ADC) in Chamber %i",iCh),  50  , 0  , 50  );
94   Add2RecPointsList(hCluSizeMip[iCh], iCh+2);
95
96    hCluSizePho[iCh]  = new TH1F(Form("CluSizePho%i",iCh ),Form("Cluster size  Phots(cluster Q < 100 ADC) in Chamber %i",iCh),  50  , 0  , 50  );
97   Add2RecPointsList(hCluSizePho[iCh], iCh+7+2);
98
99     for(Int_t iSect =0; iSect < 6; iSect++){
100       hCluQSectZoom[iCh*6+iSect] = new TH1F(Form("QClusCh%iSect%iZoom",iCh,iSect) ,Form("Zoom on Cluster charge (ADC) in Chamber %i and sector %i",iCh,iSect),100,0,100);
101       Add2RecPointsList(hCluQSectZoom[iCh*6+iSect],2+14+iCh*6+iSect);
102
103       hCluQSect[iCh*6+iSect] = new TH1F(Form("QClusCh%iSect%i",iCh,iSect) ,Form("Cluster charge (ADC) in Chamber %i and sector %i",iCh,iSect),250,0,5000);
104       Add2RecPointsList(hCluQSect[iCh*6+iSect],2+14+42+iCh*6+iSect);
105
106     }  
107   }
108 }
109 //____________________________________________________________________________
110
111 void AliHMPIDQADataMakerRec::InitRaws()
112 {
113 //
114 // Booking QA histo for Raw data
115 //
116   const Int_t kNerr = (Int_t)AliHMPIDRawStream::kSumErr+1;
117   TH1F *hSumErr[14];
118   TH2F *hDilo[14];
119
120   for(Int_t iddl =0; iddl<AliHMPIDRawStream::kNDDL; iddl++) {
121     
122     hSumErr[iddl] = new TH1F(Form("SumErrDDL%i",iddl), Form("Error summary for ddl %i",iddl), 2*kNerr,0,2*kNerr);
123     for(Int_t ilabel=0; ilabel< kNerr; ilabel++) {
124       hSumErr[iddl]->GetXaxis()->CenterLabels(kTRUE);
125       hSumErr[iddl]->GetXaxis()->SetBinLabel((2*ilabel+1),Form("%i  %s",ilabel+1,AliHMPIDRawStream::GetErrName(ilabel)));
126       }
127       
128     Add2RawsList(hSumErr[iddl],iddl);
129     
130     hDilo[iddl] = new TH2F(Form("DiloDDL%i",iddl),Form("Dilogic response at DDL",iddl),24,0,24,10,0,10);
131     
132       Add2RawsList(hDilo[iddl],14+iddl);
133  }
134 }
135 //____________________________________________________________________________
136 void AliHMPIDQADataMakerRec::InitESDs()
137 {
138   //
139   //Booking ESDs histograms
140    TH2F*  hCkovP  = new TH2F("CkovP" , "#theta_{c}, [rad];P, [GeV]"   , 150,      0,  7  ,100, 0, 1)   ;
141    TH2F*  hSigP   = new TH2F("SigP"  ,"#sigma_{#theta_c} [mrad];[GeV]", 150,      0,  7  ,100, 0, 1)   ;
142    TH2F*  hDifXY  = new TH2F("DifXY" ,"diff"                          , 200,    -10, 10  ,200,-10,10)  ;
143    TH2F*  hMvsP = new TH2F("MvsP","Reconstructed Mass vs P",60,0,6,1000,0,1) ;
144    TH1F*  hPid[5];
145    hPid[0] = new TH1F("PidE" ,"electron response"              , 101, -0.005,1.005)             ;
146    hPid[1] = new TH1F("PidMu","#mu response"                   , 101, -0.005,1.005)             ;
147    hPid[2] = new TH1F("PidPi","#pi response"                   , 101, -0.005,1.005)             ;
148    hPid[3] = new TH1F("PidK" ,"K response"                     , 101, -0.005,1.005)             ;
149    hPid[4] = new TH1F("PidP" ,"p response"                         ,101, -0.005,1.005)             ;
150
151 Add2ESDsList(hCkovP,0);
152 Add2ESDsList(hSigP ,1);
153 Add2ESDsList(hDifXY,2);
154 Add2ESDsList(hMvsP,3);
155 for(Int_t i=0; i< 5; i++) Add2ESDsList(hPid[i],i+4);
156
157
158 }
159 //____________________________________________________________________________
160 void AliHMPIDQADataMakerRec::MakeRaws(AliRawReader *rawReader)
161 {
162 //
163 // Filling Raws QA histos
164 //
165     AliHMPIDRawStream stream(rawReader);
166
167     fEvtRaw++;
168
169     while(stream.Next())
170      {
171        UInt_t ddl=stream.GetDDLNumber(); //returns 0,1,2 ... 13
172  
173        for(Int_t iErr =1; iErr<(Int_t)AliHMPIDRawStream::kSumErr; iErr++){
174  
175          Int_t numOfErr = stream.GetErrors(ddl,iErr);
176          GetRawsData(ddl)->Fill(iErr,numOfErr);
177        }
178        UInt_t word; Int_t Nddl, r, d, a;      
179        for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) {
180        AliHMPIDDigit dig(stream.GetPadArray()[iPad],stream.GetChargeArray()[iPad]);
181        dig.Raw(word,Nddl,r,d,a);
182        GetRawsData(ddl+14)->Fill(r,d);
183        }
184        
185      }
186    stream.Delete();
187    
188 }
189 //___________________________________________________________________________
190 void AliHMPIDQADataMakerRec::MakeRecPoints(TTree * clustersTree)
191 {
192   //
193   //filling QA histos for clusters
194   //
195   AliHMPIDParam *pPar =AliHMPIDParam::Instance();
196   TClonesArray *clusters = new TClonesArray("AliHMPIDCluster");
197   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
198     TBranch *branch = clustersTree->GetBranch(Form("HMPID%d",iCh));
199     branch->SetAddress(&clusters);
200     branch->GetEntry(0);
201     GetRecPointsData(0)->Fill(iCh,clusters->GetEntries());
202     TIter next(clusters);
203     AliHMPIDCluster *clu;
204     while ( (clu = dynamic_cast<AliHMPIDCluster *>(next())) ) {
205       GetRecPointsData(1)->Fill(clu->Status(),iCh);
206       Int_t sect =  pPar->InHVSector(clu->Y());
207       if(clu->Q()>100) GetRecPointsData(2+iCh)->Fill(clu->Size());
208       else {
209         GetRecPointsData(2+7+iCh)->Fill(clu->Size());
210         GetRecPointsData(2+14+iCh*6+sect)->Fill(clu->Q());
211       }    
212       GetRecPointsData(2+14+42+iCh*6+sect)->Fill(clu->Q());
213     }
214   }
215
216   clusters->Delete();
217   delete clusters;
218 }
219
220 //____________________________________________________________________________
221 void AliHMPIDQADataMakerRec::MakeESDs(AliESDEvent * esd)
222 {
223   //
224   //fills QA histos for ESD
225   //
226
227   for(Int_t iTrk = 0 ; iTrk < esd->GetNumberOfTracks() ; iTrk++){
228     AliESDtrack *pTrk = esd->GetTrack(iTrk) ;
229     GetESDsData(0)->Fill(pTrk->GetP(),pTrk->GetHMPIDsignal());
230     GetESDsData(1)->Fill( pTrk->GetP(),TMath::Sqrt(pTrk->GetHMPIDchi2()));
231     Float_t xm,ym; Int_t q,np;  
232     pTrk->GetHMPIDmip(xm,ym,q,np);                       //mip info
233     Float_t xRad,yRad,th,ph;        
234     pTrk->GetHMPIDtrk(xRad,yRad,th,ph);              //track info at the middle of the radiator
235     Float_t xPc = xRad+9.25*TMath::Tan(th)*TMath::Cos(ph); // temporar: linear extrapol (B=0!)
236     Float_t yPc = yRad+9.25*TMath::Tan(th)*TMath::Sin(ph); // temporar:          "
237     GetESDsData(2)->Fill(xm-xPc,ym-yPc); //track info
238     if(pTrk->GetHMPIDsignal()>0) {
239      Double_t a = 1.292*1.292*TMath::Cos(pTrk->GetHMPIDsignal())*TMath::Cos(pTrk->GetHMPIDsignal())-1.;
240      if(a > 0) {
241     Double_t mass = pTrk->P()*TMath::Sqrt(1.292*1.292*TMath::Cos(pTrk->GetHMPIDsignal())*TMath::Cos(pTrk->GetHMPIDsignal())-1.);
242     GetESDsData(3)->Fill( pTrk->GetP(),mass);
243      }
244     }
245    Double_t pid[5] ;      pTrk->GetHMPIDpid(pid) ;
246     for(Int_t i = 0 ; i < 5 ; i++) GetESDsData(4+i)->Fill(pid[i]) ;
247   }
248 }
249 //____________________________________________________________________________
250 void AliHMPIDQADataMakerRec::StartOfDetectorCycle()
251 {
252   //Detector specific actions at start of cycle
253   
254 }
255
256 void AliHMPIDQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray *histos)
257 {
258   //Detector specific actions at end of cycle
259   // do the QA checking
260  //  AliQAChecker::Instance()->Run(AliQA::kHMPID, task, obj);
261
262   if(task==AliQA::kRAWS) {
263     for(Int_t iddl=0;iddl<14;iddl++) {
264      TH1F *h = (TH1F*)histos->At(14+iddl); //ddl histos scaled by the number of events 
265      h->Scale(1./(Float_t)fEvtRaw);
266     }
267   }
268 }
269