]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCQADataMakerSim.cxx
Disable Evchar if no TR read.
[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 primary track in the TPC volume; Number of primary 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 (primaries); 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 (Primary tracks only); 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 (Primary tracks only); 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 (Primary tracks only); 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   // Check id histograms already created for this Event Specie
140   if ( ! GetDigitsData(kDigitsADC) )
141     InitDigits() ;
142
143   TBranch* branch = digitTree->GetBranch("Segment");
144   AliSimDigits* digArray = 0;
145   branch->SetAddress(&digArray);
146   
147   Int_t nEntries = Int_t(digitTree->GetEntries());
148   
149   for (Int_t n = 0; n < nEntries; n++) {
150     
151     digitTree->GetEvent(n);
152     
153     if (digArray->First())
154       do {
155         Float_t dig = digArray->CurrentDigit();
156         
157         GetDigitsData(kDigitsADC)->Fill(dig);
158       } while (digArray->Next());    
159   }
160 }
161
162 //____________________________________________________________________________
163 void AliTPCQADataMakerSim::MakeHits(TTree * hitTree)
164 {
165   // make QA data from Hit Tree
166  
167   // Check id histograms already created for this Event Specie
168   if ( ! GetHitsData(kNhits) )
169     InitHits() ;
170
171   const Int_t nTracks = hitTree->GetEntries();
172   TBranch* branch = hitTree->GetBranch("TPC2");
173   AliTPCv2* tpc = (AliTPCv2*)gAlice->GetDetector("TPC");
174   
175   //
176   // loop over tracks
177   //
178   for(Int_t n = 0; n < nTracks; n++){
179     Int_t nHits = 0;
180     branch->GetEntry(n);
181     
182     AliTPChit* tpcHit = (AliTPChit*)tpc->FirstHit(-1);  
183     
184     if (tpcHit) {
185       Float_t dist  = 0.;
186       Int_t   nprim = 0;
187       Float_t xold  = tpcHit->X();
188       Float_t yold  = tpcHit->Y();
189       Float_t zold  = tpcHit->Z(); 
190       Float_t radiusOld = TMath::Sqrt(xold*xold + yold*yold);
191       Float_t q     = 0.;
192       while(tpcHit) {
193         Float_t x = tpcHit->X();
194         Float_t y = tpcHit->Y();
195         Float_t z = tpcHit->Z(); 
196         Float_t radius = TMath::Sqrt(x*x + y*y);
197         
198         if(radius>50) { // Skip hits at interaction point
199           
200           nHits++;
201           
202           Int_t trackNo = tpcHit->GetTrack();
203           
204           if(trackNo==n) { // primary track
205             
206             GetHitsData(kElectrons)->Fill(tpcHit->fQ);
207             GetHitsData(kRadius)->Fill(radius);
208             
209             // find the new distance
210             dist += TMath::Sqrt((x-xold)*(x-xold) + (y-yold)*(y-yold) + 
211                                 (z-zold)*(z-zold));
212             if(dist<1.){  
213               // add data to this 1 cm step
214               nprim++;  
215               q += tpcHit->fQ;
216               
217             } else{
218               // Fill the histograms normalized to per cm 
219               
220               if(nprim==1)
221                 cout << radius << ", " << radiusOld << ", " << dist << endl; 
222               
223               GetHitsData(kPrimPerCm)->Fill((Float_t)nprim);
224               GetHitsData(kElectronsPerCm)->Fill(q);
225               
226               dist  = 0;
227               q     = 0;
228               nprim = 0;
229             }
230           }
231         }
232         
233         radiusOld = radius;
234         xold = x;
235         yold = y;
236         zold = z;
237         
238         tpcHit = (AliTPChit*) tpc->NextHit();
239         }
240       }
241       GetHitsData(kNhits)->Fill(nHits);
242     }
243   }
244