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