1 /**************************************************************************
2 * Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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
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.
36 #include "AliTPCQADataMakerSim.h"
38 // --- ROOT system ---
40 // --- Standard library ---
42 // --- AliRoot header files ---
43 #include "AliQAChecker.h"
46 #include "AliSimDigits.h"
48 ClassImp(AliTPCQADataMakerSim)
50 //____________________________________________________________________________
51 AliTPCQADataMakerSim::AliTPCQADataMakerSim() :
52 AliQADataMakerSim(AliQAv1::GetDetName(AliQAv1::kTPC),
53 "TPC Sim Quality Assurance Data Maker")
58 //____________________________________________________________________________
59 AliTPCQADataMakerSim::AliTPCQADataMakerSim(const AliTPCQADataMakerSim& qadm) :
63 SetName((const char*)qadm.GetName()) ;
64 SetTitle((const char*)qadm.GetTitle());
67 // Associate class histogram objects to the copies in the list
68 // Could also be done with the indexes
72 //__________________________________________________________________
73 AliTPCQADataMakerSim& AliTPCQADataMakerSim::operator = (const AliTPCQADataMakerSim& qadm )
76 this->~AliTPCQADataMakerSim();
77 new(this) AliTPCQADataMakerSim(qadm);
81 //____________________________________________________________________________
82 void AliTPCQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
84 //Detector specific actions at end of cycle
86 AliQAChecker::Instance()->Run(AliQAv1::kTPC, task, list) ;
89 //____________________________________________________________________________
90 void AliTPCQADataMakerSim::InitDigits()
92 const Bool_t expert = kTRUE ;
93 const Bool_t image = kTRUE ;
94 TH1F * histDigitsADC =
95 new TH1F("hDigitsADC", "Digit ADC distribution; ADC; Counts",
97 histDigitsADC->Sumw2();
98 Add2DigitsList(histDigitsADC, kDigitsADC, !expert, image);
101 //____________________________________________________________________________
102 void AliTPCQADataMakerSim::InitHits()
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",
109 histHitsNhits->Sumw2();
110 Add2HitsList(histHitsNhits, kNhits, !expert, image);
112 TH1F * histHitsElectrons =
113 new TH1F("hHitsElectrons", "Electrons per interaction; Electrons; Counts",
115 histHitsElectrons->Sumw2();
116 Add2HitsList(histHitsElectrons, kElectrons, !expert, image);
118 TH1F * histHitsRadius =
119 new TH1F("hHitsRadius", "Position of interaction; Radius; Counts",
121 histHitsRadius->Sumw2();
122 Add2HitsList(histHitsRadius, kRadius, !expert, image);
124 TH1F * histHitsPrimPerCm =
125 new TH1F("hHitsPrimPerCm", "Primaries per cm; Primaries; Counts",
127 histHitsPrimPerCm->Sumw2();
128 Add2HitsList(histHitsPrimPerCm, kPrimPerCm, !expert, image);
130 TH1F * histHitsElectronsPerCm =
131 new TH1F("hHitsElectronsPerCm", "Electrons per cm; Electrons; Counts",
133 histHitsElectronsPerCm->Sumw2();
134 Add2HitsList(histHitsElectronsPerCm, kElectronsPerCm, !expert, image);
137 //____________________________________________________________________________
138 void AliTPCQADataMakerSim::MakeDigits(TTree* digitTree)
141 TBranch* branch = digitTree->GetBranch("Segment");
142 AliSimDigits* digArray = 0;
143 branch->SetAddress(&digArray);
145 Int_t nEntries = Int_t(digitTree->GetEntries());
147 for (Int_t n = 0; n < nEntries; n++) {
149 digitTree->GetEvent(n);
151 if (digArray->First())
153 Float_t dig = digArray->CurrentDigit();
155 GetDigitsData(kDigitsADC)->Fill(dig);
156 } while (digArray->Next());
160 //____________________________________________________________________________
161 void AliTPCQADataMakerSim::MakeHits(TTree * hitTree)
163 // make QA data from Hit Tree
165 const Int_t nTracks = hitTree->GetEntries();
166 TBranch* branch = hitTree->GetBranch("TPC2");
167 AliTPCv2* tpc = (AliTPCv2*)gAlice->GetDetector("TPC");
172 for(Int_t n = 0; n < nTracks; n++){
176 AliTPChit* tpcHit = (AliTPChit*)tpc->FirstHit(-1);
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();
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);
195 if(radius>50) { // Skip hits at interaction point
199 Int_t trackNo = tpcHit->GetTrack();
201 GetHitsData(kElectrons)->Fill(tpcHit->fQ);
202 GetHitsData(kRadius)->Fill(radius);
204 if(trackNo==trackOld) { // if the same track
206 // find the new distance
207 dist += TMath::Sqrt((x-xold)*(x-xold) + (y-yold)*(y-yold) +
209 if(dist<1.){ // add data to this 1 cm step
213 } else{ // Fill the histograms normalized to per cm
216 // cout << radius << ", " << radiusOld << ", " << dist << endl;
218 GetHitsData(kPrimPerCm)->Fill((Float_t)nprim);
219 GetHitsData(kElectronsPerCm)->Fill(q);
225 } else { // reset for new track
237 trackOld = tpcHit->GetTrack();
239 tpcHit = (AliTPChit*) tpc->NextHit();
243 GetHitsData(kNhits)->Fill(nHits);