]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZEROQADataMakerRec.cxx
6c3e69bb6a816a117369219a51a9e5cf99e69362
[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 "AliVZEROQADataMakerRec.h"
35 #include "AliQAChecker.h"
36 #include "AliRawReader.h"
37 #include "AliVZERORawStream.h"
38 #include "AliVZEROReconstructor.h"
39
40
41 ClassImp(AliVZEROQADataMakerRec)
42            
43 //____________________________________________________________________________ 
44   AliVZEROQADataMakerRec::AliVZEROQADataMakerRec() : 
45   AliQADataMakerRec(AliQA::GetDetName(AliQA::kVZERO), "VZERO Quality Assurance Data Maker")
46 {
47   // constructor
48 }
49
50 //____________________________________________________________________________ 
51 AliVZEROQADataMakerRec::AliVZEROQADataMakerRec(const AliVZEROQADataMakerRec& qadm) :
52   AliQADataMakerRec()
53 {
54   //copy constructor 
55   SetName((const char*)qadm.GetName()) ; 
56   SetTitle((const char*)qadm.GetTitle()); 
57 }
58
59 //__________________________________________________________________
60 AliVZEROQADataMakerRec& AliVZEROQADataMakerRec::operator = (const AliVZEROQADataMakerRec& qadm )
61 {
62   // Equal operator
63   this->~AliVZEROQADataMakerRec();
64   new(this) AliVZEROQADataMakerRec(qadm);
65   return *this;
66 }
67  
68 //____________________________________________________________________________ 
69 void AliVZEROQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * list)
70 {
71   //Detector specific actions at end of cycle
72   // do the QA checking
73   AliQAChecker::Instance()->Run(AliQA::kVZERO, task, list) ;  
74 }
75
76 //____________________________________________________________________________ 
77 void AliVZEROQADataMakerRec::InitESDs()
78 {
79   //Create histograms to controll ESD
80  
81   TH1I * h1 = new TH1I("hVZERONbPMA", "Number of PMs fired in V0A", 80, 0, 79); 
82   h1->Sumw2() ;
83   Add2ESDsList(h1, 0)  ;  
84                                                                                                         
85   TH1I * h2 = new TH1I("hVZERONbPMC", "Number of PMs fired in V0C", 80, 0, 79); 
86   h2->Sumw2() ;
87   Add2ESDsList(h2, 1) ;
88  
89   TH1I * h3 = new TH1I("hVZEROMultA", "Multiplicity in V0A", 50, 0, 49) ; 
90   h3->Sumw2() ;
91   Add2ESDsList(h3, 2) ;
92  
93   TH1I * h4 = new TH1I("hVZEROMultC", "Multiplicity in V0C", 50, 0, 49) ; 
94   h4->Sumw2() ;
95   Add2ESDsList(h4, 3) ;
96
97   TH2F * h5 = new TH2F("fVzeroMult", "Vzero multiplicity",
98                          64, -0.5, 63.5,1000, -0.5, 99.5);
99   h5->GetXaxis()->SetTitle("Vzero PMT");
100   h5->GetYaxis()->SetTitle("Multiplicity");
101   h5->Sumw2() ;
102   Add2ESDsList(h5, 4) ;
103   TH1F * h6 = new TH1F("fBBA","BB Vzero A", 32, -0.5,31.5);
104   h6->Sumw2();
105   Add2ESDsList(h6, 5) ;
106   TH1F * h7 = new TH1F("fBGA","BG Vzero A", 32, -0.5,31.5);
107   h7->Sumw2();
108   Add2ESDsList(h7, 6) ;
109   TH1F * h8 = new TH1F("fBBC","BB Vzero C", 32, -0.5,31.5);
110   h8->Sumw2();
111   Add2ESDsList(h8, 7) ;
112   TH1F * h9 = new TH1F("fBGC","BG Vzero C", 32, -0.5,31.5);
113   h9->Sumw2();
114   Add2ESDsList(h9, 8) ;
115
116   TH2F *h10 = new TH2F("fVzeroAdc", "Vzero Adc",
117                          64, -0.5, 63.5,1024, -0.5, 1023.5);
118   h10->GetXaxis()->SetTitle("Vzero PMT");
119   h10->GetYaxis()->SetTitle("Adc counts");
120   h10->Sumw2() ;
121   Add2ESDsList(h10, 9);
122
123   TH2F *h11 = new TH2F("fVzeroTime", "Vzero Time",
124                          64, -0.5, 63.5,300, -0.5, 149.5);
125   h11->GetXaxis()->SetTitle("Vzero PMT");
126   h11->GetYaxis()->SetTitle("Time [100 ps]");
127   h11->Sumw2() ;
128   Add2ESDsList(h11,10);
129
130 }
131
132 //____________________________________________________________________________ 
133  void AliVZEROQADataMakerRec::InitRaws()
134  {
135    // create Raws histograms in Raws subdir
136
137   char ADCname[12]; 
138   char texte[40]; 
139   TH1I *fhRawADC0[64]; 
140   TH1I *fhRawADC1[64]; 
141   
142   TH2I * h0 = new TH2I("hCellADCMap0","ADC vs Cell for EVEN Integrator", 70, 0, 69, 512, 0, 1023);
143   h0->Sumw2(); 
144   Add2RawsList(h0,0) ;
145   TH2I * h1 = new TH2I("hCellADCMap1","ADC vs Cell for ODD Integrator", 70, 0, 69, 512, 0, 1023);
146   h1->Sumw2();
147   Add2RawsList(h1,1) ;
148                            
149   for (Int_t i=0; i<64; i++)
150     {
151        sprintf(ADCname,"hRaw0ADC%d",i);
152        sprintf(texte,"Raw ADC in cell %d for even integrator",i);
153        fhRawADC0[i]= new TH1I(ADCname,texte,1024,0,1023);       
154        Add2RawsList(fhRawADC0[i],i+2);
155        
156        sprintf(ADCname,"hRaw1ADC%d",i);
157        sprintf(texte,"Raw ADC in cell %d for odd integrator",i);
158        fhRawADC1[i]= new TH1I(ADCname,texte,1024,0,1023);       
159        Add2RawsList(fhRawADC1[i],i+2+64);       
160      }  
161  }
162
163 //____________________________________________________________________________
164 void AliVZEROQADataMakerRec::MakeESDs(AliESDEvent * esd)
165 {
166   // make QA data from ESDs
167
168   AliESDVZERO *esdVZERO=esd->GetVZEROData();
169    
170   if (esdVZERO) { 
171       GetESDsData(0)->Fill(esdVZERO->GetNbPMV0A());
172       GetESDsData(1)->Fill(esdVZERO->GetNbPMV0C());  
173       GetESDsData(2)->Fill(esdVZERO->GetMTotV0A());
174       GetESDsData(3)->Fill(esdVZERO->GetMTotV0C());  
175       for(Int_t i=0;i<64;i++) {
176         GetESDsData(4)->Fill((Float_t) i,(Float_t) esdVZERO->GetMultiplicity(i));
177         GetESDsData(9)->Fill((Float_t) i,(Float_t) esdVZERO->GetAdc(i));
178         GetESDsData(10)->Fill((Float_t) i,(Float_t) esdVZERO->GetTime(i));
179       }
180       for(Int_t i=0;i<32;i++) { 
181         if (esdVZERO->BBTriggerV0A(i)) 
182           GetESDsData(5)->Fill((Float_t) i);
183         if (esdVZERO->BGTriggerV0A(i)) 
184           GetESDsData(6)->Fill((Float_t) i);
185         if (esdVZERO->BBTriggerV0C(i)) 
186           GetESDsData(7)->Fill((Float_t) i);
187         if (esdVZERO->BGTriggerV0C(i)) 
188           GetESDsData(8)->Fill((Float_t) i);
189       }
190   }
191   
192 }
193
194 //____________________________________________________________________________
195  void AliVZEROQADataMakerRec::MakeRaws(AliRawReader* rawReader)
196  {
197    //Fill histograms with Raws
198        
199   AliVZERORawStream* rawStream  = new AliVZERORawStream(rawReader); 
200   rawStream->Next();
201     
202   for(Int_t i=0; i<64; i++) 
203      {
204        if(!rawStream->GetIntegratorFlag(i,10)) 
205          { 
206            // even integrator - fills index 2 to 65   
207            GetRawsData(0)->Fill(i,rawStream->GetADC(i)) ; 
208            GetRawsData(i+2)->Fill(rawStream->GetADC(i)) ; }
209        else 
210          { 
211            // odd integrator  - fills index 66 to 129   
212            GetRawsData(1)->Fill(i,rawStream->GetADC(i)) ; 
213            GetRawsData(i+2+64)->Fill(rawStream->GetADC(i)) ; }          
214       }          
215  }
216
217 //____________________________________________________________________________ 
218 void AliVZEROQADataMakerRec::StartOfDetectorCycle()
219 {
220   // Detector specific actions at start of cycle
221   
222 }