]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0QADataMakerRec.cxx
7c16fecb2a57a97a9d17530ee3cfd24cdb63ca2d
[u/mrichter/AliRoot.git] / T0 / AliT0QADataMakerRec.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 //---
20 //  Produces the data needed to calculate the quality assurance. 
21 //  Alla.Maevskaya@cern.ch
22 //---
23
24 // --- ROOT system ---
25 #include <TClonesArray.h>
26 #include <TFile.h> 
27 #include <TH1F.h> 
28 #include <TH2F.h> 
29 #include <TDirectory.h>
30 #include <TMath.h>
31 // --- Standard library ---
32
33 // --- AliRoot header files ---
34 #include "AliESDEvent.h"
35 #include "AliLog.h"
36 #include "AliT0digit.h" 
37 #include "AliT0hit.h"
38 #include "AliT0RecPoint.h"
39 #include "AliT0QADataMakerRec.h"
40 #include "AliQAChecker.h"
41 #include "AliT0RawReader.h"
42 #include "AliT0RecoParam.h"
43 #include "THnSparse.h"
44
45 #include "Riostream.h"
46 ClassImp(AliT0QADataMakerRec)
47            
48 //____________________________________________________________________________ 
49   AliT0QADataMakerRec::AliT0QADataMakerRec() : 
50 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kT0), 
51                   "T0 Quality Assurance Data Maker"),
52   fnEventCal(0),
53   fnEventPhys(0)
54 {
55   // ctor
56   for (Int_t i=0; i<6; i++) {
57     fNumTriggers[i]=0;
58     fNumTriggersCal[i]=0;
59     fTrEffCal[i] = 0;
60     fTrEffPhys[i] = 0;
61
62   }
63   for (Int_t i=0; i<24; i++)
64     {
65       feffC[i]=0;
66       feffA[i]=0;
67       feffqtc[i]=0;
68       feffPhysC[i]=0;
69       feffPhysA[i]=0;
70       feffqtcPhys[i]=0;
71
72    }
73
74   // for(Int_t ic=0; ic<24; ic++) 
75   // fhTimeDiff[ic] = new TH1F(Form("CFD1minCFD%d",ic+1),"CFD-CFD",100,-250,250);
76 }
77
78
79 //____________________________________________________________________________ 
80 AliT0QADataMakerRec::AliT0QADataMakerRec(const AliT0QADataMakerRec& qadm) :
81   AliQADataMakerRec(),
82   fnEventCal(0),
83   fnEventPhys(0)
84 {
85   //copy ctor 
86  SetName((const char*)qadm.GetName()) ; 
87   SetTitle((const char*)qadm.GetTitle()); 
88 }
89
90 //__________________________________________________________________
91 AliT0QADataMakerRec& AliT0QADataMakerRec::operator = (const AliT0QADataMakerRec& qadm )
92 {
93   // Equal operator.
94   this->~AliT0QADataMakerRec();
95   new(this) AliT0QADataMakerRec(qadm);
96   return *this;
97 }
98 //__________________________________________________________________
99 AliT0QADataMakerRec::~AliT0QADataMakerRec()
100 {
101   //destructor
102   //  for(Int_t ic=0; ic<24; ic++) {
103   ////    if (fhTimeDiff[ic]) delete fhTimeDiff[ic];
104   //  }
105 }
106 //____________________________________________________________________________
107 void AliT0QADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
108 {
109   //Detector specific actions at end of cycle
110   // do the QA checking
111   AliQAChecker::Instance()->Run(AliQAv1::kT0, task, list) ;
112   
113   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
114     if (! IsValidEventSpecie(specie, list)) 
115       continue ;
116     SetEventSpecie(AliRecoParam::ConvertIndex(specie)) ; 
117     if ( task == AliQAv1::kRAWS ) {
118       GetRawsData(0)->SetLabelSize(0.02);
119       const Char_t *triggers[6] = {"mean", "vertex","ORA","ORC","central","semi-central"};
120       for (Int_t itr=0; itr<6; itr++) {
121         if ( fnEventCal>0) 
122           fTrEffCal[itr] = Float_t (fNumTriggersCal[itr])/Float_t (fnEventCal);
123         if ( fnEventPhys>0) 
124           fTrEffPhys[itr] = Float_t (fNumTriggers[itr])/Float_t (fnEventPhys);
125
126         GetRawsData(169+250)->Fill(triggers[itr], fTrEffCal[itr]);
127         GetRawsData(169+250)->SetBinContent(itr+1, fTrEffCal[itr]);
128         GetRawsData(169)->Fill(triggers[itr], fTrEffPhys[itr]);
129         GetRawsData(169)->SetBinContent(itr+1, fTrEffPhys[itr]);
130       } 
131       Float_t effic=0;
132       for(Int_t ik=0; ik<24; ik++)
133         {  
134           effic=0;
135           if ( fnEventCal>0) effic = Float_t(feffC[ik])/Float_t(fnEventCal);
136           GetRawsData(207+250)->SetBinContent(ik+1,effic) ;
137           effic=0;
138           if ( fnEventCal>0) effic = Float_t(feffA[ik])/Float_t(fnEventCal);
139           GetRawsData(208+250)->SetBinContent(ik+1,effic );
140           effic=0;
141           if ( fnEventCal>0) effic = Float_t(feffqtc[ik])/Float_t(fnEventCal);
142           GetRawsData(209+250)->SetBinContent(ik+1, effic);
143
144           effic=0;
145           if ( fnEventPhys>0) effic = Float_t(feffPhysC[ik])/Float_t(fnEventPhys);
146           GetRawsData(207)->SetBinContent(ik+1, effic);
147
148         }
149       
150       
151     }
152     
153   }
154
155 }
156 //____________________________________________________________________________
157 void AliT0QADataMakerRec::StartOfDetectorCycle()
158 {
159   //Detector specific actions at start of cycle
160
161 }
162  
163 //____________________________________________________________________________ 
164 void AliT0QADataMakerRec::InitRaws()
165 {
166   // create Raw histograms in Raw subdir
167   const Bool_t expert   = kTRUE ; 
168   const Bool_t saveCorr = kTRUE ; 
169   const Bool_t image    = kTRUE ; 
170   Float_t low[500];
171   Float_t high[500];
172
173   
174   for (Int_t i=0; i<500; i++){
175     low[i] = 0;
176     high[i] = 30000;
177
178   }
179
180   TString timename, ampname, qtcname, ledname;
181   TString timeCalname, ampCalname, ledCalname, qtcCalname;
182   TString qt1name, qt0name, qt1Calname, qt0Calname;
183   TString nhits;
184
185   TH1F* fhRefPoint = new TH1F("hRefPoint","Ref Point", 10000, 0 ,50000);
186   Add2RawsList( fhRefPoint,0, expert, !image, !saveCorr);
187
188   TH1F* fhRefPointcal = new TH1F("hRefPointcal","Ref Point laser", 5000, 0 ,20000);
189   Add2RawsList( fhRefPointcal,250, expert, !image, !saveCorr);
190
191   TH1F *fhRawCFD[24]; 
192   TH1F * fhRawLEDamp[24];
193   TH1F *fhRawQTC[24]; TH1F * fhRawLED[24];
194   TH1F *fhRawCFDcal[24]; TH1F * fhRawLEDampcal[24]; 
195    TH1F *fhRawQTCcal[24];  TH1F * fhRawLEDcal[24];
196   TH1F *fhRawQT0cal[24]; TH1F *fhRawQT1cal[24];
197   TH1F *fhRawQT1[24]; TH1F *fhRawQT0[24];
198   TH1F* fhRawNhits[24];
199  
200   for (Int_t i=0; i<24; i++)
201     {
202       timename ="hRawCFD";
203       ledname = "hRawLED";
204       qtcname = "hRawQTC";
205       qt0name = "hRawQT0_";
206       qt1name = "hRawQT1_";
207       ampname = "hRawLEDminCFD";
208       nhits = "hRawNhits";
209       timename += i+1;
210       ampname += i+1;
211       qtcname += i+1;
212       qt0name += i+1;
213       qt1name += i+1;
214       ledname += i+1;
215       nhits   += i+1;
216       fhRawCFD[i] = new TH1F(timename.Data(), Form("%s;CFD [#channels];Counts", timename.Data()),Int_t((high[i+1]-low[i+1])/4),low[i+1],high[i+1]);
217       Add2RawsList( fhRawCFD[i],i+1, expert, !image, !saveCorr);
218       fhRawLED[i] = new TH1F(ledname.Data(),  Form("%s;LED[#channels];Counts", ledname.Data()),Int_t((high[i+25]-low[i+25])/4),low[i+25],high[i+25]);
219       Add2RawsList( fhRawLED[i],i+25, expert, !image, !saveCorr);
220       fhRawLEDamp[i] = new TH1F(ampname.Data(),  Form("%s;LED-CFD [#channels];Counts", ampname.Data()),1000,0,1000);
221       Add2RawsList( fhRawLEDamp[i],i+49, expert, !image, !saveCorr);
222       fhRawQTC[i] = new TH1F(qtcname.Data(),  Form("%s;QTC[#channels];Counts", qtcname.Data()),10000,0,10000);
223       Add2RawsList( fhRawQTC[i],i+73, expert, !image, !saveCorr);
224       fhRawQT1[i] = new TH1F(qt1name.Data(),  Form("%s;QT1[#channels];Counts", qt1name.Data()),Int_t((high[97+i]-low[97+i])/4),low[97+i],high[97+i]);
225       Add2RawsList( fhRawQT1[i],97+i, expert, !image, !saveCorr);
226       fhRawQT0[i] = new TH1F(qt0name.Data(),  Form("%s;QT0[#channels];Counts", qt0name.Data()),Int_t((high[121+i]-low[121+i])/4),low[121+i],high[121+i]);
227       Add2RawsList( fhRawQT0[i],121+i, expert, !image, !saveCorr);
228       
229       fhRawNhits[i] = new TH1F(nhits.Data(),  Form("%s;#Hits;Events", nhits.Data()),20, 0, 20);
230       Add2RawsList( fhRawNhits[i],176+i, expert, !image, !saveCorr);
231       
232       
233       timeCalname ="hRawCFDcal";
234       ledCalname = "hRawLEDcal";
235       ampCalname = "hRawLEDminCFDcal";
236       qtcCalname = "hRawQTCcal";
237       qt0Calname = "hRawQT0cal";
238       qt1Calname = "hRawQT1cal";
239       timeCalname += i+1;
240       ledCalname += i+1;
241       ampCalname += i+1;
242       qtcCalname += i+1;
243       qt0Calname += i+1;
244       qt1Calname += i+1;
245       fhRawCFDcal[i] = new TH1F(timeCalname.Data(),  Form("LASER: %s;Time ;Counts", timeCalname.Data()),Int_t((high[251+i]-low[251+i])/4),low[251+i],high[251+i]);
246       Add2RawsList( fhRawCFDcal[i],251+i, expert, !image, !saveCorr);
247       
248       fhRawLEDcal[i] = new TH1F(ledCalname.Data(),  Form("LASER: %s;Time ;Counts", ledCalname.Data()),Int_t((high[251+i+24]-low[251+i+24])/4),low[251+i+24],high[251+i+24]);
249       Add2RawsList( fhRawLEDcal[i],251+i+24, expert, !image, !saveCorr);
250       
251       fhRawLEDampcal[i] = new TH1F(ampCalname.Data(), Form("LASER: %s;Amplitude [ADC counts];Counts", ampCalname.Data()),1000,0,1000);
252       Add2RawsList( fhRawLEDampcal[i],251+i+48, expert, !image, !saveCorr);
253       
254       fhRawQTCcal[i] = new TH1F(qtcCalname.Data(), Form("LASER: %s;QTC[#channels]; ;Counts",qtcCalname.Data()),10000,0,10000);
255       Add2RawsList( fhRawQTCcal[i],251+i+72, expert, !image, !saveCorr);
256       
257       fhRawQT0cal[i] = new TH1F(qt0Calname.Data(), Form("LASER: %s;QT0[#channels] ;Counts",qt0Calname.Data()),Int_t((high[251+96+i]-low[251+96+i])/4),low[251+96+i],high[251+96+i]);
258       Add2RawsList( fhRawQT0cal[i],251+96+i, expert, !image, !saveCorr);
259       
260       fhRawQT1cal[i] = new TH1F(qt1Calname.Data(), Form("LASER: %s;QT1[#channels] ;Counts",qt1Calname.Data()),Int_t((high[i+251+120]-low[i+251+120])/4),low[i+251+120],high[i+251+120]);
261       Add2RawsList( fhRawQT1cal[i],i+251+120, expert, !image, !saveCorr);
262       
263     }
264   
265   
266   TH1F* fhRawTrigger = new TH1F("hRawTrigger"," phys triggers;Trigger ;Counts",6,0,6);
267   Add2RawsList(fhRawTrigger ,169, expert, image, !saveCorr);
268   
269   TH1F* fhRawMean = new TH1F("hRawMean","online mean signal, physics event;",Int_t((high[170]-low[170])/4),low[170],high[170]);
270   Add2RawsList( fhRawMean,170, expert, !image, !saveCorr);
271
272   TH1F* fhRawVertex = new TH1F("hRawVertex","online vertex signal; counts",Int_t((high[171]-low[171])/4),low[171],high[171]);
273   Add2RawsList( fhRawVertex,171, expert, image, !saveCorr);
274
275   TH1F* fhRawORA = new TH1F("hRawORA","online OR A; counts",Int_t((high[172]-low[172])/4),low[172],high[172]);
276   Add2RawsList( fhRawORA,172, expert, !image, !saveCorr);
277   TH1F* fhRawORC = new TH1F("hRawORC","online OR C;counts",Int_t(( high[173]-low[173])/4),low[173],high[173]);
278   Add2RawsList( fhRawORC,173, expert, !image, !saveCorr);
279   TH1F* fhMultCentr = new TH1F("hMultCentr","online trigger Central;counts ",Int_t(( high[174]-low[174])/4),low[174],high[174]);
280   Add2RawsList( fhMultCentr,174, expert, !image, !saveCorr);
281   TH1F* fhMultSeCentr = new TH1F("hMultSemiCentr","online trigger SemiCentral;counts ",Int_t(( high[175]-low[175])/4),low[175],high[175]);
282   Add2RawsList( fhMultSeCentr,175, expert, !image, !saveCorr);
283
284
285   TH1F* fhRawTriggerCal = new TH1F("hRawTriggerCal"," laser triggers",6,0,6);
286   Add2RawsList(fhRawTriggerCal ,169+250 , !expert, image, !saveCorr);
287   
288   TH1F* fhRawMeanCal = new TH1F("hRawMeanCal","online mean signal, calibration event",Int_t((high[170+250]-low[170+250])/4),low[170+250],high[170+250]);
289   Add2RawsList( fhRawMeanCal,170+250, expert, !image, !saveCorr);
290   
291   TH1F* fhRawVertexCal = new TH1F("hRawVertexCal","online vertex signal, calibration event ",Int_t((high[171+250]-low[171+250])/4),low[171+250],high[171+250] );
292   Add2RawsList( fhRawVertexCal,171+250, expert, !image, !saveCorr);
293   
294   
295   TH1F* fhRawORAcal = new TH1F("hRawORAcal","laser OR A; counts",Int_t((high[172+250]-low[172+250])/4),low[172+250],high[172+250]);
296   Add2RawsList( fhRawORAcal,172+250, expert, !image, !saveCorr );
297   
298   
299   TH1F* fhRawORCcal = new TH1F("hRawORCcal","laserOR C;counts ",Int_t(( high[173]-low[173])/4),low[173],high[173]);
300   Add2RawsList( fhRawORCcal,173+250, expert, !image, !saveCorr);
301   
302   TH1F* fhMultCentrcal = new TH1F("hMultCentrcal","laser trigger Central;counts ",Int_t(( high[174]-low[174])/4),low[174],high[174]);
303   Add2RawsList( fhMultCentrcal,174+250, expert, !image, !saveCorr);
304   TH1F* fhMultSeCentrcal = new TH1F("hMultSemiCentrcal","laser trigger SemiCentral;counts ",Int_t(( high[175]-low[175])/4),low[175],high[175]);
305   Add2RawsList( fhMultSeCentrcal,175+250, expert, !image, !saveCorr);
306   
307   //multiplicity trigger
308   //side A
309   TH1F* fhMultAcal = new TH1F("hMultAcal","laser: full mulltiplicity;Multiplicity A side;Entries",Int_t((high[201]-low[201])/4),low[201],high[201]);
310   Add2RawsList( fhMultAcal,201+250, expert, !image, !saveCorr );
311   TH1F* fhMultAScal = new TH1F("hMultASemical","laser:full multiplicity with semi-central trigger A side;Multiplicity;Entries",
312                                Int_t((high[202]-low[202])/4),low[202],high[202] );
313   Add2RawsList( fhMultAScal,202+250, expert, !image, !saveCorr);
314   TH1F* fhMultACcal = new TH1F("hMultACentrcal","laser:full multiplicity with central trigger A side;Multiplicity;Entries", 
315                                Int_t((high[203]-low[203])/4),low[203],high[203]);
316   Add2RawsList( fhMultACcal,203+250, expert, !image, !saveCorr);
317   
318   TH1F* fhMultA = new TH1F("hMultA","full mulltiplicity A side;Multiplicity;Entries", Int_t((high[201]-low[201])/4) ,low[201],high[201]);
319   Add2RawsList( fhMultA,201, expert, image, !saveCorr );
320   
321   TH1F* fhMultAS = new TH1F("hMultASemi","full multiplicity with semi-central trigger A side ;Multiplicity;Entries",
322                             Int_t((high[202]-low[202])/4),low[202],high[202] );
323   Add2RawsList( fhMultAS, 202, expert, !image, !saveCorr);
324   TH1F* fhMultAC = new TH1F("hMultACentr","full multiplicity with central trigger;Multiplicity;Entries", 
325                             Int_t((high[203]-low[203])/4),low[203],high[203]);
326   Add2RawsList( fhMultAC, 203, expert, !image, !saveCorr);
327   
328   
329   //side C
330   TH1F* fhMultCcal = new TH1F("hMultCcal","laser:full mulltiplicity C side;Multiplicity;Entries",Int_t((high[204]-low[204])/4),low[204],high[204]);
331   Add2RawsList( fhMultCcal,204+250, expert, !image, !saveCorr );
332   TH1F* fhMultCScal = new TH1F("hMultCSemical","laser:full multiplicity with semi-central trigger C side;Multiplicity;Entries",
333                                Int_t((high[205]-low[205])/4),low[205],high[205] );
334   Add2RawsList( fhMultCScal,205+250, expert, !image, !saveCorr);
335   TH1F* fhMultCCcal = new TH1F("hMultCCentrcal","laser:full multiplicity with central trigger C side;Multiplicity;Entries", 
336                                Int_t((high[206]-low[206])/4),low[206],high[206]);
337   Add2RawsList( fhMultCCcal,206+250, expert, !image, !saveCorr);
338   
339   TH1F* fhMultC = new TH1F("hMultC","full mulltiplicity C side;Multiplicity;Entries", Int_t(high[204]-low[204]/4) ,low[204],high[204]);
340   Add2RawsList( fhMultC,204, expert, image, !saveCorr );
341   TH1F* fhMultCS = new TH1F("hMultCSemi","full multiplicity with semi-central trigger C side;Multiplicity;Entries",
342                             Int_t((high[205]-low[205])/4),low[205],high[205] );
343   Add2RawsList( fhMultCS,205, expert, !image, !saveCorr);
344   TH1F* fhMultCC = new TH1F("hMultCentr","full multiplicity with central trigger C side;Multiplicity;Entries", 
345                             Int_t((high[206]-low[206])/4),low[206],high[206]);
346   Add2RawsList( fhMultCC,206, expert, !image, !saveCorr);
347   
348   
349   //efficiency
350   TH1F* fhEffCFD = new TH1F("hEffCFDcal","CFD efficiecy laser ;#PMT; #CFD counts/nEvents",24, 0 ,24); 
351   Add2RawsList( fhEffCFD,207+250, !expert, image, !saveCorr);
352   
353   TH1F* fhCFDeffpsys= new TH1F("fhCFDeffpsys"," CFD efficiency; #PMT; #CFD counts/nEvents",24, 0 ,24);  
354   // fhCFDeffpsys->SetMaximum(2);
355   Add2RawsList( fhCFDeffpsys, 207, expert, image, !saveCorr);
356   
357   TH1F* fhEffLED = new TH1F("hEffLEDcal","LEDefficiecy; #PMT; #LED counts/nEvent",24, 0 ,24);
358   Add2RawsList( fhEffLED,208+250, !expert, image, !saveCorr);
359   
360   TH1F* fhEffQTC = new TH1F("hEffQTCcal","QTC efficiecy; #PMT; QTC efficiency%s;",24, 0 ,24);
361   Add2RawsList( fhEffQTC,209+250, !expert, image, !saveCorr);
362   
363   
364  TH2F* fhCFD = new TH2F("hCFD","CFD phys; #PMT; CFD {#channnels}",25, 0 ,25,Int_t((high[210]-low[210])/4),low[210],high[210]);
365   fhCFD->SetOption("COLZ");
366   Add2RawsList( fhCFD,210, expert, image, !saveCorr);
367     
368   TH2F* fhLED = new TH2F("hLED","LED phys; #PMT; LED [#channnels]",25, 0 ,25,Int_t((high[211]-low[211])/4),low[211],high[211]);
369   fhLED->SetOption("COLZ");
370   Add2RawsList( fhLED,211, expert, image, !saveCorr);
371
372   TH2F* fhQTC = new TH2F("hQTC","QTC phys; #PMT; QTC [#channnels]",25, 0 ,25,Int_t( high[212]-low[212]),low[212],high[212]);
373   fhQTC->SetOption("COLZ");
374   Add2RawsList( fhQTC,212, expert, image, !saveCorr);
375
376    TH2F* fhCFDcal = new TH2F("hCFDcal","CFD laser; #PMT; CFD {#channnels}",25, 0 ,25,Int_t((high[210]-low[210])/4),low[210],high[210]);
377   fhCFDcal->SetOption("COLZ");
378   Add2RawsList( fhCFDcal,210+250, expert, image, !saveCorr);
379   
380   
381   TH2F* fhLEDcal = new TH2F("hLEDcal","LED laser; #PMT; LED [#channnels]",25, 0 ,25,Int_t((high[211]-low[211])/4),low[211],high[211]);
382   fhLEDcal->SetOption("COLZ");
383   Add2RawsList( fhLEDcal,211+250, expert, image, !saveCorr);
384   
385   TH2F* fhQTCcal = new TH2F("hQTCcal","QTC laser; #PMT; QTC [#channnels]",25, 0 ,25,Int_t( high[212]-low[212]),low[212],high[212]);
386   fhQTCcal->SetOption("COLZ");
387   Add2RawsList( fhQTCcal,212+250, expert, image, !saveCorr);
388   
389   
390   TH1F* fhNumPMTA= new TH1F("hNumPMTA","number of PMT hitted per event",13, 0 ,13);
391   Add2RawsList(fhNumPMTA ,213, expert, image, !saveCorr);
392   
393   TH1F* fhNumPMTC= new TH1F("hNumPMTC","number of PMT hitted per event",13, 0 ,13);
394   Add2RawsList(fhNumPMTC ,214, expert, image, !saveCorr);
395   
396   TH1F* fhHitsOrA= new TH1F("fhHitsOrA","T0_OR A hit multiplicitie",20, 0 ,20);
397   Add2RawsList( fhHitsOrA,215, expert, !image, !saveCorr);
398   
399   TH1F* fhHitsOrC= new TH1F("fhHitsOrC","T0_OR C hit multiplicitie",20, 0 ,20);
400   Add2RawsList(fhHitsOrC ,216, expert, !image, !saveCorr);
401   
402   
403   TH1F* fhOrCminOrA= new TH1F("fhOrCminOrA","T0_OR C - T0_OR A",10000,-5000,5000);
404   Add2RawsList( fhOrCminOrA,219, expert, !image, !saveCorr);
405   
406   TH1F* fhOrCminOrAcal= new TH1F("fhOrCminOrAcal","T0_OR C - T0_OR A",10000,-5000,5000);
407   Add2RawsList( fhOrCminOrAcal,219+250, expert, !image, !saveCorr);
408
409   TH1F* fhOrCminOrATvdcOn= new TH1F("fhOrCminOrATvdcOn","T0_OR C - T0_OR A TVDC on",10000,-5000,5000);
410   Add2RawsList( fhOrCminOrATvdcOn,217, expert, !image, !saveCorr);
411   
412   TH1F* fhOrCminOrATvdcOncal= new TH1F("fhOrCminOrATvdcOncal","T0_OR C - T0_OR A TVDC on laser",10000,-5000,5000);
413   Add2RawsList( fhOrCminOrATvdcOncal,217+250, expert, !image, !saveCorr);
414
415   TH1F* fhOrCminOrATvdcOff= new TH1F("fhOrCminOrATvdcOff","T0_OR C - T0_OR A TVDC off",10000,-5000,5000);
416   Add2RawsList( fhOrCminOrATvdcOff,218, expert, !image, !saveCorr);
417
418
419    TH1F* fhOrCminOrATvdcOffcal= new TH1F("fhOrCminOrATvdcOffcal","T0_OR C - T0_OR ATVDC off laser",10000,-5000,5000);
420    Add2RawsList( fhOrCminOrATvdcOffcal,218+250, expert, !image, !saveCorr);
421    
422    TH2F* fhBeam = new TH2F("fhBeam", " Mean vs Vertex ", 120, -30, 30, 120, -30, 30);
423    Add2RawsList( fhBeam,220, !expert, image, !saveCorr);
424    TH2F* fhBeamTVDCon = new TH2F("fhBeamTVDCon", " Mean vs Vertex TVDC on ", 120, -30, 30, 120, -30, 30);
425    Add2RawsList( fhBeamTVDCon,221, expert, image, !saveCorr);
426    TH2F* fhBeamTVDCoff = new TH2F("fhBeamTVDCoff", " Mean vs Vertex TVDC off", 120, -30, 30, 120, -30, 30);
427    Add2RawsList( fhBeamTVDCoff,222, expert, image, !saveCorr);
428    
429    const Char_t *triggers[6] = {"mean", "vertex","ORA","ORC","central","semi-central"};
430    for (Int_t itr=0; itr<6; itr++) {
431      GetRawsData(169)->Fill(triggers[itr], fNumTriggersCal[itr]);
432      GetRawsData(169)->SetBinContent(itr+1, fNumTriggersCal[itr]);
433      GetRawsData(169+250)->Fill(triggers[itr], fNumTriggers[itr]);
434      GetRawsData(169+250)->SetBinContent(itr+1, fNumTriggers[itr]);
435    }
436   
437     
438 }
439   
440 //____________________________________________________________________________ 
441 void AliT0QADataMakerRec::InitDigits()
442 {
443   // create Digits histograms in Digits subdir
444   const Bool_t expert   = kTRUE ; 
445   const Bool_t image    = kTRUE ; 
446   
447   TH2F * fhDigCFD = new TH2F("fhDigCFD", " CFD digits; #PMT; CFD digits[#channels]",25,-0.5,24.5,100,0,1000);
448   fhDigCFD->SetOption("COLZ");
449   Add2DigitsList( fhDigCFD,0, !expert, image);
450   TH2F *fhDigLEDamp = new TH2F("fhDigLEDamp", " LED-CFD digits; #PMT; LED-CFD amplitude ",25,-0.5,24.5,100,100,1000);
451   fhDigLEDamp->SetOption("COLZ");
452   Add2DigitsList( fhDigLEDamp,1, !expert, !image);
453   TH2F * fhDigQTC = new TH2F("fhDigQTC", " QTC digits; #PMT; QTC amplitude",25,-0.5,24.5,100,100,10000);
454   fhDigQTC->SetOption("COLZ");
455   Add2DigitsList( fhDigQTC,2, !expert, !image);
456
457
458 }
459
460 //____________________________________________________________________________ 
461
462 void AliT0QADataMakerRec::InitRecPoints()
463 {
464   // create cluster histograms in RecPoint subdir
465   const Bool_t expert   = kTRUE ; 
466   const Bool_t image    = kTRUE ; 
467
468   TH2F* fhRecCFD = new TH2F("hRecCFD"," CFD time;#PMT; CFD Time [ns];",24, 0 ,24, 
469                               100,-50,50);
470   fhRecCFD->SetOption("COLZ");
471   Add2RecPointsList ( fhRecCFD,0, !expert, image);
472
473   TH2F* fhRecAmpDiff = new TH2F("hRecAmpDiff"," LED-CFD  min QTC amplitude;#PMT; difference [MIPs];",
474                                 24, 0 ,24, 200,-10,10);
475   fhRecAmpDiff->SetOption("COLZ");
476   Add2RecPointsList (fhRecAmpDiff, 1, !expert, image);
477   
478   TH1F *fhMean = new TH1F("hMean","online - rec mean;online - rec mean[#channels];",2000, -1000, 1000);
479   Add2RecPointsList ( fhMean,2, !expert, image);
480
481 }
482
483 //____________________________________________________________________________
484 void AliT0QADataMakerRec::InitESDs()
485 {
486   //create ESDs histograms in ESDs subdir
487   const Bool_t expert   = kTRUE ; 
488   const Bool_t image    = kTRUE ; 
489   
490   TH1F *fhESDMean = new TH1F("hESDmean"," ESD mean; mean time[%channels]",1000, -5, 5);
491   Add2ESDsList(fhESDMean, 0, expert, !image) ;
492   TH1F * fhESDVertex = new TH1F("hESDvertex","ESDvertex; vertex[cm];",82,-30,30);
493   Add2ESDsList(fhESDVertex, 1, expert, !image) ;
494   
495   TH1F * fhESDResolution = new TH1F("hESDResolution","(T0A-T0C)/2 corrected by SPD vertex; ns",800,-2,2);
496   Add2ESDsList(fhESDResolution, 2, !expert, image) ;
497
498 }
499
500 //____________________________________________________________________________
501 void AliT0QADataMakerRec::MakeRaws( AliRawReader* rawReader)
502 {
503   Int_t  time[24] ;
504   for(Int_t i=0; i<24; i++) time[i] = 0;          
505   
506   
507   rawReader->Reset() ; 
508   //fills QA histos for RAW
509   Int_t shift=0;
510   // Int_t refPointParam = GetRecoParam()->GetRefPoint();
511   Int_t refpoint = 0;
512   Int_t refPointParam = 0;
513
514   AliT0RawReader *start = new AliT0RawReader(rawReader);
515   
516   if (! start->Next())
517     AliDebug(AliQAv1::GetQADebugLevel(),Form(" no raw data found!!"));
518   else
519     {  
520       UInt_t type =rawReader->GetType();
521       if (type == 8){ shift=250; fnEventCal++;} 
522       if (type == 7){ shift=0;   fnEventPhys++;}
523       //    if (type == 7){ shift=1;   fnEventPhys++;}
524       Int_t allData[110][5];
525       for (Int_t i0=0; i0<110; i0++)
526         {
527           for (Int_t j0=0; j0<5; j0++) allData[i0][j0]=0;
528         }
529       for (Int_t i=0; i<107; i++) 
530         for (Int_t iHit=0; iHit<5; iHit++)
531           allData[i][iHit]= start->GetData(i,iHit);
532       
533       if ( allData[0][0] > 0  && (type == 7))
534         GetRawsData(0) -> Fill( allData[0][0]);
535      if ( allData[0][0] > 0  && (type == 8))
536         GetRawsData(250) -> Fill( allData[0][0]);
537       
538       refpoint = allData[refPointParam][0];
539       if (refPointParam <  0 ) refpoint=0; 
540       if (refPointParam == 0 ) refpoint = allData[0][0] - 5000;
541       
542       Int_t sideshift, sideshiftqtc;
543       Int_t numPmtC=0;    
544       Int_t numPmtA=0;    
545       for (Int_t ik = 0; ik<24; ik++)
546         {
547           if(ik<12) {
548             sideshift=1;
549             sideshiftqtc=1;
550             if(allData[ik+sideshift][0]>0 && type == 7  )  numPmtC++;
551           }
552           else {
553             if(allData[ik+45][0]>0 && type == 7 ) numPmtA++;
554             sideshift=45;
555             sideshiftqtc=33;
556
557           }
558           Int_t nhitsPMT=0;
559           
560           for (Int_t iHt=0; iHt<5; iHt++){
561             //cfd
562             if(allData[ik+sideshift][iHt]>0) {
563               GetRawsData(shift+ik+1) -> Fill(allData[ik+sideshift][iHt]);
564               GetRawsData(210+shift)->Fill(ik+1, allData[ik+sideshift][iHt]);
565               if(type == 8 ) feffC[ik]++;
566               AliDebug(50,Form("%i CFD %i  data %s",ik, ik+sideshift,  GetRawsData(shift+ik+1)->GetName()));
567               if(type == 7  ) {
568                 nhitsPMT++;
569                 feffPhysC[ik]++;
570               }
571               
572             }
573             //led
574             if(allData[ik+12+sideshift][iHt] > 0) { 
575               GetRawsData(shift+ik+24+1)->  Fill(allData[ik+12+sideshift][iHt]);
576               GetRawsData(211+shift)->Fill(ik+1, allData[ik+12+sideshift][iHt]);
577               AliDebug(50,Form("%i LED %i  data %s",ik, ik+12+sideshift,  GetRawsData(shift+ik+1+24)->GetName()));
578               if(type == 8  ) {
579                 feffA[ik]++;
580               }
581             }
582             //led -cfd
583             
584             if(allData[ik+12+sideshift][iHt] > 0 && allData[ik+sideshift][iHt] >0 )
585               GetRawsData(shift+ik+48+1)->
586                 Fill(allData[ik+12+sideshift][iHt]-allData[ik+sideshift][iHt]);
587             
588             //qtc
589             if(allData[2*ik+sideshiftqtc+24][iHt] > 0 &&
590                allData[2*ik+sideshiftqtc+25][iHt] > 0) {
591               GetRawsData(shift+ik+72+1)->
592                 Fill(allData[2*ik+sideshiftqtc+24][iHt]-allData[2*ik+sideshiftqtc+25][iHt]);
593               GetRawsData(212+shift)->Fill(ik+1, allData[2*ik+sideshiftqtc+24][iHt]-allData[2*ik+sideshiftqtc+25][iHt]);
594               if(type == 8) feffqtc[ik]++;
595               AliDebug(50,Form("%i QTC %i  data %s",ik, 2*ik+sideshiftqtc+24,  GetRawsData(shift+ik+1+72)->GetName()));
596
597             }
598                 if(allData[2*ik+sideshiftqtc+24][iHt] > 0) {
599                   AliDebug(50,Form("%i QT0 %i  data %s",ik, 2*ik+sideshiftqtc+24,  GetRawsData(shift+ik+1+96)->GetName()));
600               GetRawsData(shift+ik+96+1)->Fill(allData[2*ik+sideshiftqtc+24][iHt]);
601             }
602             if(allData[2*ik+sideshiftqtc+25][iHt] > 0) {
603               AliDebug(50,Form("%i QT0 %i  data %s",ik, 2*ik+sideshiftqtc+25,  GetRawsData(shift+ik+1+120)->GetName()));
604               GetRawsData(shift+ik+120+1)->Fill(allData[2*ik+sideshiftqtc+25][iHt]);
605             }
606           }
607         
608           if(type == 7  ) {
609             GetRawsData(ik+176)->Fill(nhitsPMT);
610             GetRawsData(213)->Fill(numPmtC);
611             GetRawsData(214)->Fill(numPmtA);
612           }  
613         }   
614        
615       Int_t trChannel[6] = {49,50,51,52,55,56};  
616       Float_t ch2cm = 24.4*0.029979;     
617       Int_t nhitsOrA=0;
618       Int_t nhitsOrC=0;
619       for (Int_t iHt=0; iHt<5; iHt++) {
620         
621         //orA-orC phys tvdc 1 
622         if((allData[51][iHt]>0 && allData[52][iHt]>0) && allData[50][iHt]>0)
623           {
624             AliDebug(10,Form("orA-orC phys tvdc 1  %i  data %s", 217+shift,  GetRawsData(shift+217)->GetName()));
625
626             GetRawsData(217+shift)->Fill((allData[52][iHt]-allData[51][iHt])*ch2cm);
627             //        GetRawsData(345) ->Fill((allData[51][iHt]+allData[52][iHt])/2.);
628           }
629         //orA-orC phys tvdc 0 
630         if((allData[51][iHt]>0 && allData[52][iHt]>0) && allData[50][iHt]<=0)
631           {
632             AliDebug(10,Form("orA-orC phys tvdc 0  %i  data %s", 218+shift,  GetRawsData(shift+218)->GetName()));
633
634             GetRawsData(218+shift)->Fill((allData[52][iHt]-allData[51][iHt])*ch2cm);
635           }
636         if(allData[51][iHt]>0 && allData[52][iHt]>0) {
637           AliDebug(50,Form("orA-orC phys tvdc all  %i  data %s", 219+shift,  GetRawsData(shift+219)->GetName()));
638             GetRawsData(219+shift)->Fill((allData[52][iHt]-allData[51][iHt])*ch2cm);
639         }
640         for (Int_t itr=0; itr<6; itr++) {
641           if (allData[trChannel[itr]][iHt] >0) {
642              if(type == 7  )fNumTriggers[itr]++;
643              if(type == 8  )fNumTriggersCal[itr]++;
644              AliDebug(50,Form(" triggers %i  data %s", 170+itr+shift,  GetRawsData(170+itr+shift)->GetName()));
645
646             GetRawsData(170+itr+shift)->Fill(allData[trChannel[itr]][iHt]);
647           }
648         }
649             
650         if(type == 7) if(allData[51][iHt] >0) nhitsOrA++;
651         if(type == 7)if(allData[52][iHt] >0) nhitsOrC++;
652         
653         //mult trigger signals phys
654         //C side
655         if(allData[53][iHt]>0 && allData[54][iHt]>0) {
656           AliDebug(50,Form(" mpdA %i  data %s", 201+shift,  GetRawsData(201+shift)->GetName()));
657
658           GetRawsData(201+shift)->Fill(allData[53][iHt]-allData[54][iHt]);
659           if(allData[56][iHt]>0) GetRawsData(202+shift)->Fill(allData[53][iHt]-allData[54][iHt]);
660           if(allData[55][iHt]>0) GetRawsData(203+shift)->Fill(allData[53][iHt]-allData[54][iHt]);
661         }
662         
663         //A side 
664         if(allData[105][iHt]>0 && allData[106][iHt]>0) {
665           AliDebug(50,Form(" mpdC %i  data %s", 204+shift,  GetRawsData(204+shift)->GetName()));
666
667              GetRawsData(204+shift)->Fill(allData[105][iHt]-allData[106][iHt]);
668              if(allData[56][iHt]>0) GetRawsData(205+shift)->Fill(allData[105][iHt]-allData[106][iHt]);
669              if(allData[55][iHt]>0) GetRawsData(206+shift)->Fill(allData[105][iHt]-allData[106][iHt]);
670         }
671       }
672       
673       GetRawsData(215)->Fill(nhitsOrA);
674       GetRawsData(216)->Fill(nhitsOrC);
675
676       //draw satellite
677       Int_t besttimeA=9999999;
678       Int_t besttimeC=9999999;
679       if(type == 7){    
680         if( fnEventPhys > 2000) {
681           for (Int_t ipmt=0; ipmt<12; ipmt++){
682             if(allData[ipmt+1][0] > 1 ) {
683               time[ipmt] = allData[ipmt+1][0] - Int_t(GetRawsData(ipmt+1)->GetMean());
684               if(time[ipmt]<besttimeC)
685                 besttimeC=time[ipmt]; //timeC
686             }
687           }
688           for ( Int_t ipmt=12; ipmt<24; ipmt++){
689             if(allData[ipmt+45][0] > 0) {
690               time[ipmt] = allData[ipmt+45][0] - Int_t(GetRawsData(ipmt+1)->GetMean());
691               if(time[ipmt]<besttimeA) 
692                 besttimeA=time[ipmt]; //timeA
693             }
694           }
695           
696           if(besttimeA<99999 &&besttimeC< 99999) {
697             Float_t t0 =  24.4 * (Float_t( besttimeA+besttimeC)/2. );
698             Float_t ver = 24.4 * Float_t( besttimeA-besttimeC)/2.;
699             GetRawsData(220)->Fill(0.001*ver, 0.001*(t0));
700             if(allData[50][0] > 0)  GetRawsData(221)->Fill(0.001*ver, 0.001*(t0));
701             if(allData[50][0] <= 0) GetRawsData(222)->Fill(0.001*ver, 0.001*(t0));
702           }
703         } //event >100
704       } //type 7
705     } //next
706   
707   
708
709   delete start;
710 }
711
712
713
714
715 //____________________________________________________________________________
716 void AliT0QADataMakerRec::MakeDigits( TTree *digitsTree)
717 {
718   //fills QA histos for Digits
719
720   TArrayI *digCFD = new TArrayI(24);
721   TArrayI *digLED = new TArrayI(24);
722   TArrayI *digQT0 = new TArrayI(24);
723   TArrayI *digQT1 = new TArrayI(24);
724   Int_t refpoint=0;
725   
726   TBranch *brDigits=digitsTree->GetBranch("T0");
727   AliT0digit *fDigits = new AliT0digit() ;
728   if (brDigits) {
729     brDigits->SetAddress(&fDigits);
730   }else{
731     AliError(Form("EXEC Branch T0 digits not found"));
732     return;
733   }
734   digitsTree->GetEvent(0);
735   digitsTree->GetEntry(0);
736   brDigits->GetEntry(0);
737   fDigits->GetTimeCFD(*digCFD);
738   fDigits->GetTimeLED(*digLED);
739   fDigits->GetQT0(*digQT0);
740   fDigits->GetQT1(*digQT1);
741   refpoint = fDigits->RefPoint();
742   for (Int_t i=0; i<24; i++)
743     {
744     if (digCFD->At(i)>0) {
745       Int_t cfd=digCFD->At(i)- refpoint;
746       GetDigitsData(0) ->Fill(i,cfd);
747       GetDigitsData(1) -> Fill(i, (digLED->At(i) - digCFD->At(i)));
748       GetDigitsData(2) -> Fill(i, (digQT1->At(i) - digQT0->At(i)));
749     }
750     }  
751   
752   delete digCFD;
753   delete digLED;
754   delete digQT0;
755   delete digQT1;
756   
757 }
758
759 //____________________________________________________________________________
760 void AliT0QADataMakerRec::MakeRecPoints(TTree * clustersTree)
761 {
762   //fills QA histos for clusters
763
764   AliT0RecPoint* frecpoints= new AliT0RecPoint ();
765   if (!frecpoints) {
766     AliError(":MakeRecPoints >> no recpoints found");
767     return;
768   }
769   TBranch *brRec =clustersTree ->GetBranch("T0");
770   if (brRec) {
771     brRec->SetAddress(&frecpoints);
772   }else{
773     AliError(Form("EXEC Branch T0 rec not found "));
774     return;
775   } 
776
777   brRec->GetEntry(0);
778   
779   for ( Int_t i=0; i<24; i++) {
780     if(i<12)
781       GetRecPointsData(0) -> Fill(i, frecpoints -> GetTime(i) - frecpoints -> GetTime(0)); 
782     if(i>11)
783       GetRecPointsData(0) -> Fill(i,  frecpoints -> GetTime(i) - frecpoints -> GetTime(12)); 
784     GetRecPointsData(1) -> Fill( i, frecpoints -> GetAmp(i) - frecpoints->AmpLED(i));
785   }
786   Double_t mmm=frecpoints->GetOnlineMean()- frecpoints->GetMeanTime();
787   GetRecPointsData(2) ->Fill(mmm);
788   
789 }
790
791 //____________________________________________________________________________
792 void AliT0QADataMakerRec::MakeESDs(AliESDEvent * esd)
793 {
794   //fills QA histos for ESD
795   
796   const Double32_t  *mean;
797   mean = esd->GetT0TOF();
798   Double32_t t0time= 0.001*mean[0];
799   Double32_t orA= 0.001*mean[1];
800   Double32_t orC=0.001* mean[2];
801
802   if (t0time<99)   GetESDsData(0) -> Fill(t0time);
803   if( esd->GetT0zVertex() <99) GetESDsData(1)-> Fill(esd->GetT0zVertex());
804   if( orA<99 && orC<99) GetESDsData(2)-> Fill((orA-orC)/2.);
805   
806 }
807 //____________________________________________________________________________
808
809 /*
810 void AliT0QADataMakerRec::GetMeanAndSigma(TH1F* hist, Float_t &mean, Float_t &sigma) 
811 {
812
813   const double window = 3.;  //fit window 
814  
815   double meanEstimate, sigmaEstimate; 
816   int maxBin;
817   maxBin        =  hist->GetMaximumBin(); //position of maximum
818   meanEstimate  =  hist->GetBinCenter( maxBin); // mean of gaussian sitting in maximum
819   sigmaEstimate = hist->GetRMS();
820   TF1* fit= new TF1("fit","gaus", meanEstimate - window*sigmaEstimate, meanEstimate + window*sigmaEstimate);
821   fit->SetParameters(hist->GetBinContent(maxBin), meanEstimate, sigmaEstimate);
822   hist->Fit("fit","RQ","Q");
823
824   mean  = (Float_t) fit->GetParameter(1);
825   sigma = (Float_t) fit->GetParameter(2);
826
827   delete fit;
828 }
829 */