]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCQADataMakerSim.cxx
Fix for #83391 printings from AliTPCQADataMakerSim polluting sim.log
[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 <TTree.h>
36 #include "AliTPCQADataMakerSim.h"
37
38 // --- ROOT system ---
39
40 // --- Standard library ---
41
42 // --- AliRoot header files ---
43 #include "AliQAChecker.h"
44 #include "AliTPC.h"
45 #include "AliTPCv2.h"
46 #include "AliSimDigits.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
101 //____________________________________________________________________________ 
102 void AliTPCQADataMakerSim::InitHits()
103 {
104   const Bool_t expert   = kTRUE ; 
105   const Bool_t image    = kTRUE ; 
106   TH1F * histHitsNhits = 
107     new TH1F("hHitsNhits", "Interactions per track in the TPC volume; Number of interactions; Counts",
108              100, 0, 10000);
109   histHitsNhits->Sumw2();
110   Add2HitsList(histHitsNhits, kNhits, !expert, image);
111
112   TH1F * histHitsElectrons = 
113     new TH1F("hHitsElectrons", "Electrons per interaction; Electrons; Counts",
114              300, 0, 300);
115   histHitsElectrons->Sumw2();
116   Add2HitsList(histHitsElectrons, kElectrons, !expert, image);  
117
118   TH1F * histHitsRadius = 
119     new TH1F("hHitsRadius", "Position of interaction; Radius; Counts",
120              300, 0., 300.);  
121   histHitsRadius->Sumw2();
122   Add2HitsList(histHitsRadius, kRadius, !expert, image);  
123
124   TH1F * histHitsPrimPerCm = 
125     new TH1F("hHitsPrimPerCm", "Primaries per cm; Primaries; Counts",
126              100, 0., 100.);  
127   histHitsPrimPerCm->Sumw2();
128   Add2HitsList(histHitsPrimPerCm, kPrimPerCm, !expert, image);  
129
130   TH1F * histHitsElectronsPerCm = 
131     new TH1F("hHitsElectronsPerCm", "Electrons per cm; Electrons; Counts",
132              300, 0., 300.);  
133   histHitsElectronsPerCm->Sumw2();
134   Add2HitsList(histHitsElectronsPerCm, kElectronsPerCm, !expert, image);  
135 }
136
137 //____________________________________________________________________________
138 void AliTPCQADataMakerSim::MakeDigits(TTree* digitTree)
139 {
140
141   TBranch* branch = digitTree->GetBranch("Segment");
142   AliSimDigits* digArray = 0;
143   branch->SetAddress(&digArray);
144   
145   Int_t nEntries = Int_t(digitTree->GetEntries());
146   
147   for (Int_t n = 0; n < nEntries; n++) {
148     
149     digitTree->GetEvent(n);
150     
151     if (digArray->First())
152       do {
153         Float_t dig = digArray->CurrentDigit();
154         
155         GetDigitsData(kDigitsADC)->Fill(dig);
156       } while (digArray->Next());    
157   }
158 }
159
160 //____________________________________________________________________________
161 void AliTPCQADataMakerSim::MakeHits(TTree * hitTree)
162 {
163   // make QA data from Hit Tree
164  
165   const Int_t nTracks = hitTree->GetEntries();
166   TBranch* branch = hitTree->GetBranch("TPC2");
167   AliTPCv2* tpc = (AliTPCv2*)gAlice->GetDetector("TPC");
168   
169   //
170   // loop over tracks
171   //
172   for(Int_t n = 0; n < nTracks; n++){
173     Int_t nHits = 0;
174     branch->GetEntry(n);
175     
176     AliTPChit* tpcHit = (AliTPChit*)tpc->FirstHit(-1);  
177     
178     if (tpcHit) {
179       Float_t dist  = 0.;
180       Int_t   nprim = 0;
181       Float_t xold  = tpcHit->X();
182       Float_t yold  = tpcHit->Y();
183       Float_t zold  = tpcHit->Z(); 
184       Float_t radiusOld = TMath::Sqrt(xold*xold + yold*yold); 
185       Int_t trackOld = tpcHit->GetTrack();
186       Float_t q     = 0.;
187
188       while(tpcHit) {
189
190         Float_t x = tpcHit->X();
191         Float_t y = tpcHit->Y();
192         Float_t z = tpcHit->Z(); 
193         Float_t radius = TMath::Sqrt(x*x + y*y);
194         
195         if(radius>50) { // Skip hits at interaction point
196           
197           nHits++;
198           
199           Int_t trackNo = tpcHit->GetTrack();
200           
201           GetHitsData(kElectrons)->Fill(tpcHit->fQ);
202           GetHitsData(kRadius)->Fill(radius);
203             
204           if(trackNo==trackOld) { // if the same track
205
206             // find the new distance
207             dist += TMath::Sqrt((x-xold)*(x-xold) + (y-yold)*(y-yold) + 
208                                 (z-zold)*(z-zold));
209             if(dist<1.){ // add data to this 1 cm step
210               
211               nprim++;  
212               q += tpcHit->fQ;        
213             } else{ // Fill the histograms normalized to per cm 
214               
215               // if(nprim==1)
216               //        cout << radius << ", " << radiusOld << ", " << dist << endl; 
217               
218               GetHitsData(kPrimPerCm)->Fill((Float_t)nprim);
219               GetHitsData(kElectronsPerCm)->Fill(q);
220               
221               dist  = 0;
222               q     = 0;
223               nprim = 0;
224             }
225           } else { // reset for new track
226             
227             dist  = 0;
228             q     = 0;
229             nprim = 0;
230           }
231         }
232
233         radiusOld = radius;
234         xold = x;
235         yold = y;
236         zold = z;
237         trackOld = tpcHit->GetTrack();
238         
239         tpcHit = (AliTPChit*) tpc->NextHit();
240       }
241     }
242
243     GetHitsData(kNhits)->Fill(nHits);
244   }
245 }
246