]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZEROQADataMakerRec.cxx
Make separate, specialized geometries for RPhi and RhoZ views.
[u/mrichter/AliRoot.git] / VZERO / AliVZEROQADataMakerRec.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 /*
18   Produces the data needed to calculate the quality assurance. 
19   All data must be mergeable objects.
20 */
21
22 // --- ROOT system ---
23 #include <TClonesArray.h>
24 #include <TFile.h> 
25 #include <TH1F.h> 
26 #include <TH1I.h> 
27 #include <TH2I.h> 
28
29 // --- Standard library ---
30
31 // --- AliRoot header files ---
32 #include "AliESDEvent.h"
33 #include "AliLog.h"
34 #include "AliCDBManager.h"
35 #include "AliCDBStorage.h"
36 #include "AliCDBEntry.h"
37 #include "AliVZEROQADataMakerRec.h"
38 #include "AliQAChecker.h"
39 #include "AliRawReader.h"
40 #include "AliVZERORawStream.h"
41 #include "AliVZEROReconstructor.h"
42
43
44 ClassImp(AliVZEROQADataMakerRec)
45            
46 //____________________________________________________________________________ 
47   AliVZEROQADataMakerRec::AliVZEROQADataMakerRec() : 
48     AliQADataMakerRec(AliQA::GetDetName(AliQA::kVZERO), "VZERO Quality Assurance Data Maker"),
49     fCalibData(GetCalibData()),
50     fEvent(0)
51     
52 {
53   // constructor
54   
55    for(Int_t i=0; i<64; i++){  
56        fEven[i] = 0;   
57        fOdd[i]  = 0;  }
58   
59    for(Int_t i=0; i<128; i++){  
60        fADC_Mean[i] = 0.0;   }  
61 }
62
63 //____________________________________________________________________________ 
64   AliVZEROQADataMakerRec::AliVZEROQADataMakerRec(const AliVZEROQADataMakerRec& qadm) :
65   AliQADataMakerRec(),
66     fCalibData(GetCalibData()),
67     fEvent(0)
68   
69 {
70   //copy constructor 
71    SetName((const char*)qadm.GetName()) ; 
72    SetTitle((const char*)qadm.GetTitle()); 
73 }
74
75 //__________________________________________________________________
76 AliVZEROQADataMakerRec& AliVZEROQADataMakerRec::operator = (const AliVZEROQADataMakerRec& qadm )
77 {
78   // Equal operator
79   this->~AliVZEROQADataMakerRec();
80   new(this) AliVZEROQADataMakerRec(qadm);
81   return *this;
82 }
83
84 //____________________________________________________________________________
85 AliVZEROCalibData* AliVZEROQADataMakerRec::GetCalibData() const
86
87 {
88   AliCDBManager *man = AliCDBManager::Instance();
89
90   AliCDBEntry *entry=0;
91
92   entry = man->Get("VZERO/Calib/Data");
93
94   // Retrieval of data in directory VZERO/Calib/Data:
95
96   AliVZEROCalibData *calibdata = 0;
97
98   if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
99   if (!calibdata)  AliFatal("No calibration data from calibration database !");
100
101   return calibdata;
102 }
103  
104 //____________________________________________________________________________ 
105 void AliVZEROQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * list)
106 {
107   //Detector specific actions at end of cycle
108   // do the QA checking
109   AliQAChecker::Instance()->Run(AliQA::kVZERO, task, list) ;
110 }
111
112 //____________________________________________________________________________ 
113 void AliVZEROQADataMakerRec::InitESDs()
114 {
115   //Create histograms to controll ESD
116  
117   TH1I * h1 = new TH1I("hVZERONbPMA", "Number of PMs fired in V0A", 80, 0, 80); 
118   h1->Sumw2() ;
119   Add2ESDsList(h1, 0)  ;  
120                                                                                                         
121   TH1I * h2 = new TH1I("hVZERONbPMC", "Number of PMs fired in V0C", 80, 0, 80); 
122   h2->Sumw2() ;
123   Add2ESDsList(h2, 1) ;
124  
125   TH1I * h3 = new TH1I("hVZEROMultA", "Multiplicity in V0A", 50, 0, 50) ; 
126   h3->Sumw2() ;
127   Add2ESDsList(h3, 2) ;
128  
129   TH1I * h4 = new TH1I("hVZEROMultC", "Multiplicity in V0C", 50, 0, 50) ; 
130   h4->Sumw2() ;
131   Add2ESDsList(h4, 3) ;
132
133   TH2F * h5 = new TH2F("fVzeroMult", "Vzero multiplicity",
134                          64, -0.5, 63.5,1000, -0.5, 99.5);
135   h5->GetXaxis()->SetTitle("Vzero PMT");
136   h5->GetYaxis()->SetTitle("Multiplicity");
137   h5->Sumw2() ;
138   Add2ESDsList(h5, 4) ;
139   TH1F * h6 = new TH1F("fBBA","BB Vzero A", 32, -0.5,31.5);
140   h6->Sumw2();
141   Add2ESDsList(h6, 5) ;
142   TH1F * h7 = new TH1F("fBGA","BG Vzero A", 32, -0.5,31.5);
143   h7->Sumw2();
144   Add2ESDsList(h7, 6) ;
145   TH1F * h8 = new TH1F("fBBC","BB Vzero C", 32, -0.5,31.5);
146   h8->Sumw2();
147   Add2ESDsList(h8, 7) ;
148   TH1F * h9 = new TH1F("fBGC","BG Vzero C", 32, -0.5,31.5);
149   h9->Sumw2();
150   Add2ESDsList(h9, 8) ;
151
152   TH2F *h10 = new TH2F("fVzeroAdc", "Vzero Adc",
153                          64, -0.5, 63.5,1024, -0.5, 1023.5);
154   h10->GetXaxis()->SetTitle("Vzero PMT");
155   h10->GetYaxis()->SetTitle("Adc counts");
156   h10->Sumw2() ;
157   Add2ESDsList(h10, 9);
158
159   TH2F *h11 = new TH2F("fVzeroTime", "Vzero Time",
160                          64, -0.5, 63.5,300, -0.5, 149.5);
161   h11->GetXaxis()->SetTitle("Vzero PMT");
162   h11->GetYaxis()->SetTitle("Time [100 ps]");
163   h11->Sumw2() ;
164   Add2ESDsList(h11,10);
165
166 }
167
168 //____________________________________________________________________________ 
169  void AliVZEROQADataMakerRec::InitRaws()
170  {
171    // create Raws histograms in Raws subdir
172
173   char ADCname[12]; 
174   char texte[40]; 
175   TH1I *fhRawADC0[64]; 
176   TH1I *fhRawADC1[64];
177   
178   TH2I *h0 = new TH2I("hCellADCMap0","ADC vs Cell for EVEN Integrator", 70, 0, 70, 512, 0, 1024);
179   h0->Sumw2(); 
180   Add2RawsList(h0,0) ;
181   TH2I *h1 = new TH2I("hCellADCMap1","ADC vs Cell for ODD Integrator", 70, 0, 70, 512, 0, 1024);
182   h1->Sumw2();
183   Add2RawsList(h1,1) ;
184
185   TH2F *fhMeanADC0 = new TH2F("hCellMeanADCMap0","Mean ADC vs cell for EVEN integrator",70,-0.5,69.5,512,-0.5,511.5);       
186   Add2RawsList(fhMeanADC0,2); 
187   TH2F *fhMeanADC1 = new TH2F("hCellMeanADCMap1","Mean ADC vs cell for ODD integrator",70,-0.5,69.5,512,-0.5,511.5);       
188   Add2RawsList(fhMeanADC1,3);  
189                            
190   for (Int_t i=0; i<64; i++)
191     {
192        sprintf(ADCname,"hRaw0ADC%d",i);
193        sprintf(texte,"Raw ADC in cell %d for even integrator",i);
194        fhRawADC0[i]= new TH1I(ADCname,texte,1024,0,1024);       
195        Add2RawsList(fhRawADC0[i],i+4);
196                     
197        sprintf(ADCname,"hRaw1ADC%d",i);
198        sprintf(texte,"Raw ADC in cell %d for odd integrator",i);
199        fhRawADC1[i]= new TH1I(ADCname,texte,1024,0,1024);       
200        Add2RawsList(fhRawADC1[i],i+4+64);                     
201      }  
202  }
203
204 //____________________________________________________________________________
205 void AliVZEROQADataMakerRec::MakeESDs(AliESDEvent * esd)
206 {
207   // make QA data from ESDs
208
209   AliESDVZERO *esdVZERO=esd->GetVZEROData();
210    
211   if (esdVZERO) { 
212       GetESDsData(0)->Fill(esdVZERO->GetNbPMV0A());
213       GetESDsData(1)->Fill(esdVZERO->GetNbPMV0C());  
214       GetESDsData(2)->Fill(esdVZERO->GetMTotV0A());
215       GetESDsData(3)->Fill(esdVZERO->GetMTotV0C());  
216       for(Int_t i=0;i<64;i++) {
217         GetESDsData(4)->Fill((Float_t) i,(Float_t) esdVZERO->GetMultiplicity(i));
218         GetESDsData(9)->Fill((Float_t) i,(Float_t) esdVZERO->GetAdc(i));
219         GetESDsData(10)->Fill((Float_t) i,(Float_t) esdVZERO->GetTime(i));
220       }
221       for(Int_t i=0;i<32;i++) { 
222         if (esdVZERO->BBTriggerV0A(i)) 
223           GetESDsData(5)->Fill((Float_t) i);
224         if (esdVZERO->BGTriggerV0A(i)) 
225           GetESDsData(6)->Fill((Float_t) i);
226         if (esdVZERO->BBTriggerV0C(i)) 
227           GetESDsData(7)->Fill((Float_t) i);
228         if (esdVZERO->BGTriggerV0C(i)) 
229           GetESDsData(8)->Fill((Float_t) i);
230       }
231   }
232   
233 }
234
235 //____________________________________________________________________________
236  void AliVZEROQADataMakerRec::MakeRaws(AliRawReader* rawReader)
237  {
238   //Fill histograms with Raws, computes average ADC values dynamically (pedestal subtracted)
239                   
240         rawReader->Reset() ; 
241         AliVZERORawStream* rawStream  = new AliVZERORawStream(rawReader); 
242   rawStream->Next();
243   
244   Float_t ChargeEoI, Threshold;  // for pedestal subtraction 
245                               
246 //  for(Int_t i=0; i<128; i++) { printf(" i = %d pedestal = %f sigma = %f \n\n", 
247 //                                        i,  fCalibData->GetPedestal(i),fCalibData->GetSigma(i) );} 
248              
249   GetRawsData(2)->Reset();      // to keep only the last value of the ADC average
250   GetRawsData(3)->Reset();
251
252   for(Int_t i=0; i<64; i++) {
253       if(!rawStream->GetIntegratorFlag(i,10)) 
254            { 
255             // even integrator - fills index 4 to 67 
256             GetRawsData(0)->Fill(i,rawStream->GetADC(i)) ; 
257             GetRawsData(i+4)->Fill(rawStream->GetADC(i)) ; 
258             ChargeEoI  = (float)rawStream->GetADC(i) - fCalibData->GetPedestal(i); 
259             Threshold  = 3.0 * fCalibData->GetSigma(i);   
260             if(rawStream->GetADC(i)<1023 && ChargeEoI > Threshold) {
261                fADC_Mean[i] = ((fADC_Mean[i]*fEven[i]) + rawStream->GetADC(i)) / (fEven[i]+1); 
262                GetRawsData(2)->Fill((Float_t) i, fADC_Mean[i]);
263                fEven[i]+=1; }
264            }
265       else if(rawStream->GetIntegratorFlag(i,10)) 
266            { 
267             // odd integrator  - fills index 68 to 131  
268             GetRawsData(1)->Fill(i,rawStream->GetADC(i)) ; 
269             GetRawsData(i+4+64)->Fill(rawStream->GetADC(i)) ;
270             ChargeEoI  = (float)rawStream->GetADC(i) - fCalibData->GetPedestal(i+64) ; 
271             Threshold  = 3.0 * fCalibData->GetSigma(i+64);  
272             if(rawStream->GetADC(i)<1023 && ChargeEoI > Threshold) { 
273                fADC_Mean[i+64] = ((fADC_Mean[i+64]*fOdd[i]) + rawStream->GetADC(i)) / (fOdd[i]+1);
274                GetRawsData(3)->Fill((Float_t) i, fADC_Mean[i+64]);
275                fOdd[i]+=1; }
276            }    
277   }    
278     fEvent++;                     
279  }
280
281 //____________________________________________________________________________ 
282 void AliVZEROQADataMakerRec::StartOfDetectorCycle()
283 {
284   // Detector specific actions at start of cycle
285   
286 }