]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDQADataMakerSim.cxx
Code revisited after debugging. Optimization for functions in V3. Calculations of...
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDQADataMakerSim.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 // --- ROOT system ---
20 #include <TClonesArray.h>
21 #include <TFile.h> 
22 #include <TH1F.h> 
23 #include <TH2F.h>
24 #include <Riostream.h>
25 // --- Standard library ---
26
27 // --- AliRoot header files ---
28 #include "AliESDCaloCluster.h"
29 #include "AliESDEvent.h"
30 #include "AliQAChecker.h"
31 #include "AliLog.h"
32 #include "AliHMPIDDigit.h"
33 #include "AliHMPIDHit.h"
34 #include "AliHMPIDCluster.h"
35 #include "AliHMPIDQADataMakerSim.h"
36 #include "AliHMPIDParam.h"
37 #include "AliHMPIDRawStream.h"
38 #include "AliLog.h"
39
40 //.
41 // HMPID AliHMPIDQADataMakerSim base class
42 // for QA of simulation
43 // here also errors are calculated
44 //.
45
46 ClassImp(AliHMPIDQADataMakerSim)
47            
48 //____________________________________________________________________________ 
49   AliHMPIDQADataMakerSim::AliHMPIDQADataMakerSim() : 
50   AliQADataMakerSim(AliQA::GetDetName(AliQA::kHMPID), "HMPID Quality Assurance Data Maker")
51 {
52   // ctor
53 }
54
55 //____________________________________________________________________________ 
56 AliHMPIDQADataMakerSim::AliHMPIDQADataMakerSim(const AliHMPIDQADataMakerSim& qadm) :
57   AliQADataMakerSim() 
58 {
59   //copy ctor 
60   SetName((const char*)qadm.GetName()) ; 
61   SetTitle((const char*)qadm.GetTitle()); 
62 }
63
64 //__________________________________________________________________
65 AliHMPIDQADataMakerSim& AliHMPIDQADataMakerSim::operator = (const AliHMPIDQADataMakerSim& qadm )
66 {
67   // Equal operator.
68   this->~AliHMPIDQADataMakerSim();
69   new(this) AliHMPIDQADataMakerSim(qadm);
70   return *this;
71 }
72  
73 //____________________________________________________________________________ 
74 void AliHMPIDQADataMakerSim::InitHits()
75 {
76   // create Hits histograms in Hits subdir
77      TH1F *hHitQdc=new TH1F("HitQdc","HMPID Hit Qdc all chamber;QDC",500,0,4000);
78      Add2HitsList(hHitQdc,0);
79      TH2F *hHitMap[7];
80      for(Int_t iCh=0;iCh<7;iCh++) {
81      hHitMap[iCh]=new TH2F(Form("HMPID HitMap%i",iCh),Form("Ch%i;x_{Hit};y_{Hit}",iCh),162,-1,161,146,-1,145);   
82     Add2HitsList(hHitMap[iCh],iCh+1);
83     }
84
85 }
86
87 //____________________________________________________________________________ 
88 void AliHMPIDQADataMakerSim::InitDigits()
89 {
90   // create Digits histograms in Digits subdir
91       TH1F *hDigChEvt = new TH1F("hDigChEvt","Chamber occupancy per event",AliHMPIDParam::kMaxCh+1,AliHMPIDParam::kMinCh,AliHMPIDParam::kMaxCh+1);
92       TH1F *hDigPcEvt = new TH1F("hDigPcEvt","PC occupancy",156,-1,77);
93       TH2F *hDigMap[7];
94       TH1F *hDigQ[42];
95       for(Int_t iCh =0; iCh < 7; iCh++){
96        hDigMap[iCh] = new TH2F(Form("MapCh%i",iCh),Form("Digit Map in Chamber %i",iCh),159,0,159,143,0,143);
97        for(Int_t iPc =0; iPc < 6; iPc++ ){
98         hDigQ[iCh*6+iPc] = new TH1F(Form("QCh%iPc%i        ",iCh,iPc),Form("Charge of digits (ADC) in Chamber %i and PC %i   ",iCh,iPc),4100,0,4100);
99       }
100      }
101
102    Add2DigitsList(hDigChEvt,0);
103    Add2DigitsList(hDigPcEvt,1);
104    for(Int_t iMap=0; iMap < 7; iMap++) Add2DigitsList(hDigMap[iMap],2+iMap);
105    for(Int_t iH =0; iH < 42 ; iH++) Add2DigitsList(hDigQ[iH]    ,9+iH);
106 }
107
108 //____________________________________________________________________________ 
109 void AliHMPIDQADataMakerSim::InitSDigits()
110 {
111   // create SDigits histograms in SDigits subdir
112    TH1F   *hSDigits     = new TH1F("hHmpidSDigits",    "SDigits Q  distribution in HMPID",  500, 0., 5000.) ; 
113
114 Add2SDigitsList(hSDigits,0);
115 }
116
117 //____________________________________________________________________________ 
118
119 void AliHMPIDQADataMakerSim::MakeHits(TClonesArray * data)
120 {
121  //
122  //filling QA histos for Hits
123  //
124   TClonesArray * hits = dynamic_cast<TClonesArray *>(data) ; 
125   if (!hits){
126     AliError("Wrong type of hits container") ; 
127   } else {
128     TIter next(hits); 
129     AliHMPIDHit * hit ; 
130     while ( (hit = dynamic_cast<AliHMPIDHit *>(next())) ) {
131       if(hit->Pid()<500000) GetHitsData(0)->Fill(hit->Q()) ;
132       if(hit->Pid()<500000) GetHitsData(hit->Ch()+1)->Fill(hit->LorsX(),hit->LorsY());
133     }
134   } 
135
136 }
137 //___________________________________________________________________________
138 void AliHMPIDQADataMakerSim::MakeHits(TTree * data)
139 {
140 //
141 //Opening of the Hit TTree 
142 //
143  TClonesArray *pHits=new TClonesArray("AliHMPIDHit");  data->SetBranchAddress("HMPID",&pHits);
144   for(Int_t iEnt=0;iEnt<data->GetEntriesFast();iEnt++){//entries loop
145     data->GetEntry(iEnt);
146     MakeHits(pHits);
147   }//entries loop
148 }
149 //____________________________________________________________________________
150 void AliHMPIDQADataMakerSim::MakeDigits(TClonesArray * data)
151 {
152  //
153  //filling QA histos for Digits
154  //
155
156   TObjArray *chamber = dynamic_cast<TObjArray*>(data);
157   if ( !chamber) {
158     AliError("Wrong type of digits container") ; 
159   } else {
160     for(Int_t i =0; i< chamber->GetEntries(); i++)
161       {
162         TClonesArray * digits = dynamic_cast<TClonesArray*>(chamber->At(i)); 
163         GetDigitsData(0)->Fill(i,digits->GetEntriesFast()/(48.*80.*6.));
164         TIter next(digits); 
165         AliHMPIDDigit * digit; 
166         while ( (digit = dynamic_cast<AliHMPIDDigit *>(next())) ) {
167           GetDigitsData(1)->Fill(10.*i+digit->Pc(),1./(48.*80.));
168           GetDigitsData(2+i)->Fill(digit->PadChX(),digit->PadChY());
169           GetDigitsData(9+i*6+digit->Pc())->Fill(digit->Q());
170         }  
171       }
172   }
173 }
174 //___________________________________________________________________________
175 void AliHMPIDQADataMakerSim::MakeDigits(TTree * data)
176 {
177 //
178 //Opening the Digit Tree
179 //
180  TObjArray *pObjDig=new TObjArray(AliHMPIDParam::kMaxCh+1);
181   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
182     TClonesArray *pCA=new TClonesArray("AliHMPIDDigit");
183     pObjDig->AddAt(pCA,iCh);
184   }
185
186   pObjDig->SetOwner(kTRUE);
187
188   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
189     data->SetBranchAddress(Form("HMPID%i",iCh),&(*pObjDig)[iCh]);
190   }
191   data->GetEntry(0);
192
193    MakeDigits((TClonesArray *)pObjDig);
194 }
195 //____________________________________________________________________________
196
197 void AliHMPIDQADataMakerSim::MakeSDigits(TClonesArray * data)
198 {
199  //
200  //filling QA histos for SDigits
201  //
202   TClonesArray * sdigits = dynamic_cast<TClonesArray *>(data) ; 
203   if (!sdigits) {
204     AliError("Wrong type of sdigits container") ; 
205   } else {
206     TIter next(sdigits) ; 
207     AliHMPIDDigit * sdigit ; 
208     while ( (sdigit = dynamic_cast<AliHMPIDDigit *>(next())) ) {
209             GetSDigitsData(0)->Fill(sdigit->Q());
210     } 
211   }
212 }
213 //___________________________________________________________________________
214 void AliHMPIDQADataMakerSim::MakeSDigits(TTree * data)
215 {
216  //
217  // Opening the SDigit Tree
218  //
219  TClonesArray * sdigits = new TClonesArray("AliHMPIDDigit", 1000) ;
220
221   TBranch * branch = data->GetBranch("HMPID") ;
222   if ( ! branch ) {
223     AliError("HMPID SDigit Tree not found") ;
224     return;
225   }
226   branch->SetAddress(&sdigits) ;
227   branch->GetEntry(0) ;
228   MakeSDigits(sdigits) ;
229 }
230 //____________________________________________________________________________
231 void AliHMPIDQADataMakerSim::StartOfDetectorCycle()
232 {
233   //Detector specific actions at start of cycle
234   
235 }
236
237 void AliHMPIDQADataMakerSim::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray *obj)
238 {
239   //Detector specific actions at end of cycle
240   // do the QA checking
241   AliQAChecker::Instance()->Run(AliQA::kHMPID, task, obj) ;  
242 }
243