]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCQADataMakerSim.cxx
Fix for the problem during PbPb run of Nov 2010 (Indra)
[u/mrichter/AliRoot.git] / TPC / AliTPCQADataMakerSim.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, 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   Based on AliPHOSQADataMaker
21   Produces the data needed to calculate the quality assurance. 
22   All data must be mergeable objects.
23   P. Christiansen, Lund, January 2008
24 */
25
26 /*
27   Implementation:
28
29   We have chosen to have the histograms as non-persistent meber to
30   allow better debugging. In the copy constructor we then have to
31   assign the pointers to the existing histograms in the copied
32   list. This have been implemented but not tested.
33 */
34
35 #include "AliTPCQADataMakerSim.h"
36
37 // --- ROOT system ---
38
39 // --- Standard library ---
40
41 // --- AliRoot header files ---
42 #include "AliQAChecker.h"
43 #include "AliTPC.h"
44 #include "AliTPCv2.h"
45 #include "AliSimDigits.h"
46
47 ClassImp(AliTPCQADataMakerSim)
48
49 //____________________________________________________________________________ 
50 AliTPCQADataMakerSim::AliTPCQADataMakerSim() : 
51   AliQADataMakerSim(AliQAv1::GetDetName(AliQAv1::kTPC), 
52                     "TPC Sim Quality Assurance Data Maker")
53 {
54   // ctor
55 }
56
57 //____________________________________________________________________________ 
58 AliTPCQADataMakerSim::AliTPCQADataMakerSim(const AliTPCQADataMakerSim& qadm) :
59   AliQADataMakerSim()
60 {
61   //copy ctor 
62   SetName((const char*)qadm.GetName()) ; 
63   SetTitle((const char*)qadm.GetTitle()); 
64   
65   //
66   // Associate class histogram objects to the copies in the list
67   // Could also be done with the indexes
68   //
69  }
70
71 //__________________________________________________________________
72 AliTPCQADataMakerSim& AliTPCQADataMakerSim::operator = (const AliTPCQADataMakerSim& qadm )
73 {
74   // Equal operator.
75   this->~AliTPCQADataMakerSim();
76   new(this) AliTPCQADataMakerSim(qadm);
77   return *this;
78 }
79  
80 //____________________________________________________________________________ 
81 void AliTPCQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
82 {
83   //Detector specific actions at end of cycle
84   // do the QA checking
85   AliQAChecker::Instance()->Run(AliQAv1::kTPC, task, list) ;  
86 }
87
88 //____________________________________________________________________________ 
89 void AliTPCQADataMakerSim::InitDigits()
90 {
91   const Bool_t expert   = kTRUE ; 
92   const Bool_t image    = kTRUE ; 
93   TH1F * histDigitsADC = 
94     new TH1F("hDigitsADC", "Digit ADC distribution; ADC; Counts",
95              1000, 0, 1000);
96   histDigitsADC->Sumw2();
97   Add2DigitsList(histDigitsADC, kDigitsADC, !expert, image);
98 }
99
100 //____________________________________________________________________________ 
101 void AliTPCQADataMakerSim::InitHits()
102 {
103   const Bool_t expert   = kTRUE ; 
104   const Bool_t image    = kTRUE ; 
105   TH1F * histHitsNhits = 
106     new TH1F("hHitsNhits", "Interactions per track in the TPC volume; Number of interactions; Counts",
107              100, 0, 10000);
108   histHitsNhits->Sumw2();
109   Add2HitsList(histHitsNhits, kNhits, !expert, image);
110
111   TH1F * histHitsElectrons = 
112     new TH1F("hHitsElectrons", "Electrons per interaction; Electrons; Counts",
113              300, 0, 300);
114   histHitsElectrons->Sumw2();
115   Add2HitsList(histHitsElectrons, kElectrons, !expert, image);  
116
117   TH1F * histHitsRadius = 
118     new TH1F("hHitsRadius", "Position of interaction; Radius; Counts",
119              300, 0., 300.);  
120   histHitsRadius->Sumw2();
121   Add2HitsList(histHitsRadius, kRadius, !expert, image);  
122
123   TH1F * histHitsPrimPerCm = 
124     new TH1F("hHitsPrimPerCm", "Primaries per cm; Primaries; Counts",
125              100, 0., 100.);  
126   histHitsPrimPerCm->Sumw2();
127   Add2HitsList(histHitsPrimPerCm, kPrimPerCm, !expert, image);  
128
129   TH1F * histHitsElectronsPerCm = 
130     new TH1F("hHitsElectronsPerCm", "Electrons per cm; Electrons; Counts",
131              300, 0., 300.);  
132   histHitsElectronsPerCm->Sumw2();
133   Add2HitsList(histHitsElectronsPerCm, kElectronsPerCm, !expert, image);  
134 }
135
136 //____________________________________________________________________________
137 void AliTPCQADataMakerSim::MakeDigits(TTree* digitTree)
138 {
139
140   TBranch* branch = digitTree->GetBranch("Segment");
141   AliSimDigits* digArray = 0;
142   branch->SetAddress(&digArray);
143   
144   Int_t nEntries = Int_t(digitTree->GetEntries());
145   
146   for (Int_t n = 0; n < nEntries; n++) {
147     
148     digitTree->GetEvent(n);
149     
150     if (digArray->First())
151       do {
152         Float_t dig = digArray->CurrentDigit();
153         
154         GetDigitsData(kDigitsADC)->Fill(dig);
155       } while (digArray->Next());    
156   }
157 }
158
159 //____________________________________________________________________________
160 void AliTPCQADataMakerSim::MakeHits(TTree * hitTree)
161 {
162   // make QA data from Hit Tree
163  
164   const Int_t nTracks = hitTree->GetEntries();
165   TBranch* branch = hitTree->GetBranch("TPC2");
166   AliTPCv2* tpc = (AliTPCv2*)gAlice->GetDetector("TPC");
167   
168   //
169   // loop over tracks
170   //
171   for(Int_t n = 0; n < nTracks; n++){
172     Int_t nHits = 0;
173     branch->GetEntry(n);
174     
175     AliTPChit* tpcHit = (AliTPChit*)tpc->FirstHit(-1);  
176     
177     if (tpcHit) {
178       Float_t dist  = 0.;
179       Int_t   nprim = 0;
180       Float_t xold  = tpcHit->X();
181       Float_t yold  = tpcHit->Y();
182       Float_t zold  = tpcHit->Z(); 
183       Float_t radiusOld = TMath::Sqrt(xold*xold + yold*yold); 
184       Int_t trackOld = tpcHit->GetTrack();
185       Float_t q     = 0.;
186
187       while(tpcHit) {
188
189         Float_t x = tpcHit->X();
190         Float_t y = tpcHit->Y();
191         Float_t z = tpcHit->Z(); 
192         Float_t radius = TMath::Sqrt(x*x + y*y);
193         
194         if(radius>50) { // Skip hits at interaction point
195           
196           nHits++;
197           
198           Int_t trackNo = tpcHit->GetTrack();
199           
200           GetHitsData(kElectrons)->Fill(tpcHit->fQ);
201           GetHitsData(kRadius)->Fill(radius);
202             
203           if(trackNo==trackOld) { // if the same track
204
205             // find the new distance
206             dist += TMath::Sqrt((x-xold)*(x-xold) + (y-yold)*(y-yold) + 
207                                 (z-zold)*(z-zold));
208             if(dist<1.){ // add data to this 1 cm step
209               
210               nprim++;  
211               q += tpcHit->fQ;        
212             } else{ // Fill the histograms normalized to per cm 
213               
214               if(nprim==1)
215                 cout << radius << ", " << radiusOld << ", " << dist << endl; 
216               
217               GetHitsData(kPrimPerCm)->Fill((Float_t)nprim);
218               GetHitsData(kElectronsPerCm)->Fill(q);
219               
220               dist  = 0;
221               q     = 0;
222               nprim = 0;
223             }
224           } else { // reset for new track
225             
226             dist  = 0;
227             q     = 0;
228             nprim = 0;
229           }
230         }
231
232         radiusOld = radius;
233         xold = x;
234         yold = y;
235         zold = z;
236         trackOld = tpcHit->GetTrack();
237         
238         tpcHit = (AliTPChit*) tpc->NextHit();
239       }
240     }
241
242     GetHitsData(kNhits)->Fill(nHits);
243   }
244 }
245