]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCQADataMakerSim.cxx
Removing obsolete mapping macro
[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(AliQA::GetDetName(AliQA::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(AliQA::TASKINDEX_t task, TObjArray ** list)
82 {
83   //Detector specific actions at end of cycle
84   // do the QA checking
85   AliQAChecker::Instance()->Run(AliQA::kTPC, task, list) ;  
86 }
87
88 //____________________________________________________________________________ 
89 void AliTPCQADataMakerSim::InitDigits()
90 {
91   TH1F * histDigitsADC = 
92     new TH1F("hDigitsADC", "Digit ADC distribution; ADC; Counts",
93              1000, 0, 1000);
94   histDigitsADC->Sumw2();
95   Add2DigitsList(histDigitsADC, kDigitsADC);
96 }
97
98 //____________________________________________________________________________ 
99 void AliTPCQADataMakerSim::InitHits()
100 {
101   TH1F * histHitsNhits = 
102     new TH1F("hHitsNhits", "Interactions per primary track in the TPC volume; Number of primary interactions; Counts",
103              100, 0, 10000);
104   histHitsNhits->Sumw2();
105   Add2HitsList(histHitsNhits, kNhits);
106
107   TH1F * histHitsElectrons = 
108     new TH1F("hHitsElectrons", "Electrons per interaction (primaries); Electrons; Counts",
109              300, 0, 300);
110   histHitsElectrons->Sumw2();
111   Add2HitsList(histHitsElectrons, kElectrons);  
112
113   TH1F * histHitsRadius = 
114     new TH1F("hHitsRadius", "Position of interaction (Primary tracks only); Radius; Counts",
115              300, 0., 300.);  
116   histHitsRadius->Sumw2();
117   Add2HitsList(histHitsRadius, kRadius);  
118
119   TH1F * histHitsPrimPerCm = 
120     new TH1F("hHitsPrimPerCm", "Primaries per cm (Primary tracks only); Primaries; Counts",
121              100, 0., 100.);  
122   histHitsPrimPerCm->Sumw2();
123   Add2HitsList(histHitsPrimPerCm, kPrimPerCm);  
124
125   TH1F * histHitsElectronsPerCm = 
126     new TH1F("hHitsElectronsPerCm", "Electrons per cm (Primary tracks only); Electrons; Counts",
127              300, 0., 300.);  
128   histHitsElectronsPerCm->Sumw2();
129   Add2HitsList(histHitsElectronsPerCm, kElectronsPerCm);  
130 }
131
132 //____________________________________________________________________________
133 void AliTPCQADataMakerSim::MakeDigits(TTree* digitTree)
134 {
135   TBranch* branch = digitTree->GetBranch("Segment");
136   AliSimDigits* digArray = 0;
137   branch->SetAddress(&digArray);
138   
139   Int_t nEntries = Int_t(digitTree->GetEntries());
140   
141   for (Int_t n = 0; n < nEntries; n++) {
142     
143     digitTree->GetEvent(n);
144     
145     if (digArray->First())
146       do {
147         Float_t dig = digArray->CurrentDigit();
148         
149         GetDigitsData(kDigitsADC)->Fill(dig);
150       } while (digArray->Next());    
151   }
152 }
153
154 //____________________________________________________________________________
155 void AliTPCQADataMakerSim::MakeHits(TTree * hitTree)
156 {
157   // make QA data from Hit Tree
158   const Int_t nTracks = hitTree->GetEntries();
159   TBranch* branch = hitTree->GetBranch("TPC2");
160   AliTPCv2* tpc = (AliTPCv2*)gAlice->GetDetector("TPC");
161   
162   //
163   // loop over tracks
164   //
165   for(Int_t n = 0; n < nTracks; n++){
166     Int_t nHits = 0;
167     branch->GetEntry(n);
168     
169     AliTPChit* tpcHit = (AliTPChit*)tpc->FirstHit(-1);  
170     
171     if (tpcHit) {
172       Float_t dist  = 0.;
173       Int_t   nprim = 0;
174       Float_t xold  = tpcHit->X();
175       Float_t yold  = tpcHit->Y();
176       Float_t zold  = tpcHit->Z(); 
177       Float_t radiusOld = TMath::Sqrt(xold*xold + yold*yold);
178       Float_t q     = 0.;
179       while(tpcHit) {
180         Float_t x = tpcHit->X();
181         Float_t y = tpcHit->Y();
182         Float_t z = tpcHit->Z(); 
183         Float_t radius = TMath::Sqrt(x*x + y*y);
184         
185         if(radius>50) { // Skip hits at interaction point
186           
187           nHits++;
188           
189           Int_t trackNo = tpcHit->GetTrack();
190           
191           if(trackNo==n) { // primary track
192             
193             GetHitsData(kElectrons)->Fill(tpcHit->fQ);
194             GetHitsData(kRadius)->Fill(radius);
195             
196             // find the new distance
197             dist += TMath::Sqrt((x-xold)*(x-xold) + (y-yold)*(y-yold) + 
198                                 (z-zold)*(z-zold));
199             if(dist<1.){  
200               // add data to this 1 cm step
201               nprim++;  
202               q += tpcHit->fQ;
203               
204             } else{
205               // Fill the histograms normalized to per cm 
206               
207               if(nprim==1)
208                 cout << radius << ", " << radiusOld << ", " << dist << endl; 
209               
210               GetHitsData(kPrimPerCm)->Fill((Float_t)nprim);
211               GetHitsData(kElectronsPerCm)->Fill(q);
212               
213               dist  = 0;
214               q     = 0;
215               nprim = 0;
216             }
217           }
218         }
219         
220         radiusOld = radius;
221         xold = x;
222         yold = y;
223         zold = z;
224         
225         tpcHit = (AliTPChit*) tpc->NextHit();
226         }
227       }
228       GetHitsData(kNhits)->Fill(nHits);
229     }
230   }
231