]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCQADataMakerSim.cxx
Adding macro to create Pad response function for the GEM readout
[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 #include <TTree.h>
47
48 ClassImp(AliTPCQADataMakerSim)
49
50 //____________________________________________________________________________ 
51 AliTPCQADataMakerSim::AliTPCQADataMakerSim() : 
52   AliQADataMakerSim(AliQAv1::GetDetName(AliQAv1::kTPC), 
53                     "TPC Sim Quality Assurance Data Maker")
54 {
55   // ctor
56 }
57
58 //____________________________________________________________________________ 
59 AliTPCQADataMakerSim::AliTPCQADataMakerSim(const AliTPCQADataMakerSim& qadm) :
60   AliQADataMakerSim()
61 {
62   //copy ctor 
63   SetName((const char*)qadm.GetName()) ; 
64   SetTitle((const char*)qadm.GetTitle()); 
65   
66   //
67   // Associate class histogram objects to the copies in the list
68   // Could also be done with the indexes
69   //
70  }
71
72 //__________________________________________________________________
73 AliTPCQADataMakerSim& AliTPCQADataMakerSim::operator = (const AliTPCQADataMakerSim& qadm )
74 {
75   // Equal operator.
76   this->~AliTPCQADataMakerSim();
77   new(this) AliTPCQADataMakerSim(qadm);
78   return *this;
79 }
80  
81 //____________________________________________________________________________ 
82 void AliTPCQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
83 {
84   //Detector specific actions at end of cycle
85   // do the QA checking
86   AliQAChecker::Instance()->Run(AliQAv1::kTPC, task, list) ;  
87 }
88
89 //____________________________________________________________________________ 
90 void AliTPCQADataMakerSim::InitDigits()
91 {
92   const Bool_t expert   = kTRUE ; 
93   const Bool_t image    = kTRUE ; 
94   TH1F * histDigitsADC = 
95     new TH1F("hDigitsADC", "Digit ADC distribution; ADC; Counts",
96              1000, 0, 1000);
97   histDigitsADC->Sumw2();
98   Add2DigitsList(histDigitsADC, kDigitsADC, !expert, image);
99   //
100   ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
101 }
102
103 //____________________________________________________________________________ 
104 void AliTPCQADataMakerSim::InitHits()
105 {
106   const Bool_t expert   = kTRUE ; 
107   const Bool_t image    = kTRUE ; 
108   TH1F * histHitsNhits = 
109     new TH1F("hHitsNhits", "Interactions per track in the TPC volume; Number of interactions; Counts",
110              100, 0, 10000);
111   histHitsNhits->Sumw2();
112   Add2HitsList(histHitsNhits, kNhits, !expert, image);
113
114   TH1F * histHitsElectrons = 
115     new TH1F("hHitsElectrons", "Electrons per interaction; Electrons; Counts",
116              300, 0, 300);
117   histHitsElectrons->Sumw2();
118   Add2HitsList(histHitsElectrons, kElectrons, !expert, image);  
119
120   TH1F * histHitsRadius = 
121     new TH1F("hHitsRadius", "Position of interaction; Radius; Counts",
122              300, 0., 300.);  
123   histHitsRadius->Sumw2();
124   Add2HitsList(histHitsRadius, kRadius, !expert, image);  
125
126   TH1F * histHitsPrimPerCm = 
127     new TH1F("hHitsPrimPerCm", "Primaries per cm; Primaries; Counts",
128              100, 0., 100.);  
129   histHitsPrimPerCm->Sumw2();
130   Add2HitsList(histHitsPrimPerCm, kPrimPerCm, !expert, image);  
131
132   TH1F * histHitsElectronsPerCm = 
133     new TH1F("hHitsElectronsPerCm", "Electrons per cm; Electrons; Counts",
134              300, 0., 300.);  
135   histHitsElectronsPerCm->Sumw2();
136   Add2HitsList(histHitsElectronsPerCm, kElectronsPerCm, !expert, image);  
137   //
138   ClonePerTrigClass(AliQAv1::kHITS); // this should be the last line
139 }
140
141 //____________________________________________________________________________
142 void AliTPCQADataMakerSim::MakeDigits(TTree* digitTree)
143 {
144
145   TBranch* branch = digitTree->GetBranch("Segment");
146   AliSimDigits* digArray = 0;
147   branch->SetAddress(&digArray);
148   
149   Int_t nEntries = Int_t(digitTree->GetEntries());
150   
151   for (Int_t n = 0; n < nEntries; n++) {
152     
153     digitTree->GetEvent(n);
154     
155     if (digArray->First())
156       do {
157         Float_t dig = digArray->CurrentDigit();
158         
159         FillDigitsData(kDigitsADC,dig);
160       } while (digArray->Next());    
161   }
162   //
163   IncEvCountCycleDigits();
164   IncEvCountTotalDigits();  
165   //
166 }
167
168 //____________________________________________________________________________
169 void AliTPCQADataMakerSim::MakeHits(TTree * hitTree)
170 {
171   // make QA data from Hit Tree
172  
173   const Int_t nTracks = hitTree->GetEntries();
174   TBranch* branch = hitTree->GetBranch("TPC2");
175   AliTPCv2* tpc = (AliTPCv2*)gAlice->GetDetector("TPC");
176   
177   //
178   // loop over tracks
179   //
180   for(Int_t n = 0; n < nTracks; n++){
181     Int_t nHits = 0;
182     branch->GetEntry(n);
183     
184     AliTPChit* tpcHit = (AliTPChit*)tpc->FirstHit(-1);  
185     
186     if (tpcHit) {
187       Float_t dist  = 0.;
188       Int_t   nprim = 0;
189       Float_t xold  = tpcHit->X();
190       Float_t yold  = tpcHit->Y();
191       Float_t zold  = tpcHit->Z(); 
192       Float_t radiusOld = TMath::Sqrt(xold*xold + yold*yold); 
193       Int_t trackOld = tpcHit->GetTrack();
194       Float_t q     = 0.;
195
196       while(tpcHit) {
197
198         Float_t x = tpcHit->X();
199         Float_t y = tpcHit->Y();
200         Float_t z = tpcHit->Z(); 
201         Float_t radius = TMath::Sqrt(x*x + y*y);
202         
203         if(radius>50) { // Skip hits at interaction point
204           
205           nHits++;
206           
207           Int_t trackNo = tpcHit->GetTrack();
208           
209           FillHitsData(kElectrons,tpcHit->fQ);
210           FillHitsData(kRadius,radius);
211             
212           if(trackNo==trackOld) { // if the same track
213
214             // find the new distance
215             dist += TMath::Sqrt((x-xold)*(x-xold) + (y-yold)*(y-yold) + 
216                                 (z-zold)*(z-zold));
217             if(dist<1.){ // add data to this 1 cm step
218               
219               nprim++;  
220               q += tpcHit->fQ;        
221             } else{ // Fill the histograms normalized to per cm 
222               
223               // if(nprim==1)
224               //        cout << radius << ", " << radiusOld << ", " << dist << endl; 
225               
226               FillHitsData(kPrimPerCm,(Float_t)nprim);
227               FillHitsData(kElectronsPerCm,q);
228               
229               dist  = 0;
230               q     = 0;
231               nprim = 0;
232             }
233           } else { // reset for new track
234             
235             dist  = 0;
236             q     = 0;
237             nprim = 0;
238           }
239         }
240
241         radiusOld = radius;
242         xold = x;
243         yold = y;
244         zold = z;
245         trackOld = tpcHit->GetTrack();
246         
247         tpcHit = (AliTPChit*) tpc->NextHit();
248       }
249     }
250
251     FillHitsData(kNhits,nHits);
252   }
253   //
254   IncEvCountCycleHits();
255   IncEvCountTotalHits();
256   //
257 }
258