]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCQADataMakerRec.cxx
Changes in QA to be able to process separately different triggers (Ruben)
[u/mrichter/AliRoot.git] / ZDC / AliZDCQADataMakerRec.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 ZDC quality assurance. //
19 //  QA objects are 1 & 2 Dimensional histograms.                     //
20 //  author: C. Oppedisano                                            //
21 //                                                                   //
22 ///////////////////////////////////////////////////////////////////////
23
24 // --- ROOT system ---
25 #include <TClonesArray.h>
26 #include <TFile.h> 
27 #include <TH1F.h> 
28 #include <TH2F.h>
29 #include <TLine.h>
30 #include <TProfile.h>
31 #include <Riostream.h>
32 // --- Standard library ---
33
34 // --- AliRoot header files ---
35 #include "AliLog.h"
36 #include "AliQAChecker.h"
37 #include "AliZDCReco.h"
38 #include "AliRawReader.h"
39 #include "AliZDCQADataMakerRec.h"
40 #include "AliZDCPedestals.h"
41 #include "AliZDCRawStream.h"
42 #include "AliZDCDigit.h"
43 #include "AliESDZDC.h"
44 #include "AliESDEvent.h"
45
46 ClassImp(AliZDCQADataMakerRec)
47            
48 //____________________________________________________________________________ 
49   AliZDCQADataMakerRec::AliZDCQADataMakerRec() : 
50   AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kZDC), "ZDC Quality Assurance Data Maker"), 
51   fPedCalibData(0x0)
52 {
53   // ctor
54 }
55
56 //____________________________________________________________________________ 
57 AliZDCQADataMakerRec::AliZDCQADataMakerRec(const AliZDCQADataMakerRec& qadm) :
58   AliQADataMakerRec(),      
59   fPedCalibData(qadm.fPedCalibData)
60
61 {
62   //copy ctor 
63   SetName((const char*)qadm.GetName()); 
64   SetTitle((const char*)qadm.GetTitle()); 
65 }
66
67 //__________________________________________________________________
68 AliZDCQADataMakerRec& AliZDCQADataMakerRec::operator = (const AliZDCQADataMakerRec& qadm )
69 {
70   // Equal operator.
71   this->~AliZDCQADataMakerRec();
72   new(this) AliZDCQADataMakerRec(qadm);
73   return *this;
74 }
75
76 //____________________________________________________________________________ 
77 AliZDCQADataMakerRec::~AliZDCQADataMakerRec()
78 {
79   if(fPedCalibData && !(AliCDBManager::Instance()->GetCacheFlag())){
80     delete fPedCalibData;
81     fPedCalibData=0;
82   } 
83 }
84
85 //____________________________________________________________________________ 
86 AliZDCPedestals* AliZDCQADataMakerRec::GetPedCalibData() const
87 {
88
89   // Retrieving pedestal calibration object from OCDB
90   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
91   if(!entry) AliWarning("No calibration data loaded!");  
92
93   AliZDCPedestals *calibdata = (AliZDCPedestals*)  (entry->GetObject());
94   if(!calibdata) AliFatal("Wrong calibration object in calibration  file!");
95
96   return calibdata;
97
98 }
99
100 //____________________________________________________________________________ 
101 void AliZDCQADataMakerRec::InitDigits()
102 {
103   // create Digits histograms in Digits subdir
104   //
105   const Bool_t expert   = kTRUE ; 
106   const Bool_t image    = kTRUE ; 
107   
108   // ------------------- HIGH GAIN CHAIN ---------------------------
109   TH1F * hDigZNCTot = new TH1F("hDigZNCTot", "Signal in ZNC;Amplitude [ADC counts];Counts", 100, 0., 6000.);
110   TH1F * hDigZNATot = new TH1F("hDigZNATot", "Signal in ZNA;Amplitude [ADC counts];Counts", 100, 0., 6000.);
111   TH1F * hDigZPCTot = new TH1F("hDigZPCTot", "Signal in ZPC;Amplitude [ADC counts];Counts", 100, 0., 6000.);
112   TH1F * hDigZPATot = new TH1F("hDigZPATot", "Signal in ZPA;Amplitude [ADC counts];Counts", 100, 0., 6000.);
113   Add2DigitsList(hDigZNCTot, 0, !expert, image);
114   Add2DigitsList(hDigZNATot, 1, !expert, image);
115   Add2DigitsList(hDigZPCTot, 2, !expert, image);
116   Add2DigitsList(hDigZPATot, 3, !expert, image);
117   //
118   TH1F * hDigSumQZNC = new TH1F("hDigSumQZNC", "Signal in 4 ZNC PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
119   TH1F * hDigSumQZNA = new TH1F("hDigSumQZNA", "Signal in 4 ZNA PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
120   TH1F * hDigSumQZPC = new TH1F("hDigSumQZPC", "Signal in 4 ZPC PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
121   TH1F * hDigSumQZPA = new TH1F("hDigSumQZPA", "Signal in 4 ZPA PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
122   Add2DigitsList(hDigSumQZNC, 4, expert, !image);
123   Add2DigitsList(hDigSumQZNA, 5, expert, !image);
124   Add2DigitsList(hDigSumQZPC, 6, expert, !image);
125   Add2DigitsList(hDigSumQZPA, 7, expert, !image);
126   //
127   TH1F * hDigPMCZNC = new TH1F("hDigPMCZNC", "Signal in ZNC PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
128   TH1F * hDigPMCZNA = new TH1F("hDigPMCZNA", "Signal in ZNA PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
129   TH1F * hDigPMCZPC = new TH1F("hDigPMCZPC", "Signal in ZPC PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
130   TH1F * hDigPMCZPA = new TH1F("hDigPMCZPA", "Signal in ZPA PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
131   Add2DigitsList(hDigPMCZNC, 8, expert, !image);
132   Add2DigitsList(hDigPMCZNA, 9, expert, !image);
133   Add2DigitsList(hDigPMCZPC, 10, expert, !image);
134   Add2DigitsList(hDigPMCZPA, 11, expert, !image);
135   // 
136   ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
137 }
138
139 //____________________________________________________________________________
140 void AliZDCQADataMakerRec::InitRaws()
141 {
142   // create Digits histograms in Digits subdir
143   const Bool_t expert   = kTRUE ; 
144   const Bool_t image    = kTRUE ; 
145
146   TH1F * hZNCSpectrum = new TH1F("hZNCSpectrum","ZNC spectrum;Amplitude [ADC counts];Counts",100,0.,1200.);
147   TH1F * hZPCSpectrum = new TH1F("hZPCSpectrum","ZPC spectrum;Amplitude [ADC counts];Counts",100,0.,1200.);
148   TH1F * hZNASpectrum = new TH1F("hZNASpectrum","ZNA spectrum;Amplitude [ADC counts];Counts",100,0.,1200.);
149   TH1F * hZPASpectrum = new TH1F("hZPASpectrum","ZPA spectrum;Amplitude [ADC counts];Counts",100,0.,1200.);
150   TH1F * hZEM1Spectrum = new TH1F("hZEM1Spectrum","ZEM1 spectrum;Amplitude [ADC counts];Counts",100,8., 1208.);
151   TH1F * hZEM2Spectrum = new TH1F("hZEM2Spectrum","ZEM2 spectrum;Amplitude [ADC counts];Counts",100,8., 1208.);
152   Add2RawsList(hZNCSpectrum, 0, expert, !image);
153   Add2RawsList(hZNASpectrum, 1, expert, !image);
154   Add2RawsList(hZPCSpectrum, 2, expert, !image);
155   Add2RawsList(hZPASpectrum, 3, expert, !image);
156   Add2RawsList(hZEM1Spectrum, 4, expert, !image);
157   Add2RawsList(hZEM2Spectrum, 5, expert, !image);
158     
159   TH1F * hRawPMCZNC = new TH1F("hRawPMCZNC", "Raw ZNC PMC;Amplitude [ADC counts];Counts",100, 8., 1208.);
160   TH1F * hRawPMCZNA = new TH1F("hRawPMCZNA", "Raw ZNA PMC;Amplitude [ADC counts];Counts",100, 8., 1208.);
161   TH1F * hRawPMCZPC = new TH1F("hRawPMCZPC", "Raw ZPC PMC;Amplitude [ADC counts];Counts",100, 8., 1208.);
162   TH1F * hRawPMCZPA = new TH1F("hRawPMCZPA", "Raw ZPA PMC;Amplitude [ADC counts];Counts",100, 8., 1208.);
163   Add2RawsList(hRawPMCZNC, 6, expert, !image);
164   Add2RawsList(hRawPMCZNA, 7, expert, !image);
165   Add2RawsList(hRawPMCZPC, 8, expert, !image);
166   Add2RawsList(hRawPMCZPA, 9, expert, !image);
167   TH1F * hRawSumQZNC = new TH1F("hRawSumQZNC", "Raw sumQ ZNC;Amplitude [ADC counts];Counts",100, 0., 1200.);
168   TH1F * hRawSumQZNA = new TH1F("hRawSumQZNA", "Raw sumQ ZNA;Amplitude [ADC counts];Counts",100, 0., 1200.);
169   TH1F * hRawSumQZPC = new TH1F("hRawSumQZPC", "Raw sumQ ZPC;Amplitude [ADC counts];Counts",100, 0., 1200.);
170   TH1F * hRawSumQZPA = new TH1F("hRawSumQZPA", "Raw sumQ ZPA;Amplitude [ADC counts];Counts",100, 0., 1200.);
171   Add2RawsList(hRawSumQZNC, 10, expert, !image);
172   Add2RawsList(hRawSumQZNA, 11, expert, !image);
173   Add2RawsList(hRawSumQZPC, 12, expert, !image);
174   Add2RawsList(hRawSumQZPA, 13, expert, !image);
175   
176   TH1F * hRawTDCZEM1 = new TH1F("hRawTDCZEM1", "Raw TDC ZEM1;TDC [ns]",160, -350., -310.);
177   Add2RawsList(hRawTDCZEM1, 14, expert, !image);
178   TH1F * hRawTDCZPC = new TH1F("hRawTDCZPC", "Raw TDC ZPC;TDC [ns]",160, -350., -310.);
179   Add2RawsList(hRawTDCZPC, 15, expert, !image);
180   
181   TProfile * hRawADCProfs = new TProfile("hRawADCProfs", "ADC profiles;ADC id;Mean ADC values",22,-0.5,21.5,10.,1210.,"");
182   Add2RawsList(hRawADCProfs, 16, expert, !image);
183   TProfile * hRawTDCProfs = new TProfile("hRawTDCProfs", "TDC profiles;TDC id;Mean TDC values",6,0.5,6.5,-340.,-300.,"S");
184   Add2RawsList(hRawTDCProfs, 17, expert, !image);
185   
186   TH1F * hRawADCs = new TH1F("hRawADCs", "ADCs;ADC id;Mean ADC values",22,-0.5,21.5);
187   Add2RawsList(hRawADCs, 18, !expert, image);
188  
189   TH1F * hRawTDCs = new TH1F("hRawTDCs", "TDCs;TDC id;Mean TDC values",6,0.5,6.5);
190 //  hRawTDCs->SetMaximum(-300); hRawTDCs->SetMinimum(-340);
191   Add2RawsList(hRawTDCs, 19, !expert, image);
192   
193   TH2F *hZNCrawCentr  = new TH2F("hZNCrawCentr", "Centroid in ZNC;X (cm);Y(cm)", 100, -5.,5.,100,-5.,5.);
194   Add2RawsList(hZNCrawCentr, 20, expert, !image);
195   TH2F *hZNArawCentr  = new TH2F("hZNArawCentr", "Centroid in ZNA;X (cm);Y(cm)", 100, -5.,5.,100,-5.,5.);
196   Add2RawsList(hZNArawCentr, 21, expert, !image);
197   
198   TH2F *hTimeZDC = new TH2F("fTimeZDC", "ZDC;TDC timing A+C vs C-A; events", 20,-30,30,120,-100,-40);
199   Add2RawsList(hTimeZDC, 22, expert, !image);
200   //
201   ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
202 }
203
204 //____________________________________________________________________________
205 void AliZDCQADataMakerRec::InitRecPoints()
206 {
207   // create Digits histograms in Digits subdir
208   const Bool_t expert = kTRUE ; 
209   const Bool_t image  = kTRUE ; 
210   //
211   // ------------------- HIGH GAIN CHAIN ---------------------------
212   TH1F * hRecZNCTot = new TH1F("hRecZNCTot", "Rec signal in ZNC;Amplitude [ADC counts];Counts", 100, 0., 2000.);
213   TH1F * hRecZNATot = new TH1F("hRecZNATot", "Rec signal in ZNA;Amplitude [ADC counts];Counts", 100, 0., 2000.);
214   TH1F * hRecZPCTot = new TH1F("hRecZPCTot", "Rec signal in ZPC;Amplitude [ADC counts];Counts", 100, 0., 8000.);
215   TH1F * hRecZPATot = new TH1F("hRecZPATot", "Rec signal in ZPA;Amplitude [ADC counts];Counts", 100, 0., 8000.);
216   Add2RecPointsList(hRecZNCTot, 0, !expert, image);
217   Add2RecPointsList(hRecZNATot, 1, !expert, image);
218   Add2RecPointsList(hRecZPCTot, 2, !expert, image);
219   Add2RecPointsList(hRecZPATot, 3, !expert, image);
220   //
221   TH1F * hRecSumQZNC = new TH1F("hRecSumQZNC", "Rec summed 4 ZNC quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
222   TH1F * hRecSumQZNA = new TH1F("hRecSumQZNA", "Rec summed 4 ZNA quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
223   TH1F * hRecSumQZPC = new TH1F("hRecSumQZPC", "Rec summed 4 ZPC quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
224   TH1F * hRecSumQZPA = new TH1F("hRecSumQZPA", "Rec summed 4 ZPA quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
225   Add2RecPointsList(hRecSumQZNC, 4, expert, !image);
226   Add2RecPointsList(hRecSumQZNA, 5, expert, !image);
227   Add2RecPointsList(hRecSumQZPC, 6, expert, !image);
228   Add2RecPointsList(hRecSumQZPA, 7, expert, !image);
229   //
230   TH1F * hRecPMCZNC = new TH1F("hRecPMCZNC", "Rec common ZNC PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
231   TH1F * hRecPMCZNA = new TH1F("hRecPMCZNA", "Rec common ZNA PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
232   TH1F * hRecPMCZPC = new TH1F("hRecPMCZPC", "Rec common ZPC PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
233   TH1F * hRecPMCZPA = new TH1F("hRecPMCZPA", "Rec common ZPA PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
234   Add2RecPointsList(hRecPMCZNC, 8 , expert, !image);
235   Add2RecPointsList(hRecPMCZNA, 9 , expert, !image);
236   Add2RecPointsList(hRecPMCZPC, 10, expert, !image);
237   Add2RecPointsList(hRecPMCZPA, 11, expert, !image); 
238   //
239   ClonePerTrigClass(AliQAv1::kRECPOINTS); // this should be the last line
240 }
241
242 //____________________________________________________________________________
243 void AliZDCQADataMakerRec::InitESDs()
244 {
245   //Booking ESDs histograms
246   //
247   const Bool_t expert = kTRUE ; 
248   const Bool_t image  = kTRUE ; 
249   
250   // ------------------- HIGH GAIN CHAIN ---------------------------
251   TH1F * hESDZNCTot = new TH1F("hESDZNCTot", "Energy in ZNC", 100, 0., 4000.);
252   TH1F * hESDZNATot = new TH1F("hESDZNATot", "Energy in ZNA", 100, 0., 4000.);
253   TH1F * hESDZPCTot = new TH1F("hESDZPCTot", "Energy in ZPC", 100, 0., 4000.);
254   TH1F * hESDZPATot = new TH1F("hESDZPATot", "Energy in ZPA", 100, 0., 4000.);
255   Add2ESDsList(hESDZNCTot, 0, !expert, image);
256   Add2ESDsList(hESDZNATot, 1, !expert, image);
257   Add2ESDsList(hESDZPCTot, 2, !expert, image);
258   Add2ESDsList(hESDZPATot, 3, !expert, image);
259   //
260   TH1F * hESDZEM1 = new TH1F("hESDZEM1", "Energy in ZEM1", 100, 0., 2000.);
261   TH1F * hESDZEM2 = new TH1F("hESDZEM2", "Energy in ZEM2", 100, 0., 2000.);
262   Add2ESDsList(hESDZEM1,4, !expert, image);
263   Add2ESDsList(hESDZEM2,5, !expert, image);
264   //
265   TH1F * hESDSumQZNC = new TH1F("hESDSumQZNC", "Sum of 4 ZNC energy",100, 0., 2000.);
266   TH1F * hESDSumQZNA = new TH1F("hESDSumQZNA", "Sum of 4 ZNA energy",100, 0., 2000.);
267   TH1F * hESDSumQZPC = new TH1F("hESDSumQZPC", "Sum of 4 ZPC energy",100, 0., 2000.);
268   TH1F * hESDSumQZPA = new TH1F("hESDSumQZPA", "Sum of 4 ZPA energy",100, 0., 2000.);
269   Add2ESDsList(hESDSumQZNC, 6, expert, !image);
270   Add2ESDsList(hESDSumQZNA, 7, expert, !image);
271   Add2ESDsList(hESDSumQZPC, 8, expert, !image);
272   Add2ESDsList(hESDSumQZPA, 9, expert, !image);
273   //
274   TH1F * hESDPMCZNC = new TH1F("hESDPMCZNC", "Energy in ZNC PMC",100, 0., 2000.);
275   TH1F * hESDPMCZNA = new TH1F("hESDPMCZNA", "Energy in ZNA PMC",100, 0., 2000.);
276   TH1F * hESDPMCZPC = new TH1F("hESDPMCZPC", "Energy in ZPC PMC",100, 0., 2000.);
277   TH1F * hESDPMCZPA = new TH1F("hESDPMCZPA", "Energy in ZPA PMC",100, 0., 2000.);
278   Add2ESDsList(hESDPMCZNC, 10, expert, !image);
279   Add2ESDsList(hESDPMCZNA, 11, expert, !image);
280   Add2ESDsList(hESDPMCZPC, 12, expert, !image);
281   Add2ESDsList(hESDPMCZPA, 13, expert, !image);
282   // 
283   ClonePerTrigClass(AliQAv1::kESDS); // this should be the last line
284 }
285
286 //___________________________________________________________________________
287 void AliZDCQADataMakerRec::MakeDigits(TTree *digitTree)
288 {
289   // makes data from Digit Tree
290   if(!GetDigitsData(0)) InitDigits();
291
292   if(!digitTree){
293     AliError("Can't get ZDC digit tree!!");
294     return; 
295   }     
296    
297   TBranch * branch = digitTree->GetBranch("ZDC");
298   if(!branch){
299     AliError("ZDC branch in digit tree not found"); 
300     return;
301   } 
302     
303   AliZDCDigit *digit = 0x0;
304   branch->SetAddress(&digit);
305      
306   Float_t adcSum_ZNC=0., adcSum_ZNA=0., adcSum_ZPC=0., adcSum_ZPA=0.;
307   Float_t adcSumQ_ZNC=0., adcSumQ_ZNA=0., adcSumQ_ZPC=0., adcSumQ_ZPA=0.;
308   //Float_t adcSum_ZNC_lg=0., adcSum_ZNA_lg=0., adcSum_ZPC_lg=0., adcSum_ZPA_lg=0.;
309   //Float_t adcSumQ_ZNC_lg=0., adcSumQ_ZNA_lg=0., adcSumQ_ZPC_lg=0., adcSumQ_ZPA_lg=0.;
310   
311   Int_t ndig = digitTree->GetEntries();
312   for(Int_t i=0; i<ndig; i++){
313       branch->GetEntry(i);
314       
315       if(digit->GetSector(0)==1 && digit->GetSector(1)!=5){
316           adcSum_ZNC += digit->GetADCValue(0);
317           //adcSum_ZNC_lg += digit->GetADCValue(1);
318           //
319           if(digit->GetSector(1)!=0){
320               adcSumQ_ZNC += digit->GetADCValue(0);
321               //adcSumQ_ZNC_lg+= digit->GetADCValue(1);
322           }
323           else{
324               FillDigitsData(8,digit->GetADCValue(0));
325               //FillDigitsData(20,digit->GetADCValue(1));
326           }
327       }
328       else if(digit->GetSector(0)==2){
329           adcSum_ZPC += digit->GetADCValue(0);
330           //adcSum_ZPC_lg += digit->GetADCValue(1);
331           //
332           if(digit->GetSector(1)!=0){
333               adcSumQ_ZPC += digit->GetADCValue(0);
334               //adcSumQ_ZPC_lg+= digit->GetADCValue(1);
335           }
336           else{
337               FillDigitsData(10,digit->GetADCValue(0));
338               //FillDigitsData(22,digit->GetADCValue(1));
339           }
340       }
341       else if(digit->GetSector(0)==4 && digit->GetSector(1)!=5){
342           adcSum_ZNA += digit->GetADCValue(0);
343           //adcSum_ZNA_lg += digit->GetADCValue(1);
344           //
345           if(digit->GetSector(1)!=0){
346               adcSumQ_ZNA += digit->GetADCValue(0);
347               //adcSumQ_ZNA_lg+= digit->GetADCValue(1);
348           }
349           else{
350               FillDigitsData(9,digit->GetADCValue(0));
351               //FillDigitsData(21,digit->GetADCValue(1));
352           }
353       }
354       else if(digit->GetSector(0)==5){
355           adcSum_ZPA += digit->GetADCValue(0);
356           //adcSum_ZPA_lg += digit->GetADCValue(1);
357           //
358           if(digit->GetSector(1)!=0){
359               adcSumQ_ZPA += digit->GetADCValue(0);
360               //adcSumQ_ZPA_lg+= digit->GetADCValue(1);
361           }
362           else{
363               FillDigitsData(11,digit->GetADCValue(0));
364               //FillDigitsData(23,digit->GetADCValue(1));
365           }
366       }
367   }
368   //
369   FillDigitsData(0,adcSum_ZNC);
370   FillDigitsData(1,adcSum_ZNA);
371   FillDigitsData(2,adcSum_ZPC);
372   FillDigitsData(3,adcSum_ZPA);
373   //
374   FillDigitsData(4,adcSumQ_ZNC);
375   FillDigitsData(5,adcSumQ_ZNA);
376   FillDigitsData(6,adcSumQ_ZPC);
377   FillDigitsData(7,adcSumQ_ZPA);
378   
379   delete digit;
380   digit=0;
381   //
382   IncEvCountCycleDigits();
383   IncEvCountTotalDigits();
384   //
385 }
386
387
388 //____________________________________________________________________________
389 void AliZDCQADataMakerRec::MakeRaws(AliRawReader *rawReader)
390 {
391   // Filling Raws QA histos
392   //
393   // Checking the event type 
394 //  if (rawReader->GetType()!=7){
395   
396     // Check if histograms already created for this Event Specie
397     if(!GetRawsData(0)) InitRaws();
398     // Parameters for mean value pedestal subtraction
399     int const kNch = 24;
400     Float_t meanPed[2*kNch];    
401     for(Int_t jj=0; jj<2*kNch; jj++) meanPed[jj] = fPedCalibData->GetMeanPed(jj);
402     
403     Float_t zncSignal=0., znaSignal=0., zpcSignal=0., zpaSignal=0.;
404     Float_t zncSumQ=0., znaSumQ=0., zpcSumQ=0., zpaSumQ=0.;
405     Float_t zncpmC=0., znapmC=0., zpcpmC=0., zpapmC=0.;
406     Bool_t isZNCFired=kFALSE, isZPCFired=kFALSE, isZNAFired=kFALSE, isZPAFired=kFALSE;
407     Int_t  indZNC=0, indZNA=0, indZPC=0, indZPA=0;
408     Float_t zncTDC[10], zpcTDC[10], zem1TDC[10], zem2TDC[10], znaTDC[10], zpaTDC[10];
409     Float_t zncSumTDC[10], znaSumTDC[10];
410     for(Int_t i=0; i<10; i++){
411        zncTDC[i]=zpcTDC[i]=zem1TDC[i]=zem2TDC[i]=znaTDC[i]=zpaTDC[i]=zncSumTDC[i]=znaSumTDC[i]=-999.;
412     }
413     Float_t tdcGate=-999., l0=-999.;
414     Int_t iMultZNCTDC=0, iMultZPCTDC=0, iMultZEM1TDC=0, iMultZEM2TDC=0, iMultZNATDC=0, iMultZPATDC=0;
415     Int_t iMultTDCC=0, iMultTDCA=0;
416     
417     const Float_t x[4] = {-1.75, 1.75, -1.75, 1.75};
418     const Float_t y[4] = {-1.75, -1.75, 1.75, 1.75};
419     const Float_t alpha=0.5;
420     Float_t numXZNC=0., numYZNC=0., denZNC=0., wZNC=0.; 
421     Float_t numXZNA=0., numYZNA=0., denZNA=0., wZNA=0.; 
422     
423     rawReader->Reset();
424     AliZDCRawStream stream(rawReader);
425     while(stream.Next()){
426
427       if(stream.IsADCDataWord() && 
428          (stream.GetADCModule()==0 || stream.GetADCModule()==1)){
429        
430          Int_t det = stream.GetSector(0);
431          Int_t quad = stream.GetSector(1);
432          Int_t gain = stream.GetADCGain();
433          Int_t pedindex=0;
434          
435          // Stuff for pedestal subtraction
436          if(quad != 5){ // ZDCs (not reference PTMs)
437           Float_t pedSubVal=-99.;
438           if(det == 1){    
439             pedindex = quad;
440             if(gain == 0){
441               pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
442               zncSignal  += pedSubVal; 
443               isZNCFired = kTRUE;
444               if(quad!=0){
445                 zncSumQ += pedSubVal;
446                 if(pedSubVal>0.&& zncpmC>50.){
447                   wZNC = TMath::Power(pedSubVal, alpha);
448                   numXZNC += x[quad-1]*wZNC;
449                   numYZNC += y[quad-1]*wZNC;
450                   denZNC += wZNC;
451                 }
452               }
453               else{
454                 zncpmC = pedSubVal;
455                 FillRawsData(6,zncpmC);
456               }
457               indZNC++;
458               
459               FillRawsData(16, pedindex, pedSubVal);
460             }
461           }
462           else if(det == 2){ 
463             pedindex = quad+5;
464             if(gain == 0){
465               pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
466               zpcSignal += pedSubVal; 
467               isZPCFired = kTRUE;
468               if(quad!=0) zpcSumQ += pedSubVal;
469               else{
470                 zpcpmC = pedSubVal;
471                 FillRawsData(8,zpcpmC);
472               }
473               indZPC++;
474               
475               FillRawsData(16, pedindex, pedSubVal);
476             }
477           }
478           else if(det == 3){ 
479             pedindex = quad+9;
480             if(quad==1){     
481               if(gain == 0){
482                 pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
483                 FillRawsData(4,pedSubVal);
484                 FillRawsData(16,pedindex, pedSubVal);
485               }
486             }
487             else if(quad==2){ 
488               if(gain == 0){
489                 pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
490                 FillRawsData(5,pedSubVal); 
491                 FillRawsData(16,pedindex, pedSubVal);
492               }
493             }
494           }
495           else if(det == 4){       
496             pedindex = quad+12;
497             if(gain == 0){
498               pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
499               znaSignal  += pedSubVal; 
500               isZNAFired = kTRUE;
501               if(quad!=0){
502                 znaSumQ += pedSubVal;
503                 if(pedSubVal>0.&& znapmC>50.) {
504                   wZNA = TMath::Power(pedSubVal, alpha);
505                   numXZNA += x[quad-1]*wZNA;
506                   numYZNA += y[quad-1]*wZNA;
507                   denZNA += wZNA;
508                 }
509               }
510               else{
511                 znapmC = pedSubVal;
512                 FillRawsData(7,znapmC);
513               }
514               indZNA++;
515               
516               FillRawsData(16,pedindex, pedSubVal);
517               // RS: invalid
518               //GetRawsData(22)->SetBinContent(pedindex+1, GetRawsData(20)->GetBinContent(pedindex+1));
519               //GetRawsData(22)->SetBinError(pedindex+1, GetRawsData(20)->GetBinError(pedindex+1));
520             }
521           }
522           else if(det == 5){
523             pedindex = quad+17;
524             if(gain == 0){
525               pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
526               zpaSignal  += pedSubVal; 
527               isZPAFired = kTRUE;
528               if(quad!=0) zpaSumQ += pedSubVal;
529               else{
530                 zpapmC = pedSubVal;
531                 FillRawsData(9,zpapmC);
532               }
533               indZPA++;
534               
535               FillRawsData(16,pedindex, pedSubVal);
536               // RS: invalid
537               //GetRawsData(22)->SetBinContent(pedindex+1, GetRawsData(20)->GetBinContent(pedindex+1));
538               //GetRawsData(22)->SetBinError(pedindex+1, GetRawsData(20)->GetBinError(pedindex+1));
539             }
540           }
541                          
542          }
543
544          if(isZNCFired && indZNC==5){
545            FillRawsData(0,zncSignal);
546            FillRawsData(10,zncSumQ); 
547            //
548            Float_t xZNC, yZNC;           
549            if(denZNC!=0){
550              xZNC = numXZNC/denZNC;
551              yZNC = numYZNC/denZNC;
552            } 
553            else xZNC = yZNC = 999.;
554            FillRawsData(20,xZNC, yZNC);
555          }
556          if(isZPCFired && indZPC==5){
557            FillRawsData(2,zpcSignal);
558            FillRawsData(12,zpcSumQ); 
559          }
560          if(isZNAFired && indZNA==5){ 
561            FillRawsData(1,znaSignal);
562            FillRawsData(11,znaSumQ); 
563            //
564            Float_t xZNA, yZNA;
565            if(denZNA!=0){
566              xZNA = numXZNA/denZNA;
567              yZNA = numYZNA/denZNA;
568            } 
569            else xZNA = yZNA = 999.;
570            FillRawsData(21,xZNA, yZNA);
571          }
572          if(isZPAFired && indZPA==5){ 
573            FillRawsData(3,zpaSignal);
574            FillRawsData(13,zpaSumQ); 
575          }
576          
577          if(indZNC==5){
578            zncSignal = zncSumQ = zncpmC = 0;
579            isZNCFired=kFALSE; indZNC=0;
580          }
581          if(indZPC==5){
582            zpcSignal = zpcSumQ = zpcpmC = 0;
583            isZPCFired=kFALSE; indZPC=0;
584          }
585          if(indZNA==5){
586            znaSignal = znaSumQ = znapmC = 0;
587            isZNAFired=kFALSE; indZNA=0;
588          }
589          if(indZPA==5){
590            zpaSignal = zpaSumQ = zpapmC = 0;
591            isZPAFired=kFALSE; indZPA=0;
592          }
593          
594       } //IsADCDataWord && signal ADCs
595       else if(stream.IsZDCTDCDatum()){
596          if(stream.GetChannel()==1){
597             zncTDC[iMultZNCTDC] = (0.025*stream.GetZDCTDCDatum());
598             iMultZNCTDC++;
599          }
600          else if(stream.GetChannel()==3){
601             zpcTDC[iMultZPCTDC] = (0.025*stream.GetZDCTDCDatum());
602             iMultZPCTDC++;
603          }
604          else if(stream.GetChannel()==5){
605             znaTDC[iMultZNATDC] = (0.025*stream.GetZDCTDCDatum());
606             iMultZNATDC++;
607          }
608          else if(stream.GetChannel()==7){
609             zpaTDC[iMultZPATDC] = (0.025*stream.GetZDCTDCDatum());
610             iMultZPATDC++;
611          }
612          else if(stream.GetChannel()==8){
613             zem1TDC[iMultZEM1TDC] = (0.025*stream.GetZDCTDCDatum());
614             iMultZEM1TDC++;
615          }
616          else if(stream.GetChannel()==9){
617             zem2TDC[iMultZEM2TDC] = (0.025*stream.GetZDCTDCDatum());
618             iMultZEM2TDC++;
619          }
620          else if(stream.GetChannel()==10){
621             zncSumTDC[iMultZEM2TDC] = (0.025*stream.GetZDCTDCDatum());
622             iMultTDCC++;
623          }
624          else if(stream.GetChannel()==12){
625             znaSumTDC[iMultZEM2TDC] = (0.025*stream.GetZDCTDCDatum());
626             iMultTDCA++;
627          }
628          else if(stream.GetChannel()==14) tdcGate = (0.025*stream.GetZDCTDCDatum());
629          else if(stream.GetChannel()==15) l0 = (0.025*stream.GetZDCTDCDatum());
630          
631          if(stream.GetChannel()==16 && tdcGate!=-999.){
632            for(Int_t iHit=0; iHit<10; iHit++){
633               if(zncTDC[iHit]!=-999.){
634                 if(zncTDC[iHit]-tdcGate>-340. && zncTDC[iHit]-tdcGate<-300.) 
635                    FillRawsData(17,1, zncTDC[iHit]-tdcGate);
636               }
637               if(zpcTDC[iHit]!=-999.){
638                 Float_t diffZPC = zpcTDC[iHit]-tdcGate;
639                 FillRawsData(15,diffZPC);
640                 if(diffZPC>-340. && diffZPC<-300.) FillRawsData(17,2, diffZPC);
641               }
642               if(znaTDC[iHit]!=-999.){
643                 if(znaTDC[iHit]-tdcGate>-340. && znaTDC[iHit]-tdcGate<-300.) 
644                   FillRawsData(17,3, znaTDC[iHit]-tdcGate);
645               }
646               if(zpaTDC[iHit]!=-999.){
647                 if(zpaTDC[iHit]-tdcGate>-340. && zpaTDC[iHit]-tdcGate<-300.) 
648                   FillRawsData(17,4, zpaTDC[iHit]-tdcGate);
649               }
650               if(zem1TDC[iHit]!=-999.){
651                 Float_t diffZEM1 = zem1TDC[iHit]-tdcGate;
652                 FillRawsData(14,diffZEM1);
653                 if(diffZEM1>-340. && diffZEM1<-300.) FillRawsData(17,5, diffZEM1);
654               }
655               if(zem2TDC[iHit]!=-999.){
656                 if(zem2TDC[iHit]-tdcGate>-340. && zem2TDC[iHit]-tdcGate<-300.) 
657                   FillRawsData(17,6, zem2TDC[iHit]-tdcGate);
658               }
659               if(zncSumTDC[iHit]!=-999.){
660                  Float_t tdcC = zncSumTDC[iHit]-l0;
661                  if(znaSumTDC[iHit]!=-999.){
662                     Float_t tdcA = znaSumTDC[iHit]-l0;
663                     //if (((tdcC-tdcA-refDelta)*(tdcC-tdcA-refDelta)/(sigmaDelta*sigmaDelta) +
664                         //(tdcC+tdcA-refSum)*(tdcC+tdcA-refSum)/(sigmaSum*sigmaSum))< 1.0)
665                         FillRawsData(22,tdcC-tdcA,tdcC+tdcA);
666                     
667                  }
668               }
669            }
670            //
671            tdcGate = -999.;
672            for(Int_t i=0; i<10; i++){
673               zncTDC[i] = zpcTDC[i] = zem1TDC[i] = zem2TDC[i] = znaTDC[i] = zpaTDC[i] = -999.;
674               zncSumTDC[i] = znaSumTDC[i] = -999.;
675            } 
676          }
677       }
678     
679     } //stream.Next()
680 //  } // check on event type
681 //  else{
682 //    AliDebug(1,Form("Skipping non-physics event for QA -> event type %d \n", rawReader->GetType())); 
683 //  }
684 //
685     IncEvCountCycleRaws();
686     IncEvCountTotalRaws();
687     //
688 }
689
690 //____________________________________________________________________________
691 void AliZDCQADataMakerRec::MakeRecPoints(TTree * clustersTree)
692 {
693   // Filling QA histos from RecPoints
694
695   TBranch *branch = clustersTree->GetBranch("ZDC");
696   if(!branch){ 
697     AliError("Can't get the ZDC branch for rec points!");
698     return;
699   }
700   
701   if(!GetRecPointsData(0)) InitRecPoints() ;
702
703   Float_t sum_ZNC=0., sum_ZNA=0., sum_ZPC=0., sum_ZPA=0.;
704   Float_t sumQ_ZNC=0., sumQ_ZNA=0., sumQ_ZPC=0., sumQ_ZPA=0.;
705
706   AliZDCReco reco;
707   AliZDCReco* preco = &reco;
708   clustersTree->SetBranchAddress("ZDC", &preco);
709
710   clustersTree->GetEntry(0);
711   for(Int_t i=0; i<5; i++){
712     sum_ZNC += reco.GetZN1HREnTow(i);
713     sum_ZPC += reco.GetZN2HREnTow(i);
714     sum_ZNA += reco.GetZP1HREnTow(i);
715     sum_ZPA += reco.GetZP2HREnTow(i);
716     if(i==0){
717       FillRecPointsData(8,reco.GetZN1HREnTow(i));
718       FillRecPointsData(9,reco.GetZN2HREnTow(i));
719       FillRecPointsData(10,reco.GetZP1HREnTow(i));
720       FillRecPointsData(11,reco.GetZP2HREnTow(i));
721     }
722     else{
723       sumQ_ZNC += reco.GetZN1HREnTow(i);
724       sumQ_ZPC += reco.GetZN2HREnTow(i);
725       sumQ_ZNA += reco.GetZP1HREnTow(i);
726       sumQ_ZPA += reco.GetZP2HREnTow(i);
727     }
728   }
729   
730   FillRecPointsData(0,sum_ZNC);
731   FillRecPointsData(1,sum_ZNA);
732   FillRecPointsData(2,sum_ZPC);
733   FillRecPointsData(3,sum_ZPA);
734   //
735   FillRecPointsData(4,sumQ_ZNC);
736   FillRecPointsData(5,sumQ_ZNA);
737   FillRecPointsData(6,sumQ_ZPC);
738   FillRecPointsData(7,sumQ_ZPA);
739   //
740   IncEvCountCycleRecPoints();
741   IncEvCountTotalRecPoints();
742   //  
743 }  
744
745 //____________________________________________________________________________
746 void AliZDCQADataMakerRec::MakeESDs(AliESDEvent * esd)
747 {
748   // make QA data from ESDs
749   //
750   
751   // Check id histograms already created for this Event Specie
752   if(!GetESDsData(0)) InitESDs() ;
753
754   AliESDZDC * zdcESD =  esd->GetESDZDC();
755   //
756   /*TString beamType = esd->GetBeamType();
757   Double_t centr_ZNC[2]={999.,999}, centr_ZNA[2]={999.,999};
758   if(((beamType.CompareTo("pp"))==0) || ((beamType.CompareTo("p-p"))==0)
759      ||((beamType.CompareTo("PP"))==0) || ((beamType.CompareTo("P-P"))==0)){
760     zdcESD->GetZNCentroidInpp(centr_ZNC, centr_ZNA);
761   }
762   else if((beamType.CompareTo("A-A")) == 0 || (beamType.CompareTo("Pb-Pb")) == 0){
763     Float_t beamEne = esd->GetBeamEnergy();
764     zdcESD->GetZNCentroidInPbPb(beamEne, centr_ZNC, centr_ZNA);
765   }
766   else printf("\n WARNING!!! AliZDCQADataMakerRec::MakeESDs: can't calculate centroids for beam type: %s\n\n",beamType.Data());
767   FillESDsData(0,centr_ZNC[0], centr_ZNC[1]);
768   FillESDsData(1,centr_ZNA[0], centr_ZNA[1]);*/
769
770   FillESDsData(0,esd->GetZDCN1Energy());
771   FillESDsData(1,esd->GetZDCN2Energy());
772   FillESDsData(2,esd->GetZDCP1Energy());
773   FillESDsData(3,esd->GetZDCP2Energy());
774   FillESDsData(4,esd->GetZDCEMEnergy(0));
775   FillESDsData(5,esd->GetZDCEMEnergy(1));
776   //
777   Double_t sumQZNC=0., sumQZPC=0., sumQZNA=0., sumQZPA=0.;
778   //Double_t sumQZNC_lg=0., sumQZPC_lg=0., sumQZNA_lg=0., sumQZPA_lg=0.;
779   //
780   const Double_t *towZNC, *towZPC, *towZNA, *towZPA;
781   //const Double_t *towZNC_lg, *towZPC_lg, *towZNA_lg, *towZPA_lg;
782   //
783   towZNC = zdcESD->GetZN1TowerEnergy();
784   towZPC = zdcESD->GetZP1TowerEnergy();
785   towZNA = zdcESD->GetZN2TowerEnergy();
786   towZPA = zdcESD->GetZP2TowerEnergy();
787   //
788   /*towZNC_lg = zdcESD->GetZN1TowerEnergyLR();
789   towZPC_lg = zdcESD->GetZP1TowerEnergyLR();
790   towZNA_lg = zdcESD->GetZN2TowerEnergyLR();
791   towZPA_lg = zdcESD->GetZP2TowerEnergyLR();*/
792   //
793   for(Int_t i=0; i<5; i++){
794      if(i==0){
795        FillESDsData(10,towZNC[i]);
796        FillESDsData(11,towZNA[i]);
797        FillESDsData(12,towZPC[i]);
798        FillESDsData(13,towZPA[i]);
799      }
800      else{
801        sumQZNC += towZNC[i];
802        sumQZPC += towZPC[i];
803        sumQZNA += towZNA[i];
804        sumQZPA += towZPA[i];
805      }
806   }
807   FillESDsData(6,sumQZNC);
808   FillESDsData(7,sumQZNA);
809   FillESDsData(8,sumQZPC);
810   FillESDsData(9,sumQZPA);
811   //
812   IncEvCountCycleESDs();
813   IncEvCountTotalESDs();
814   //
815 }
816
817 //____________________________________________________________________________
818 void AliZDCQADataMakerRec::StartOfDetectorCycle()
819 {
820   //Detector specific actions at start of cycle
821
822   fPedCalibData = GetPedCalibData();
823   
824 }
825
826 //____________________________________________________________________________ 
827 void AliZDCQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
828 {
829   //Detector specific actions at end of cycle
830   // do the QA checking
831   //
832   ResetEventTrigClasses();
833   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {  // RS: loop over event types
834     //
835     if (!IsValidEventSpecie(specie, list)) continue;
836     SetEventSpecie(AliRecoParam::ConvertIndex(specie));
837     //
838     for (int itc=-1;itc<GetNTrigClasses();itc++) { // RS: loop over eventual clones per trigger class
839       //
840       if( task == AliQAv1::kRAWS) {
841         TH1* h16 = GetRawsData(16, itc);
842         TH1* h17 = GetRawsData(17, itc);
843         TH1* h18 = GetRawsData(18, itc);
844         TH1* h19 = GetRawsData(19, itc);
845         //      TH1* h20 = GetRawsData(20, itc);
846         //      TH1* h21 = GetRawsData(21, itc);
847         if (!h16 || !h17 || !h18 || !h19) {
848          AliWarning("AliZDCQADataMaker Rec -> Not all histogram for DQM found!"); 
849          AliWarning(Form("for specie %s and trigger class %s",
850                          AliRecoParam::GetEventSpecieName(specie), AliQADataMaker::GetTrigClassName(itc)));
851         }
852         else{
853           h16->Draw("");
854           for(Int_t ibin=1; ibin<=h16->GetNbinsX(); ibin++){
855             h18->SetBinContent(ibin, h16->GetBinContent(ibin)); 
856             h18->SetBinError(ibin, h16->GetBinError(ibin));
857           }
858           for(Int_t ibin=1; ibin<=h17->GetNbinsX(); ibin++){
859             h19->SetBinContent(ibin, h17->GetBinContent(ibin)); 
860             h19->SetBinError(ibin, h17->GetBinError(ibin));
861           }
862           h18->SetLineColor(kBlue); h18->SetLineWidth(2);
863           h19->SetLineColor(kAzure-3); h19->SetLineWidth(2);
864         }
865       }
866     } // loop over t
867   } //  loop over species
868         
869   AliQAChecker::Instance()->Run(AliQAv1::kZDC, task, list) ;  
870 }