]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDQualAssDataMaker.cxx
Changing name QualAss to QA
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDQualAssDataMaker.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 //---
20 //  Produces the data needed to calculate the quality assurance. 
21 //  All data must be mergeable objects.
22 //  A. Mastroserio
23 //---
24
25 // --- ROOT system ---
26 #include <TClonesArray.h>
27 #include <TFile.h> 
28 #include <TH1F.h> 
29 #include <TH2F.h>
30 #include <TH1I.h> 
31 #include <TDirectory.h>
32 #include <Riostream.h>
33 // --- Standard library ---
34
35 // --- AliRoot header files ---
36 #include "AliESDCaloCluster.h"
37 #include "AliESDEvent.h"
38 #include "AliLog.h"
39 #include "AliHMPIDDigit.h"
40 #include "AliHMPIDHit.h"
41 #include "AliHMPIDCluster.h"
42 #include "AliHMPIDQualAssDataMaker.h"
43
44 ClassImp(AliHMPIDQualAssDataMaker)
45            
46 //____________________________________________________________________________ 
47   AliHMPIDQualAssDataMaker::AliHMPIDQualAssDataMaker() : 
48   AliQualAssDataMaker(AliQualAss::GetDetName(AliQualAss::kHMPID), "HMPID Quality Assurance Data Maker"),
49   fhHitQdc(0x0), 
50   fhSDigits(0x0),
51   fhDigPcEvt(0x0),
52   fhDigChEvt(0x0),
53   fhDigQ(0x0),
54   fhCluEvt(0x0),
55   fhCluChi2(0x0),
56   fhCluQ(0x0),
57   fhCluFlg(0x0), 
58   fhCluSize(0x0),  
59   fhMipCluSize(0x0),
60   fhCkovP(0x0),
61   fhSigP(0x0),
62   fhMipXY(0x0),
63   fhDifXY(0x0)
64 {
65   // ctor
66   for(Int_t i=0; i<7; i++) fhHitMap[i]=0x0;
67   for(Int_t j=0; j<5; j++) fhPid[j]=0x0;
68 //   fDetectorDir = fOutput->GetDirectory(GetName()) ;  
69 //   if (!fDetectorDir) 
70 //     fDetectorDir = fOutput->mkdir(GetName()) ;  
71 }
72
73 //____________________________________________________________________________ 
74 AliHMPIDQualAssDataMaker::AliHMPIDQualAssDataMaker(const AliHMPIDQualAssDataMaker& qadm) :
75   AliQualAssDataMaker(), 
76   fhHitQdc(qadm.fhHitQdc), 
77   fhSDigits(qadm.fhSDigits),
78   fhDigPcEvt(qadm.fhDigPcEvt),
79   fhDigChEvt(qadm.fhDigChEvt),
80   fhDigQ(qadm.fhDigQ),
81   fhCluEvt(qadm.fhCluEvt),
82   fhCluChi2(qadm.fhCluChi2),
83   fhCluQ(qadm.fhCluQ),
84   fhCluFlg(qadm.fhCluFlg),
85   fhCluSize(qadm.fhCluSize),
86   fhMipCluSize(qadm.fhMipCluSize),
87   fhCkovP(qadm.fhCkovP),
88   fhSigP(qadm.fhSigP),
89   fhMipXY(qadm.fhMipXY),
90   fhDifXY(qadm.fhDifXY)
91 {
92   //copy ctor 
93   for(Int_t i=0; i<7; i++) fhHitMap[i]=qadm.fhHitMap[i];
94   for(Int_t j=0; j<5; j++) fhPid[j]=qadm.fhPid[j];
95
96   SetName((const char*)qadm.GetName()) ; 
97   SetTitle((const char*)qadm.GetTitle()); 
98 }
99
100 //__________________________________________________________________
101 AliHMPIDQualAssDataMaker& AliHMPIDQualAssDataMaker::operator = (const AliHMPIDQualAssDataMaker& qadm )
102 {
103   // Equal operator.
104   this->~AliHMPIDQualAssDataMaker();
105   new(this) AliHMPIDQualAssDataMaker(qadm);
106   return *this;
107 }
108  
109 //____________________________________________________________________________ 
110 void AliHMPIDQualAssDataMaker::InitHits()
111 {
112   // create Hits histograms in Hits subdir
113      fhHitQdc=new TH1F("HitQdc","HMPID Hit Qdc all chamber;QDC",500,0,4000);
114      for(Int_t iCh=0;iCh<7;iCh++) fhHitMap[iCh]=new TH2F(Form("HMPID HitMap%i",iCh),Form("Ch%i;x_{Hit};y_{Hit}",iCh),162,-1,161,146,-1,145);   
115 }
116
117 //____________________________________________________________________________ 
118 void AliHMPIDQualAssDataMaker::InitDigits()
119 {
120   // create Digits histograms in Digits subdir
121       fhDigPcEvt=new TH1F("hDigPcEvt","PC occupancy",156,-1,77);
122       fhDigChEvt=new TH1F("hDigChEvt","Chamber occupancy",32,-1,7);
123       fhDigQ  =new TH1F("Q        "," digit charge               ",3000,0,3000);
124 }
125
126 //____________________________________________________________________________ 
127 void AliHMPIDQualAssDataMaker::InitSDigits()
128 {
129   // create SDigits histograms in SDigits subdir
130       fhSDigits     = new TH1F("hHmpidSDigits",    "SDigits Q  distribution in HMPID",  500, 0., 5000.) ; 
131 }
132
133 //____________________________________________________________________________ 
134
135 void AliHMPIDQualAssDataMaker::InitRecPoints()
136 {
137   // create cluster histograms in RecPoint subdir
138       fhCluEvt=new TH1F("CluPerEvt","# clusters per chamber",16,-1,7);
139       fhCluChi2  =new TH1F("CluChi2"  ,"Chi2 "               ,1000,0,100);
140       fhCluQ   =new TH1F("CluQ"   ,"Cluster charge"        ,3000,0,3000);
141       fhCluFlg   =new TH1F("CluFlg"   ,"Cluster flag"        ,14,-1.5,12.5);
142       fhCluSize  =new TH1F("CluSize"  ,"Raw cluster size    ",100,0,100);
143       fhMipCluSize =new TH1F("MipCluSize"  ,"Mip cluster size    ",100,0,100);
144 }
145 //____________________________________________________________________________
146 void AliHMPIDQualAssDataMaker::InitESDs()
147 {
148   //create ESDs histograms in ESDs subdir
149      fhCkovP  = new TH2F("CkovP" , "#theta_{c}, [rad];P, [GeV]"   , 150,   0,  7  ,100, 0, 1)   ;
150      fhSigP   = new TH2F("SigP"  ,"#sigma_{#theta_c} [mrad];[GeV]", 150,   0,  7  ,100, 0, 1)   ;
151      fhMipXY  = new TH2F("MipXY" ,"mip position"                  , 260,   0,130  ,252, 0,126)  ;
152      fhDifXY  = new TH2F("DifXY" ,"diff"                          , 200, -10, 10  ,200,-10,10)  ;
153      fhPid[0] = new TH1F("PidE" ,"PID: e yellow #mu magenta"  ,100,0,1)                         ;
154      fhPid[1] = new TH1F("PidMu","pid of #mu"                 ,100,0,1)                         ;
155      fhPid[2] = new TH1F("PidPi","PID: #pi red K green p blue",100,0,1)                         ;
156      fhPid[3] = new TH1F("PidK" ,"pid of K"                   ,100,0,1)                         ;
157      fhPid[4] = new TH1F("PidP" ,"pid of p"                   ,100,0,1)                         ;
158 }
159
160 //____________________________________________________________________________
161 void AliHMPIDQualAssDataMaker::MakeHits(TObject * data)
162 {
163   //fills QA histos for Hits
164   TClonesArray * hits = dynamic_cast<TClonesArray *>(data) ; 
165   if (!hits){
166     AliError("Wrong type of hits container") ; 
167   } else {
168     TIter next(hits); 
169     AliHMPIDHit * hit ; 
170     while ( (hit = dynamic_cast<AliHMPIDHit *>(next())) ) {
171       if(hit->Pid()<500000) fhHitQdc->Fill(hit->Q()) ;
172       if(hit->Pid()<500000) fhHitMap[hit->Ch()]->Fill(hit->LorsX(),hit->LorsY());
173     }
174   } 
175 }
176
177 //____________________________________________________________________________
178 void AliHMPIDQualAssDataMaker::MakeDigits( TObject * data)
179 {
180   //fills QA histos for Digits
181   TObjArray *chambers = dynamic_cast<TObjArray*>(data);
182   if ( !chambers) {
183     AliError("Wrong type of digits container") ; 
184   } else {
185     for(Int_t i =0; i< chambers->GetEntries(); i++)
186       {
187         TClonesArray * digits = dynamic_cast<TClonesArray*>(chambers->At(i)); 
188         fhDigChEvt->Fill(i,digits->GetEntriesFast()/(48.*80.*6.));
189         TIter next(digits); 
190         AliHMPIDDigit * digit; 
191         while ( (digit = dynamic_cast<AliHMPIDDigit *>(next())) ) {
192           fhDigPcEvt->Fill(10.*i+digit->Pc(),1./(48.*80.));
193           fhDigQ->Fill(digit->Q());
194         }  
195       }
196   }
197 }
198
199 //____________________________________________________________________________
200 void AliHMPIDQualAssDataMaker::MakeSDigits( TObject * data)
201 {
202   //fills QA histos for SDigits
203   TClonesArray * sdigits = dynamic_cast<TClonesArray *>(data) ; 
204   if (!sdigits) {
205     AliError("Wrong type of sdigits container") ; 
206   } else {
207     AliHMPIDDigit *ref = (AliHMPIDDigit *)sdigits->At(0);
208     Float_t zero = ref->GetTrack(0); 
209     TIter next(sdigits) ; 
210     AliHMPIDDigit * sdigit ; 
211     while ( (sdigit = dynamic_cast<AliHMPIDDigit *>(next())) ) {
212       fhSDigits->Fill(sdigit->Q()) ;
213       if(zero == sdigit->GetTrack(0)) continue;
214       else zero = sdigit->GetTrack(0);
215     } 
216   }
217 }
218
219 //____________________________________________________________________________
220 void AliHMPIDQualAssDataMaker::MakeRecPoints(TTree * clustersTree)
221 {
222   //fills QA histos for clusters
223
224   TClonesArray *clusters = new TClonesArray("AliHMPIDCluster");
225   for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){
226     TBranch *branch = clustersTree->GetBranch(Form("HMPID%d",i));
227     branch->SetAddress(&clusters);
228     branch->GetEntry(0);
229
230     fhCluEvt->Fill(i,clusters->GetEntries());
231     TIter next(clusters);
232     AliHMPIDCluster *clu;
233     while ( (clu = dynamic_cast<AliHMPIDCluster *>(next())) ) {;
234       fhCluFlg->Fill(clu->Status());  fhCluChi2->Fill(clu->Chi2());  fhCluSize->Fill(clu->Size());
235       fhCluQ->Fill(clu->Q()); 
236       Int_t qCut=100;
237       if(clu->Q()>qCut) {
238         fhMipCluSize->SetTitle(Form("Mip cluster size at a Qcut = %i ADC",qCut));
239         fhMipCluSize->Fill(clu->Size());
240       }
241     }
242   }
243
244   clusters->Delete();
245   delete clusters;
246 }
247
248 //____________________________________________________________________________
249 void AliHMPIDQualAssDataMaker::MakeESDs(AliESDEvent * esd)
250 {
251   //fills QA histos for ESD
252   for(Int_t iTrk = 0 ; iTrk < esd->GetNumberOfTracks() ; iTrk++){
253     AliESDtrack *pTrk = esd->GetTrack(iTrk) ;
254     fhCkovP->Fill(pTrk->GetP(),pTrk->GetHMPIDsignal());
255     fhSigP->Fill( pTrk->GetP(),TMath::Sqrt(pTrk->GetHMPIDchi2()));
256     Float_t xm,ym; Int_t q,np;  
257     pTrk->GetHMPIDmip(xm,ym,q,np);                       //mip info
258     fhMipXY->Fill(xm,ym);
259     Float_t xRad,yRad,th,ph;        
260     pTrk->GetHMPIDtrk(xRad,yRad,th,ph);              //track info at the middle of the radiator
261     Float_t xPc = xRad+9.25*TMath::Tan(th)*TMath::Cos(ph); // temporar: linear extrapol (B=0!)
262     Float_t yPc = yRad+9.25*TMath::Tan(th)*TMath::Sin(ph); // temporar:          "
263     fhDifXY->Fill(xm-xPc,ym-yPc); //track info
264     Double_t pid[5] ;      pTrk->GetHMPIDpid(pid) ;
265     for(Int_t i = 0 ; i < 5 ; i++) fhPid[i]->Fill(pid[i]) ;
266   }
267 }
268