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