]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALQADataMakerRec.cxx
Correct setting of non linearity correction in Tender
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALQADataMakerRec.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 Based on the QA code for PHOS written by Yves Schutz July 2007
17
18 Authors:  J.Klay (Cal Poly) May 2008
19           S. Salur LBL April 2008
20  
21 Created one histogram for QA shifter;-- Yaxian Mao: 11/2009
22 The idea:average counts for all the towers should be flat 
23 Change all existing histograms as experts
24
25 Change histograms for DQM shifter: --  Yaxian Mao 04/2010
26 Calcuate the amplitude ratio from current run and the LED reference, for QAChecker use
27 Also calculate the ratio of amplitude from LED Monitor system (current/Reference), to check LED system  
28  
29 */
30
31 // --- ROOT system ---
32 #include <TClonesArray.h>
33 #include <TFile.h> 
34 #include <TH1F.h> 
35 #include <TH1I.h> 
36 #include <TH2F.h> 
37 #include <TLine.h>
38 #include <TText.h>
39 #include <TProfile.h> 
40 #include <TProfile2D.h> 
41 #include <TStyle.h>
42 // --- Standard library ---
43
44
45 // --- AliRoot header files ---
46 #include "AliDAQ.h"
47 #include "AliESDCaloCluster.h"
48 #include "AliESDCaloCells.h"
49 #include "AliESDEvent.h"
50 #include "AliLog.h"
51 #include "AliEMCALQADataMakerRec.h"
52 #include "AliQAChecker.h"
53 #include "AliEMCALDigit.h" 
54 #include "AliEMCALRecPoint.h" 
55 #include "AliEMCALRawUtils.h"
56 #include "AliEMCALReconstructor.h"
57 #include "AliEMCALRecParam.h"
58 #include "AliRawReader.h"
59 #include "AliCaloRawStreamV3.h"
60 #include "AliEMCALGeoParams.h"
61 #include "AliRawEventHeaderBase.h"
62 #include "AliQAManager.h"
63 #include "AliCDBEntry.h"
64
65 #include "AliCaloBunchInfo.h"
66 #include "AliCaloFitResults.h"
67 #include "AliCaloRawAnalyzerFastFit.h"
68 #include "AliCaloRawAnalyzerNN.h"
69 //#include "AliCaloRawAnalyzerLMS.h"
70 #include "AliCaloRawAnalyzerKStandard.h"
71 #include "AliCaloRawAnalyzerPeakFinder.h"
72 #include "AliCaloRawAnalyzerCrude.h"
73 #include "AliEMCALGeometry.h"
74 #include "AliEMCALTriggerSTURawStream.h"
75
76 #include "AliCaloRawAnalyzerFactory.h"
77
78 using namespace std;
79
80 ClassImp(AliEMCALQADataMakerRec)
81            
82 //____________________________________________________________________________ 
83 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(fitAlgorithm fitAlgo) : 
84   AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kEMCAL), "EMCAL Quality Assurance Data Maker"),
85   fFittingAlgorithm(0),
86   fRawAnalyzer(0),
87   fRawAnalyzerTRU(0),
88         fGeom(0),
89   fSuperModules(10), // FIXME!!! number of SuperModules; 10 for 2011; update default for later runs 
90   fFirstPedestalSample(0),
91   fLastPedestalSample(3),
92   fFirstPedestalSampleTRU(0),
93   fLastPedestalSampleTRU(3),
94   fMinSignalLG(0),
95   fMaxSignalLG(AliEMCALGeoParams::fgkSampleMax),
96   fMinSignalHG(0),
97   fMaxSignalHG(AliEMCALGeoParams::fgkSampleMax),
98   fMinSignalTRU(0),
99   fMaxSignalTRU(AliEMCALGeoParams::fgkSampleMax),
100   fMinSignalLGLEDMon(0),
101   fMaxSignalLGLEDMon(AliEMCALGeoParams::fgkSampleMax),
102   fMinSignalHGLEDMon(0),
103   fMaxSignalHGLEDMon(AliEMCALGeoParams::fgkSampleMax),
104   fCalibRefHistoPro(NULL),
105   fCalibRefHistoH2F(NULL),
106   fLEDMonRefHistoPro(NULL),
107   fHighEmcHistoH2F(NULL)
108 //  fTextSM(new TText*[fSuperModules]) ,
109 //  fLineCol(NULL),
110 //  fLineRow(NULL)
111
112 {
113   // ctor
114   SetFittingAlgorithm(fitAlgo);
115   
116   //fRawAnalyzerTRU = new AliCaloRawAnalyzerLMS();
117   
118   fRawAnalyzerTRU =  ( AliCaloRawAnalyzerKStandard*)AliCaloRawAnalyzerFactory::CreateAnalyzer(kLMS);
119   
120   fRawAnalyzerTRU->SetFixTau(kTRUE); 
121   fRawAnalyzerTRU->SetTau(2.5); // default for TRU shaper
122
123         fGeom = new AliEMCALGeometry("EMCAL_COMPLETEV1", "EMCAL");
124 //  for (Int_t sm = 0 ; sm < fSuperModules ; sm++){
125 //    fTextSM[sm] = NULL ;
126 //  }
127 }
128
129 //____________________________________________________________________________ 
130 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(const AliEMCALQADataMakerRec& qadm) :
131   AliQADataMakerRec(), 
132   fFittingAlgorithm(0),
133   fRawAnalyzer(0),
134   fRawAnalyzerTRU(0),
135         fGeom(0),
136   fSuperModules(qadm.GetSuperModules()), 
137   fFirstPedestalSample(qadm.GetFirstPedestalSample()), 
138   fLastPedestalSample(qadm.GetLastPedestalSample()),  
139   fFirstPedestalSampleTRU(qadm.GetFirstPedestalSampleTRU()), 
140   fLastPedestalSampleTRU(qadm.GetLastPedestalSampleTRU()),  
141   fMinSignalLG(qadm.GetMinSignalLG()),
142   fMaxSignalLG(qadm.GetMaxSignalLG()),
143   fMinSignalHG(qadm.GetMinSignalHG()),
144   fMaxSignalHG(qadm.GetMaxSignalHG()),
145   fMinSignalTRU(qadm.GetMinSignalTRU()),
146   fMaxSignalTRU(qadm.GetMaxSignalTRU()),
147   fMinSignalLGLEDMon(qadm.GetMinSignalLGLEDMon()),
148   fMaxSignalLGLEDMon(qadm.GetMaxSignalLGLEDMon()),
149   fMinSignalHGLEDMon(qadm.GetMinSignalHGLEDMon()),
150   fMaxSignalHGLEDMon(qadm.GetMaxSignalHGLEDMon()),
151   fCalibRefHistoPro(NULL),
152   fCalibRefHistoH2F(NULL),
153   fLEDMonRefHistoPro(NULL),
154   fHighEmcHistoH2F(NULL)
155 //  fTextSM(new TText*[fSuperModules]) ,
156 //  fLineCol(NULL),
157 //  fLineRow(NULL)
158 {
159   //copy ctor 
160   SetName((const char*)qadm.GetName()) ; 
161   SetTitle((const char*)qadm.GetTitle()); 
162   SetFittingAlgorithm(qadm.GetFittingAlgorithm());
163   
164   //fRawAnalyzerTRU = new AliCaloRawAnalyzerLMS();
165   fRawAnalyzerTRU = (AliCaloRawAnalyzerKStandard*)AliCaloRawAnalyzerFactory::CreateAnalyzer(kLMS);
166   fRawAnalyzerTRU->SetFixTau(kTRUE); 
167   fRawAnalyzerTRU->SetTau(2.5); // default for TRU shaper
168 //  for (Int_t sm = 0 ; sm < fSuperModules ; sm++){
169 //    fTextSM[sm] = qadm.fTextSM[sm] ;
170 //  }  
171 }
172
173 //__________________________________________________________________
174 AliEMCALQADataMakerRec& AliEMCALQADataMakerRec::operator = (const AliEMCALQADataMakerRec& qadm )
175 {
176   // Equal operator.
177   this->~AliEMCALQADataMakerRec();
178   new(this) AliEMCALQADataMakerRec(qadm);
179 //  fLineCol = NULL;
180 //  fLineRow = NULL;
181 //  for (Int_t sm = 0 ; sm < fSuperModules ; sm++){
182 //    fTextSM[sm] = qadm.fTextSM[sm] ;
183 //  }    
184   return *this;
185 }
186  
187 //____________________________________________________________________________ 
188 void AliEMCALQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
189 {
190   //Detector specific actions at end of cycle
191         
192 //  if(fCycleCounter)
193 //        GetRawsData(kNEventsPerTower)->Scale(1./fCycleCounter);
194
195   // do the QA checking
196   ResetEventTrigClasses(); // reset triggers list to select all histos
197   AliQAChecker::Instance()->Run(AliQAv1::kEMCAL, task, list) ;  
198 }
199
200 //____________________________________________________________________________ 
201 void AliEMCALQADataMakerRec::GetCalibRefFromOCDB()
202 {
203   //Get the reference histogram from OCDB
204   TString sName1("hHighEmcalRawMaxMinusMin") ;
205   TString sName2("hLowLEDMonEmcalRawMaxMinusMin") ;
206   sName1.Prepend(Form("%s_", AliRecoParam::GetEventSpecieName(AliRecoParam::kCalib))) ; 
207   sName2.Prepend(Form("%s_", AliRecoParam::GetEventSpecieName(AliRecoParam::kCalib))) ; 
208   
209   TString refStorage(AliQAv1::GetQARefStorage()) ;
210   if (!refStorage.Contains(AliQAv1::GetLabLocalOCDB()) && !refStorage.Contains(AliQAv1::GetLabAliEnOCDB())) {
211     AliFatal(Form("%s is not a valid location for reference data", refStorage.Data())) ; 
212   } else {
213     AliQAManager* manQA = AliQAManager::QAManager(AliQAv1::kRAWS) ;    
214     AliQAv1::SetQARefDataDirName(AliRecoParam::kCalib) ;
215     if ( ! manQA->GetLock() ) { 
216       manQA->SetDefaultStorage(AliQAv1::GetQARefStorage()) ; 
217       manQA->SetSpecificStorage("*", AliQAv1::GetQARefStorage()) ;
218       manQA->SetRun(AliCDBManager::Instance()->GetRun()) ; 
219       manQA->SetLock() ; 
220     }
221     char * detOCDBDir = Form("%s/%s/%s", GetName(), AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ; 
222     AliCDBEntry * entry = manQA->Get(detOCDBDir, manQA->GetRun()) ;
223     if (entry) {
224       TList * listDetQAD =static_cast<TList *>(entry->GetObject()) ;
225       if ( strcmp(listDetQAD->ClassName(), "TList") != 0 ) {
226         AliError(Form("Expected a Tlist and found a %s for detector %s", listDetQAD->ClassName(), GetName())) ; 
227         listDetQAD = NULL ; 
228       }
229       TObjArray * dirOCDB= NULL ; 
230       if ( listDetQAD )
231         dirOCDB = static_cast<TObjArray *>(listDetQAD->FindObject(Form("%s/%s", AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), AliRecoParam::GetEventSpecieName(AliRecoParam::kCalib)))) ;       
232       if (dirOCDB){
233         fCalibRefHistoPro = dynamic_cast<TProfile *>(dirOCDB->FindObject(sName1.Data())) ; 
234         fLEDMonRefHistoPro = dynamic_cast<TProfile *>(dirOCDB->FindObject(sName2.Data())) ; 
235       }
236     }
237   }
238  
239   if(fCalibRefHistoPro && fLEDMonRefHistoPro){
240     
241     //Defining histograms binning, each 2D histogram covers all SMs
242     Int_t nSMSectors = fSuperModules / 2; // 2 SMs per sector
243     Int_t nbinsZ = 2*AliEMCALGeoParams::fgkEMCALCols;
244     Int_t nbinsPhi = nSMSectors * AliEMCALGeoParams::fgkEMCALRows;
245     
246     if(!fCalibRefHistoH2F)
247       fCalibRefHistoH2F =  new TH2F("hCalibRefHisto", "hCalibRefHisto", nbinsZ, -0.5, nbinsZ - 0.5, nbinsPhi, -0.5, nbinsPhi -0.5);
248     ConvertProfile2H(fCalibRefHistoPro,fCalibRefHistoH2F) ; 
249   } else {
250     AliFatal(Form("No reference object with name %s or %s found", sName1.Data(), sName2.Data())) ; 
251   }
252 }
253 //____________________________________________________________________________ 
254 void AliEMCALQADataMakerRec::InitESDs()
255 {
256   //Create histograms to controll ESD
257   const Bool_t expert   = kTRUE ; 
258   const Bool_t image    = kTRUE ; 
259   
260   TH1F * h1 = new TH1F("hESDCaloClusterE",  "ESDs CaloCluster energy in EMCAL;Energy [GeV];Counts",    200, 0., 100.) ; 
261   h1->Sumw2() ;
262   Add2ESDsList(h1, kESDCaloClusE, !expert, image)  ;                                                     
263
264   TH1I * h2 = new TH1I("hESDCaloClusterM", "ESDs CaloCluster multiplicity in EMCAL;# of Clusters;Entries", 100, 0,  100) ; 
265   h2->Sumw2() ;
266   Add2ESDsList(h2, kESDCaloClusM, !expert, image)  ;
267
268   TH1F * h3 = new TH1F("hESDCaloCellA",  "ESDs CaloCell amplitude in EMCAL;Energy [GeV];Counts",    500, 0., 50.) ; 
269   h3->Sumw2() ;
270   Add2ESDsList(h3, kESDCaloCellA, !expert, image)  ;  
271  
272   TH1I * h4 = new TH1I("hESDCaloCellM", "ESDs CaloCell multiplicity in EMCAL;# of Clusters;Entries", 200, 0,  1000) ; 
273   h4->Sumw2() ;
274   Add2ESDsList(h4, kESDCaloCellM, !expert, image) ;
275   //
276   ClonePerTrigClass(AliQAv1::kESDS); // this should be the last line    
277 }
278
279 //____________________________________________________________________________ 
280 void AliEMCALQADataMakerRec::InitDigits()
281 {
282   // create Digits histograms in Digits subdir
283   const Bool_t expert   = kTRUE ; 
284   const Bool_t image    = kTRUE ; 
285   
286   TH1I * h0 = new TH1I("hEmcalDigits",    "Digits amplitude distribution in EMCAL;Amplitude [ADC counts];Counts",    500, 0, 500) ; 
287   h0->Sumw2() ;
288   Add2DigitsList(h0, 0, !expert, image) ;
289   TH1I * h1 = new TH1I("hEmcalDigitsMul", "Digits multiplicity distribution in EMCAL;# of Digits;Entries", 200, 0, 2000) ; 
290   h1->Sumw2() ;
291   Add2DigitsList(h1, 1, !expert, image) ;
292   //
293   ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
294 }
295
296 //____________________________________________________________________________ 
297 void AliEMCALQADataMakerRec::InitRecPoints()
298 {
299   // create Reconstructed PoInt_ts histograms in RecPoints subdir
300   const Bool_t expert   = kTRUE ; 
301   const Bool_t image    = kTRUE ; 
302   
303   TH1F* h0 = new TH1F("hEMCALRpE","EMCAL RecPoint energies;Energy [GeV];Counts",200, 0.,20.); //GeV
304   h0->Sumw2();
305   Add2RecPointsList(h0,kRecPE, !expert, image);
306
307   TH1I* h1 = new TH1I("hEMCALRpM","EMCAL RecPoint multiplicities;# of Clusters;Entries",100,0,100);
308   h1->Sumw2();
309   Add2RecPointsList(h1,kRecPM, !expert, image);
310
311   TH1I* h2 = new TH1I("hEMCALRpDigM","EMCAL RecPoint Digit Multiplicities;# of Digits;Entries",20,0,20);
312   h2->Sumw2();
313   Add2RecPointsList(h2,kRecPDigM, !expert, image);
314   //
315   ClonePerTrigClass(AliQAv1::kRECPOINTS); // this should be the last line
316 }
317
318 //____________________________________________________________________________ 
319 void AliEMCALQADataMakerRec::InitRaws()
320 {
321   // create Raws histograms in Raws subdir
322   const Bool_t expert   = kTRUE ; 
323   const Bool_t saveCorr = kTRUE ; 
324   const Bool_t image    = kTRUE ; 
325   const Option_t *profileOption = "s";
326
327   Int_t nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
328   Int_t nTot = fSuperModules * nTowersPerSM; // max number of towers in all SuperModules
329     
330   //Defining histograms binning, each 2D histogram covers all SMs
331   Int_t nSMSectors = fSuperModules / 2; // 2 SMs per sector
332   Int_t nbinsZ = 2*AliEMCALGeoParams::fgkEMCALCols;
333   Int_t nbinsPhi = nSMSectors * AliEMCALGeoParams::fgkEMCALRows;
334         
335   Int_t nTRUCols = 2*AliEMCALGeoParams::fgkEMCALTRUCols; //total TRU columns for 2D TRU histos
336   Int_t nTRURows = nSMSectors*AliEMCALGeoParams::fgkEMCALTRUsPerSM*AliEMCALGeoParams::fgkEMCALTRURows; //total TRU rows for 2D TRU histos
337    // counter info: number of channels per event (bins are SM index)
338   TProfile * h0 = new TProfile("hLowEmcalSupermodules", "Low Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
339                                fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
340   Add2RawsList(h0, kNsmodLG, expert, !image, !saveCorr) ;
341   TProfile * h1 = new TProfile("hHighEmcalSupermodules", "High Gain EMC: # of towers vs SuperMod;SM Id;# of towers",  
342                                fSuperModules, -0.5, fSuperModules-0.5, profileOption) ; 
343   Add2RawsList(h1, kNsmodHG, expert, !image, !saveCorr) ;
344
345   // where did max sample occur? (bins are towers)
346   TProfile * h2 = new TProfile("hLowEmcalRawtime", "Low Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]", 
347                                nTot, -0.5, nTot-0.5, profileOption) ;
348   Add2RawsList(h2, kTimeLG, expert, !image, !saveCorr) ;
349   TProfile * h3 = new TProfile("hHighEmcalRawtime", "High Gain EMC: Time at Max vs towerId;Tower Id;Time [ticks]", 
350                                nTot, -0.5, nTot-0.5, profileOption) ;
351   Add2RawsList(h3, kTimeHG, expert, !image, !saveCorr) ;
352
353   // how much above pedestal was the max sample?  (bins are towers)
354   TProfile * h4 = new TProfile("hLowEmcalRawMaxMinusMin", "Low Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]", 
355                                nTot, -0.5, nTot-0.5, profileOption) ;
356   Add2RawsList(h4, kSigLG, expert, !image, !saveCorr) ;
357   TProfile * h5 = new TProfile("hHighEmcalRawMaxMinusMin", "High Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
358                                nTot, -0.5, nTot-0.5, profileOption) ;
359   Add2RawsList(h5, kSigHG, expert, !image, !saveCorr) ;
360
361   // total counter: channels per event
362   TH1I * h6 = new TH1I("hLowNtot", "Low Gain EMC: Total Number of found towers;# of Towers;Counts", 200, 0, nTot) ;
363   h6->Sumw2() ;
364   Add2RawsList(h6, kNtotLG, expert, !image, !saveCorr) ;
365   TH1I * h7 = new TH1I("hHighNtot", "High Gain EMC: Total Number of found towers;# of Towers;Counts", 200,0, nTot) ;
366   h7->Sumw2() ;
367   Add2RawsList(h7, kNtotHG, expert, !image, !saveCorr) ;
368
369   // pedestal (bins are towers)
370   TProfile * h8 = new TProfile("hLowEmcalRawPed", "Low Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]", 
371                                nTot, -0.5, nTot-0.5, profileOption) ;
372   Add2RawsList(h8, kPedLG, expert, !image, !saveCorr) ;
373   TProfile * h9 = new TProfile("hHighEmcalRawPed", "High Gain EMC: Pedestal vs towerId;Tower Id;Pedestal [ADC counts]",
374                                nTot, -0.5, nTot-0.5, profileOption) ;
375   Add2RawsList(h9, kPedHG, expert, !image, !saveCorr) ;
376         
377  
378   // now repeat the same for TRU and LEDMon data
379   Int_t nTot2x2 = fSuperModules * AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // max number of TRU channels for all SuperModules
380
381   // counter info: number of channels per event (bins are SM index)
382   TProfile * hT0 = new TProfile("hTRUEmcalSupermodules", "TRU EMC: # of TRU channels vs SuperMod;SM Id;# of TRU channels",
383                                 fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
384   Add2RawsList(hT0, kNsmodTRU, expert, !image, !saveCorr) ;
385
386   // how much above pedestal was the max sample?  (bins are TRU channels)
387   TProfile * hT1 = new TProfile("hTRUEmcalRawMaxMinusMin", "TRU EMC: Max - Min vs 2x2Id;2x2 Id;Max-Min [ADC counts]", 
388                                 nTot2x2, -0.5, nTot2x2-0.5, profileOption) ;
389   Add2RawsList(hT1, kSigTRU, expert, !image, !saveCorr) ;
390
391   // total counter: channels per event
392   TH1I * hT2 = new TH1I("hTRUNtot", "TRU EMC: Total Number of found TRU channels;# of TRU Channels;Counts", 200, 0, nTot2x2) ;
393   hT2->Sumw2() ;
394   Add2RawsList(hT2, kNtotTRU, expert, !image, !saveCorr) ;
395
396   // L0 trigger hits: # of hits (bins are TRU channels)
397   TH2I * hT3 = new TH2I("hTRUEmcalL0hits", "L0 trigger hits: Total number of 2x2 L0 generated",  nTRUCols, -0.5, nTRUCols - 0.5, nTRURows, -0.5, nTRURows-0.5);
398   hT3->SetOption("COLZ");
399   //hT3->Sumw2();
400   Add2RawsList(hT3, kNL0TRU, !expert, image, !saveCorr);
401
402   // L0 trigger hits: average time (bins are TRU channels)
403   TProfile2D * hT4 = new TProfile2D("hTRUEmcalL0hitsAvgTime", "L0 trigger hits: average time bin", nTRUCols, -0.5, nTRUCols - 0.5, nTRURows, -0.5, nTRURows-0.5, profileOption);
404   hT4->SetOption("COLZ");
405   Add2RawsList(hT4, kTimeL0TRU, !expert, image, !saveCorr);
406
407   // L0 trigger hits: first in the event (bins are TRU channels)
408   TH1I * hT5 = new TH1I("hTRUEmcalL0hitsFirst", "L0 trigger hits: First hit in the event", nTot2x2, -0.5, nTot2x2);
409   hT5->Sumw2();
410   Add2RawsList(hT5, kNL0FirstTRU, expert, !image, !saveCorr);
411         
412   // L0 trigger hits: average time of first hit in the event (bins are TRU channels)
413   TProfile * hT6 = new TProfile("hTRUEmcalL0hitsFirstAvgTime", "L0 trigger hits: average time of first hit", nTot2x2, -0.5, nTot2x2, profileOption); 
414   Add2RawsList(hT6, kTimeL0FirstTRU, expert, !image, !saveCorr);
415
416   // and also LED Mon..
417   // LEDMon has both high and low gain channels, just as regular FEE/towers
418   Int_t nTotLEDMon = fSuperModules * AliEMCALGeoParams::fgkEMCALLEDRefs; // max number of LEDMon channels for all SuperModules
419
420   // counter info: number of channels per event (bins are SM index)
421   TProfile * hL0 = new TProfile("hLowLEDMonEmcalSupermodules", "LowLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",
422                                fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
423   Add2RawsList(hL0, kNsmodLGLEDMon, expert, !image, !saveCorr) ;
424   TProfile * hL1 = new TProfile("hHighLEDMonEmcalSupermodules", "HighLEDMon Gain EMC: # of strips vs SuperMod;SM Id;# of strips",  
425                                fSuperModules, -0.5, fSuperModules-0.5, profileOption) ; 
426   Add2RawsList(hL1, kNsmodHGLEDMon, expert, !image, !saveCorr) ;
427
428   // where did max sample occur? (bins are strips)
429   TProfile * hL2 = new TProfile("hLowLEDMonEmcalRawtime", "LowLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]", 
430                                nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
431   Add2RawsList(hL2, kTimeLGLEDMon, expert, !image, !saveCorr) ;
432   TProfile * hL3 = new TProfile("hHighLEDMonEmcalRawtime", "HighLEDMon Gain EMC: Time at Max vs stripId;Strip Id;Time [ticks]", 
433                                nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
434   Add2RawsList(hL3, kTimeHGLEDMon, expert, !image, !saveCorr) ;
435
436   // how much above pedestal was the max sample?  (bins are strips)
437   TProfile * hL4 = new TProfile("hLowLEDMonEmcalRawMaxMinusMin", "LowLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]", 
438                                nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
439   Add2RawsList(hL4, kSigLGLEDMon, expert, !image, !saveCorr) ;
440   TProfile * hL5 = new TProfile("hHighLEDMonEmcalRawMaxMinusMin", "HighLEDMon Gain EMC: Max - Min vs stripId;Strip Id;Max-Min [ADC counts]",
441                                nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
442   Add2RawsList(hL5, kSigHGLEDMon, expert, !image, !saveCorr) ;
443   
444     // total counter: channels per event
445   TH1I * hL6 = new TH1I("hLowLEDMonNtot", "LowLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200, 0, nTotLEDMon) ;
446   hL6->Sumw2() ;
447   Add2RawsList(hL6, kNtotLGLEDMon, expert, !image, !saveCorr) ;
448   TH1I * hL7 = new TH1I("hHighLEDMonNtot", "HighLEDMon Gain EMC: Total Number of found strips;# of Strips;Counts", 200,0, nTotLEDMon) ;
449   hL7->Sumw2() ;
450   Add2RawsList(hL7, kNtotHGLEDMon, expert, !image, !saveCorr) ;
451
452   // pedestal (bins are strips)
453   TProfile * hL8 = new TProfile("hLowLEDMonEmcalRawPed", "LowLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]", 
454                                nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
455   Add2RawsList(hL8, kPedLGLEDMon, expert, !image, !saveCorr) ;
456   TProfile * hL9 = new TProfile("hHighLEDMonEmcalRawPed", "HighLEDMon Gain EMC: Pedestal vs stripId;Strip Id;Pedestal [ADC counts]",
457                                nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
458   Add2RawsList(hL9, kPedHGLEDMon, expert, !image, !saveCorr) ;
459   
460   //temp 2D amplitude histogram for the current run
461   fHighEmcHistoH2F = new TH2F("h2DHighEC2", "High Gain EMC:Max - Min [ADC counts]", nbinsZ, -0.5 , nbinsZ-0.5, nbinsPhi, -0.5, nbinsPhi-0.5);
462    fHighEmcHistoH2F->SetDirectory(0) ; // this histo must be memory resident
463   //add ratio histograms: to comapre the current run with the reference data 
464   TH2F * h15 = new TH2F("h2DRatioAmp", "High Gain Ratio to Reference:Amplitude_{current run}/Amplitude_{reference run}", nbinsZ, -0.5 , nbinsZ-0.5, 
465                         nbinsPhi, -0.5, nbinsPhi-0.5);
466   //settings for display in amore
467   h15->SetTitle("Amplitude_{current run}/Amplitude_{reference run}"); 
468   h15->SetMaximum(2.0);
469   h15->SetMinimum(0.1);
470   h15->SetOption("COLZ");
471   gStyle->SetOptStat(0);
472   Int_t color[] = {4,3,2} ;
473   gStyle->SetPalette(3,color);
474   h15->GetZaxis()->SetNdivisions(3);
475   h15->UseCurrentStyle();
476   h15->SetDirectory(0);
477   Add2RawsList(h15, k2DRatioAmp, !expert, image, !saveCorr) ;
478
479   TH1F * h16 = new TH1F("hRatioDist", "Amplitude_{current run}/Amplitude_{reference run} ratio distribution", nTot, 0., 2.);
480   h16->SetMinimum(0.1); 
481   h16->SetMaximum(100.);
482   gStyle->SetOptStat(0);
483   h16->UseCurrentStyle();
484   h16->SetDirectory(0);
485   Add2RawsList(h16, kRatioDist, !expert, image, !saveCorr) ;
486
487   //add two histograms for shifter from the LED monitor system: comapre LED monitor with the reference run
488   //to be used for decision whether we need to change reference data
489   TH1F * hL10 = new TH1F("hMaxMinusMinLEDMonRatio", "LEDMon amplitude, Ratio to reference run", nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
490   //settings for display in amore
491   hL10->SetTitle("Amplitude_{LEDMon current}/Amplitude_{LEDMon reference}"); 
492   hL10->SetMaximum(2.0);
493   hL10->SetMinimum(0.1); 
494   gStyle->SetOptStat(0);
495   hL10->UseCurrentStyle();
496   hL10->SetDirectory(0);
497 //  hL10->SetOption("E");
498   Add2RawsList(hL10, kLEDMonRatio, !expert, image, !saveCorr) ;
499
500   TH1F * hL11 = new TH1F("hMaxMinusMinLEDMonRatioDist", "LEDMon amplitude, Ratio distribution", nTotLEDMon, 0, 2);
501   hL11->SetMinimum(0.1) ;
502   gStyle->SetOptStat(0);
503   hL11->UseCurrentStyle();
504   hL11->SetDirectory(0);
505   Add2RawsList(hL11, kLEDMonRatioDist, !expert, image, !saveCorr) ;
506   
507   GetCalibRefFromOCDB();   
508
509
510         //STU histgrams
511
512  //histos
513  Int_t nSTUCols = AliEMCALGeoParams::fgkEMCALSTUCols;
514  Int_t nSTURows = AliEMCALGeoParams::fgkEMCALSTURows;
515 //              kAmpL1, kGL1, kJL1,
516 //              kGL1V0, kJL1V0, kSTUTRU  
517         
518  TProfile2D *hS0 = new TProfile2D("hL1Amp", "Mean STU signal per Row and Column", nSTUCols, -0.5, nSTUCols-0.5, nSTURows, -0.5, nSTURows-0.5);
519  Add2RawsList(hS0, kAmpL1, expert, !image, !saveCorr) ;
520         
521  TH2F *hS1 = new TH2F("hL1Gamma", "L1 Gamma patch position (FastOR top-left)", nSTUCols, -0.50, nSTUCols-0.5, nSTURows, -0.5, nSTURows-0.5);
522  Add2RawsList(hS1, kGL1, expert, image, !saveCorr) ;
523         
524  TH2F *hS2 = new TH2F("hL1Jet", "L1 Jet patch position (FastOR top-left)", 12, -0.5, nSTUCols-0.5, 16, 0, nSTURows-0.5);
525  Add2RawsList(hS2, kJL1, expert, image, !saveCorr) ;
526         
527  TH2I *hS3 = new TH2I("hL1GV0", "L1 Gamma patch amplitude versus V0 signal", 500, 0, 50000, 1500, 0, 1500);
528  Add2RawsList(hS3, kGL1V0, expert, !image, !saveCorr) ;
529         
530  TH2I *hS4 = new TH2I("hL1JV0", "L1 Jet patch amplitude versus V0 signal", 500, 0, 50000, 1000, 0, 1000);
531  Add2RawsList(hS4, kJL1V0, expert, !image, !saveCorr) ;
532
533  TH1I *hS5 = new TH1I("hFrameR","Link between TRU and STU", 32, 0, 32);
534  Add2RawsList(hS5, kSTUTRU, expert, !image, !saveCorr) ;
535
536
537
538   //
539   ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
540 }
541
542 //____________________________________________________________________________
543 void AliEMCALQADataMakerRec::MakeESDs(AliESDEvent * esd)
544 {
545   // make QA data from ESDs
546
547   Int_t nTot = 0 ; 
548   for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
549     AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
550     if( clu->IsEMCAL() ) {
551       FillESDsData(kESDCaloClusE,clu->E()) ;
552       nTot++ ;
553     } 
554   }
555   FillESDsData(kESDCaloClusM,nTot) ;
556
557   //fill calo cells
558   AliESDCaloCells* cells = esd->GetEMCALCells();
559   FillESDsData(kESDCaloCellM,cells->GetNumberOfCells()) ;
560
561   for ( Int_t index = 0; index < cells->GetNumberOfCells() ; index++ ) {
562     FillESDsData(kESDCaloCellA,cells->GetAmplitude(index)) ;
563   }
564   //
565   IncEvCountCycleESDs();
566   IncEvCountTotalESDs();
567 }
568
569 //____________________________________________________________________________
570 void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
571 {
572   // Check that all the reference histograms exist before we try to use them - otherwise call InitRaws
573   // RS: Attention: the counters are increments after custom modification of eventSpecie
574   if (!fCalibRefHistoPro || !fCalibRefHistoH2F || !fLEDMonRefHistoPro || !fHighEmcHistoH2F) {
575     InitRaws();
576   }
577
578   // make sure EMCal was readout during the event
579   Int_t emcID = AliDAQ::DetectorID("EMCAL"); // bit 18..
580   const UInt_t *detPattern = rawReader->GetDetectorPattern(); 
581   UInt_t emcInReadout = ( ((1 << emcID) & detPattern[0]) >> emcID);
582   if (! emcInReadout) return; // no poInt_t in looking at this event, if no EMCal data
583
584   // setup
585   rawReader->Reset() ;
586   AliCaloRawStreamV3 in(rawReader,"EMCAL"); 
587   rawReader->Select("EMCAL", 0, AliEMCALGeoParams::fgkLastAltroDDL) ; //select EMCAL DDL's 
588
589   AliRecoParam::EventSpecie_t saveSpecie = fEventSpecie ;
590   if (rawReader->GetType() == AliRawEventHeaderBase::kCalibrationEvent) { 
591     SetEventSpecie(AliRecoParam::kCalib) ;      
592   }
593
594   const Int_t nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
595   const Int_t nRows        = AliEMCALGeoParams::fgkEMCALRows; // number of rows per SuperModule
596   const Int_t nStripsPerSM = AliEMCALGeoParams::fgkEMCALLEDRefs; // number of strips per SuperModule
597   const Int_t n2x2PerSM    = AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // number of TRU 2x2's per SuperModule
598   const Int_t n2x2PerTRU   = AliEMCALGeoParams::fgkEMCAL2x2PerTRU;
599   const Int_t nTot2x2      = fSuperModules * n2x2PerSM; // total TRU channel
600
601   // SM counters; decl. should be safe, assuming we don't get more than expected SuperModules..
602   Int_t nTotalSMLG[AliEMCALGeoParams::fgkEMCALModules]       = {0};
603   Int_t nTotalSMHG[AliEMCALGeoParams::fgkEMCALModules]       = {0};
604   Int_t nTotalSMTRU[AliEMCALGeoParams::fgkEMCALModules]      = {0};
605   Int_t nTotalSMLGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
606   Int_t nTotalSMHGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
607
608   const Int_t nTRUL0ChannelBits = 10; // used for L0 trigger bits checks
609         int firstL0TimeBin = 999;
610   int triggers[nTot2x2][24]; //auxiliary array for L0 trigger - TODO remove hardcoded 24
611         memset(triggers, 0, sizeof(int) * 24 * nTot2x2);
612
613   Int_t iSM = 0; // SuperModule index 
614   // start loop over input stream  
615   while (in.NextDDL()) {
616     Int_t iRCU = in.GetDDLNumber() % 2; // RCU0 or RCU1, within SuperModule
617     Int_t iDDL = in.GetDDLNumber();
618     fRawAnalyzer->SetIsZeroSuppressed( in.GetZeroSupp() ); 
619     
620     while (in.NextChannel()) {
621       Int_t iBranch = in.GetBranch();
622       
623       iSM = in.GetModule(); // SuperModule
624       //prInt_tf("iSM %d DDL %d", iSM, in.GetDDLNumber()); 
625       if (iSM>=0 && iSM<fSuperModules) { // valid module reading
626
627         Int_t nsamples = 0;
628         vector<AliCaloBunchInfo> bunchlist; 
629         while (in.NextBunch()) {
630           nsamples += in.GetBunchLength();
631           bunchlist.push_back( AliCaloBunchInfo(in.GetStartTimeBin(), in.GetBunchLength(), in.GetSignals() ) );
632         } 
633         
634         if (nsamples > 0) { // this check is needed for when we have zero-supp. on, but not sparse readout
635           Float_t time = 0.; 
636           Float_t amp  = 0.; 
637           // indices for pedestal calc.
638           Int_t firstPedSample = 0;
639           Int_t lastPedSample  = 0;
640           bool isTRUL0IdData   = false;
641
642           if (! in.IsTRUData() ) { // high gain, low gain, LED Mon data - all have the same shaper/sampling 
643             AliCaloFitResults fitResults = fRawAnalyzer->Evaluate( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2()); 
644             amp  = fitResults.GetAmp();
645             time = fitResults.GetTof(); 
646             firstPedSample = fFirstPedestalSample;
647             lastPedSample  = fLastPedestalSample;
648           }
649           else { // TRU data is special, needs its own analyzer
650             AliCaloFitResults fitResults = fRawAnalyzerTRU->Evaluate( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2()); 
651             amp  = fitResults.GetAmp();
652             time = fitResults.GetTof(); 
653             firstPedSample = fFirstPedestalSampleTRU;
654             lastPedSample  = fLastPedestalSampleTRU;
655             if (in.GetColumn() >= n2x2PerTRU) {
656               isTRUL0IdData = true;
657             }
658           }
659   
660           // pedestal samples
661           Int_t nPed = 0;
662           vector<Int_t> pedSamples; 
663         
664           // select earliest bunch 
665           unsigned int bunchIndex = 0;
666           unsigned int startBin = bunchlist.at(0).GetStartBin();
667           if (bunchlist.size() > 0) {
668             for(unsigned int ui=1; ui < bunchlist.size(); ui++ ) {
669               if (startBin > bunchlist.at(ui).GetStartBin() ) {
670                 startBin = bunchlist.at(ui).GetStartBin();
671                 bunchIndex = ui;
672               }
673             }
674           }
675
676           // check bunch for entries in the pedestal sample range
677           Int_t bunchLength = bunchlist.at(bunchIndex).GetLength(); 
678           const UShort_t *sig = bunchlist.at(bunchIndex).GetData();
679           Int_t timebin = 0;
680
681           if (! isTRUL0IdData) { // regular data, can look at pedestals
682             for (Int_t i = 0; i<bunchLength; i++) {
683               timebin = startBin--;
684               if ( firstPedSample<=timebin && timebin<=lastPedSample ) {
685                 pedSamples.push_back( sig[i] );
686                 nPed++;
687               }     
688             } // i
689           }
690           else { // TRU L0 Id Data
691             // which TRU the channel belongs to?
692             Int_t iTRUId = in.GetModule()*3 + (iRCU*in.GetBranch() + iRCU);
693
694             for (Int_t i = 0; i< bunchLength; i++) {
695               for( Int_t j = 0; j < nTRUL0ChannelBits; j++ ){
696                 // check if the bit j is 1
697                 if( (sig[i] & ( 1 << j )) > 0 ){
698                   Int_t iTRUIdInSM = (in.GetColumn() - n2x2PerTRU)*nTRUL0ChannelBits+j;
699                   if(iTRUIdInSM < n2x2PerTRU) {
700                     Int_t iTRUAbsId = iTRUIdInSM + n2x2PerTRU * iTRUId;
701                     // Fill the histograms
702                     Int_t globTRUCol, globTRURow;
703                     GetTruChannelPosition(globTRURow, globTRUCol, iSM, iDDL, iBranch, iTRUIdInSM );
704                     
705                     FillRawsData(kNL0TRU, globTRUCol, globTRURow);
706                     FillRawsData(kTimeL0TRU, globTRUCol, globTRURow, startBin);
707                     triggers[iTRUAbsId][startBin] = 1;
708                     
709                     if((int)startBin < firstL0TimeBin) firstL0TimeBin = startBin;
710                   }
711                 }
712               }
713               startBin--;
714             } // i      
715           } // TRU L0 Id data                   
716           
717           // fill histograms
718           if ( in.IsLowGain() || in.IsHighGain() ) { // regular towers
719             Int_t towerId = iSM*nTowersPerSM + in.GetColumn()*nRows + in.GetRow();
720             if ( in.IsLowGain() ) { 
721               nTotalSMLG[iSM]++; 
722               if ( (amp > fMinSignalLG) && (amp < fMaxSignalLG) ) { 
723                 FillRawsData(kSigLG,towerId, amp);
724                 FillRawsData(kTimeLG,towerId, time);
725               }
726               if (nPed > 0) {
727                 for (Int_t i=0; i<nPed; i++) {
728                   FillRawsData(kPedLG,towerId, pedSamples[i]);
729                 }
730               }
731             } // gain==0
732             else if ( in.IsHighGain() ) {               
733               nTotalSMHG[iSM]++; 
734               if ( (amp > fMinSignalHG) && (amp < fMaxSignalHG) ) { 
735                 FillRawsData(kSigHG,towerId, amp);
736                 FillRawsData(kTimeHG,towerId, time);
737               } 
738               if (nPed > 0) {
739                 for (Int_t i=0; i<nPed; i++) {
740                   FillRawsData(kPedHG,towerId, pedSamples[i]);
741                 }
742               }
743             } // gain==1
744           } // low or high gain
745           // TRU
746           else if ( in.IsTRUData() && in.GetColumn()<AliEMCALGeoParams::fgkEMCAL2x2PerTRU) {
747             // for TRU data, the mapping class holds the TRU Int_ternal 2x2 number (0..95) in the Column var..
748             Int_t iTRU = (iRCU*in.GetBranch() + iRCU); //TRU0 is from RCU0, TRU1 from RCU1, TRU2 is from branch B on RCU1
749             Int_t iTRU2x2Id = iSM*n2x2PerSM + iTRU*AliEMCALGeoParams::fgkEMCAL2x2PerTRU 
750               + in.GetColumn();
751             nTotalSMTRU[iSM]++; 
752             if ( (amp > fMinSignalTRU) && (amp < fMaxSignalTRU) ) { 
753               FillRawsData(kSigTRU,iTRU2x2Id, amp);
754               //FillRawsData(kTimeTRU,iTRU2x2Id, time);
755             }
756             //if (nPed > 0) {
757               //for (Int_t i=0; i<nPed; i++) {
758                 //FillRawsData(kPedTRU,iTRU2x2Id, pedSamples[i]);
759               //}
760             //}
761           }
762           // LED Mon
763           else if ( in.IsLEDMonData() ) {
764             // for LED Mon data, the mapping class holds the gain info in the Row variable
765             // and the Strip number in the Column..
766             Int_t gain = in.GetRow(); 
767             Int_t stripId = iSM*nStripsPerSM + in.GetColumn();
768           
769             if ( gain == 0 ) { 
770               nTotalSMLGLEDMon[iSM]++; 
771               if ( (amp > fMinSignalLGLEDMon) && (amp < fMaxSignalLGLEDMon) ) {
772                 FillRawsData(kSigLGLEDMon,stripId, amp);
773                 FillRawsData(kTimeLGLEDMon,stripId, time);
774               }
775               if (nPed > 0) {
776                 for (Int_t i=0; i<nPed; i++) {
777                   FillRawsData(kPedLGLEDMon,stripId, pedSamples[i]);
778                 }
779               }
780             } // gain==0
781             else if ( gain == 1 ) {             
782               nTotalSMHGLEDMon[iSM]++; 
783               if ( (amp > fMinSignalHGLEDMon) && (amp < fMaxSignalHGLEDMon) ) { 
784                 FillRawsData(kSigHGLEDMon,stripId, amp);
785                 FillRawsData(kTimeHGLEDMon,stripId, time);
786               }
787               if (nPed > 0) {
788                 for (Int_t i=0; i<nPed; i++) {
789                   FillRawsData(kPedHGLEDMon,stripId, pedSamples[i]);
790                 }
791               }
792             } // low or high gain
793           } // LEDMon
794
795         } // SM index OK
796
797       } // nsamples>0 check, some data found for this channel; not only trailer/header
798     }// end while over channel 
799    
800   }//end while over DDL's, of input stream 
801   //filling some L0 trigger histos
802   if( firstL0TimeBin < 999 ){
803     for(Int_t i = 0; i < nTot2x2; i++) {        
804       if( triggers[i][firstL0TimeBin] > 0 ) {
805         //histo->Fill(i,j);
806         FillRawsData(kNL0FirstTRU, i);
807         FillRawsData(kTimeL0FirstTRU, i, firstL0TimeBin);
808       }
809     }
810   }
811   
812   //calculate the ratio of the amplitude and fill the histograms, only if the events type is Calib
813   // RS: operation on the group of histos kSigHG,k2DRatioAmp,kRatioDist,kLEDMonRatio,kLEDMonRatio,kSigLGLEDMon
814   const int hGrp[] = {kSigHG,k2DRatioAmp,kRatioDist,kLEDMonRatio,kLEDMonRatioDist,kSigLGLEDMon};
815   if ( rawReader->GetType() == AliRawEventHeaderBase::kCalibrationEvent &&
816        CheckCloningConsistency(fRawsQAList, hGrp, sizeof(hGrp)/sizeof(int)) ) {  // RS converting original code to loop over all matching triggers
817     int nTrig =IsClonedPerTrigClass(kSigHG,fRawsQAList) ? GetNEventTrigClasses() : 0; // loop over triggers only if histos were cloned
818     //
819     for (int itr=-1;itr<nTrig;itr++) { // start from -1 to acknowledge original histos if they were kept
820       TObjArray* trArr = GetMatchingRawsHistosSet(hGrp, sizeof(hGrp)/sizeof(int) ,itr);
821       if (!trArr) continue;  // no histos for current trigger
822       //
823       Double_t binContent = 0.;
824       TProfile* prSigHG      = (TProfile *)trArr->At(0); //kSigHG
825       TH1* th2DRatioAmp      = (TH1*) trArr->At(1); //k2DRatioAmp
826       TH1* thRatioDist       = (TH1*) trArr->At(2); //kRatioDist
827       TH1* thLEDMonRatio     = (TH1*) trArr->At(3); //kLEDMonRatio
828       TH1* thLEDMonRatioDist = (TH1*) trArr->At(4); //kLEDMonRatio
829       TH1* hSigLGLEDMon      = (TH1*) trArr->At(5); //kSigLGLEDMon
830       th2DRatioAmp->Reset("ICE");
831       thRatioDist->Reset("ICE");
832       thLEDMonRatio->Reset("ICE");
833       thLEDMonRatioDist->Reset("ICE");
834       th2DRatioAmp->ResetStats();
835       thRatioDist->ResetStats();
836       thLEDMonRatio->ResetStats();
837       thLEDMonRatioDist->ResetStats();
838       ConvertProfile2H(prSigHG, fHighEmcHistoH2F);  
839       //
840       for(Int_t ix = 1; ix <= fHighEmcHistoH2F->GetNbinsX(); ix++) {
841         for(Int_t iy = 1; iy <= fHighEmcHistoH2F->GetNbinsY(); iy++) { 
842           if(fCalibRefHistoH2F->GetBinContent(ix, iy)) 
843             binContent = fHighEmcHistoH2F->GetBinContent(ix, iy)/fCalibRefHistoH2F->GetBinContent(ix, iy);
844           th2DRatioAmp->SetBinContent(ix, iy, binContent);
845           thRatioDist->Fill(binContent);
846         }
847       } 
848       //
849       //Now for LED monitor system, to calculate the ratio as well
850       Double_t binError = 0. ;
851       // for the binError, we add the relative errors, squared
852       Double_t relativeErrorSqr = 0. ;
853       //
854       for(int ib = 1; ib <= fLEDMonRefHistoPro->GetNbinsX(); ib++) {
855         //
856         if(fLEDMonRefHistoPro->GetBinContent(ib) != 0) {
857           binContent = hSigLGLEDMon->GetBinContent(ib) / fLEDMonRefHistoPro->GetBinContent(ib);
858
859           relativeErrorSqr = TMath::Power( (fLEDMonRefHistoPro->GetBinError(ib) / fLEDMonRefHistoPro->GetBinContent(ib)), 2);
860           if( hSigLGLEDMon->GetBinContent(ib) != 0) {
861             relativeErrorSqr += TMath::Power( (hSigLGLEDMon->GetBinError(ib)/hSigLGLEDMon->GetBinContent(ib)), 2);
862           }
863         }
864         else { // ref. run info is zero
865           binContent = -1;
866           relativeErrorSqr = 1;
867         }
868         thLEDMonRatio->SetBinContent(ib, binContent);
869         
870         binError = sqrt(relativeErrorSqr) * binContent;
871         thLEDMonRatio->SetBinError(ib, binError);
872         thLEDMonRatioDist->Fill(thLEDMonRatio->GetBinContent(ib));
873       }
874     } // loop over eventual trigger clones
875   } 
876   // let's also fill the SM and event counter histograms
877   Int_t nTotalHG = 0;
878   Int_t nTotalLG = 0;
879   Int_t nTotalTRU = 0;
880   Int_t nTotalHGLEDMon = 0;
881   Int_t nTotalLGLEDMon = 0;
882   for (iSM=0; iSM<fSuperModules; iSM++) {  
883     nTotalLG += nTotalSMLG[iSM]; 
884     nTotalHG += nTotalSMHG[iSM]; 
885     nTotalTRU += nTotalSMTRU[iSM]; 
886     nTotalLGLEDMon += nTotalSMLGLEDMon[iSM]; 
887     nTotalHGLEDMon += nTotalSMHGLEDMon[iSM]; 
888     FillRawsData(kNsmodLG,iSM, nTotalSMLG[iSM]); 
889     FillRawsData(kNsmodHG,iSM, nTotalSMHG[iSM]); 
890     FillRawsData(kNsmodTRU,iSM, nTotalSMTRU[iSM]); 
891     FillRawsData(kNsmodLGLEDMon,iSM, nTotalSMLGLEDMon[iSM]); 
892     FillRawsData(kNsmodHGLEDMon,iSM, nTotalSMHGLEDMon[iSM]); 
893   }
894  
895   FillRawsData(kNtotLG,nTotalLG);
896   FillRawsData(kNtotHG,nTotalHG);
897   FillRawsData(kNtotTRU,nTotalTRU);
898   FillRawsData(kNtotLGLEDMon,nTotalLGLEDMon);
899   FillRawsData(kNtotHGLEDMon,nTotalHGLEDMon);
900  
901   IncEvCountCycleESDs();
902   IncEvCountTotalESDs();
903   SetEventSpecie(saveSpecie) ; 
904   
905         MakeRawsSTU(rawReader);
906
907   // just in case the next rawreader consumer forgets to reset; let's do it here again..
908   rawReader->Reset() ;
909   return;
910 }
911
912 //____________________________________________________________________________
913 void AliEMCALQADataMakerRec::MakeDigits()
914 {
915   // makes data from Digits
916   FillDigitsData(1,fDigitsArray->GetEntriesFast()) ; 
917   TIter next(fDigitsArray) ; 
918   AliEMCALDigit * digit ; 
919   while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) {
920     FillDigitsData(0, digit->GetAmplitude()) ;
921   }  
922   //
923 }
924
925 //____________________________________________________________________________
926 void AliEMCALQADataMakerRec::MakeDigits(TTree * digitTree)
927 {
928   // makes data from Digit Tree
929   // RS: Attention: the counters are increments in the MakeDigits()
930   if (fDigitsArray) 
931     fDigitsArray->Clear("C") ; 
932   else
933     fDigitsArray = new TClonesArray("AliEMCALDigit", 1000) ; 
934   
935   TBranch * branch = digitTree->GetBranch("EMCAL") ;
936   if ( ! branch ) { AliWarning("EMCAL branch in Digit Tree not found"); return; }
937   //
938   branch->SetAddress(&fDigitsArray) ;
939   branch->GetEntry(0) ; 
940   MakeDigits() ; 
941   //
942   IncEvCountCycleDigits();
943   IncEvCountTotalDigits();  
944   //  
945 }
946
947 //____________________________________________________________________________
948 void AliEMCALQADataMakerRec::MakeRecPoints(TTree * clustersTree)
949 {
950   // makes data from RecPoints
951   TBranch *emcbranch = clustersTree->GetBranch("EMCALECARP");
952   if (!emcbranch) { 
953     AliError("can't get the branch with the EMCAL clusters !");
954     return;
955   }
956   
957   TObjArray * emcRecPoints = new TObjArray(100) ;
958   emcbranch->SetAddress(&emcRecPoints);
959   emcbranch->GetEntry(0);
960   
961   FillRecPointsData(kRecPM,emcRecPoints->GetEntriesFast()) ; 
962   TIter next(emcRecPoints) ; 
963   AliEMCALRecPoint * rp ; 
964   while ( (rp = dynamic_cast<AliEMCALRecPoint *>(next())) ) {
965     FillRecPointsData(kRecPE,rp->GetEnergy()) ;
966     FillRecPointsData(kRecPDigM,rp->GetMultiplicity());
967   }
968   emcRecPoints->Delete();
969   delete emcRecPoints;
970   IncEvCountCycleRecPoints();
971   IncEvCountTotalRecPoints();
972 }
973
974 //____________________________________________________________________________ 
975 void AliEMCALQADataMakerRec::StartOfDetectorCycle()
976 {
977   //Detector specific actions at start of cycle
978   
979 }
980
981 //____________________________________________________________________________ 
982 void AliEMCALQADataMakerRec::SetFittingAlgorithm(Int_t fitAlgo)              
983 {
984   //Set fitting algorithm and initialize it if this same algorithm was not set before.
985   //printf("**** Set Algorithm , number %d ****\n",fitAlgo);
986
987   
988   fRawAnalyzer =  AliCaloRawAnalyzerFactory::CreateAnalyzer(fitAlgo);
989   fFittingAlgorithm = fitAlgo; 
990
991   /*
992   if(fitAlgo == fFittingAlgorithm && fRawAnalyzer) {
993     //Do nothing, this same algorithm already set before.
994     //printf("**** Algorithm already set before, number %d, %s ****\n",fitAlgo, fRawAnalyzer->GetName());
995     return;
996   }
997   //Initialize the requested algorithm
998   if(fitAlgo != fFittingAlgorithm || !fRawAnalyzer) {
999     //printf("**** Init Algorithm , number %d ****\n",fitAlgo);
1000                 
1001     fFittingAlgorithm = fitAlgo; 
1002     if (fRawAnalyzer) delete fRawAnalyzer;  // delete prev. analyzer if existed.
1003                 
1004     if (fitAlgo == kFastFit) {
1005       fRawAnalyzer = new AliCaloRawAnalyzerFastFit();
1006     }
1007     else if (fitAlgo == kNeuralNet) {
1008       fRawAnalyzer = new AliCaloRawAnalyzerNN();
1009     }
1010     else if (fitAlgo == kLMS) {
1011       fRawAnalyzer = new AliCaloRawAnalyzerLMS();
1012     }
1013     else if (fitAlgo == kPeakFinder) {
1014       fRawAnalyzer = new AliCaloRawAnalyzerPeakFinder();
1015     }
1016     else if (fitAlgo == kCrude) {
1017       fRawAnalyzer = new AliCaloRawAnalyzerCrude();
1018     }
1019     else {
1020       AliWarning("EMCAL QA invalid fit algorithm choice") ; 
1021     }
1022
1023   }
1024   return;
1025   */
1026 }
1027
1028 //_____________________________________________________________________________________
1029 void AliEMCALQADataMakerRec::ConvertProfile2H(TProfile * p, TH2 * histo) 
1030 {  
1031   // reset histogram
1032   histo->Reset("ICE") ; 
1033   histo->ResetStats(); 
1034
1035   Int_t nbinsProf = p->GetNbinsX();
1036   
1037   // loop through the TProfile p and fill the TH2F histo 
1038   Int_t row = 0;
1039   Int_t col = 0;
1040   Double_t binContent = 0;
1041   Int_t towerNum = 0; // global tower Id
1042   //  i = 0; // tower Id within SuperModule
1043   Int_t iSM = 0; // SuperModule index 
1044   Int_t iSMSide = 0; // 0=A, 1=C side
1045   Int_t iSMSector = 0; // 2 SM's per sector
1046   
1047   // indices for 2D plots
1048   Int_t col2d = 0;
1049   Int_t row2d = 0;
1050   
1051   for (Int_t ibin = 1; ibin <= nbinsProf; ibin++) {
1052     towerNum = (Int_t) p->GetBinCenter(ibin);
1053     binContent = p->GetBinContent(ibin);
1054     
1055     // figure out what the tower indices are: col, row within a SuperModule
1056     iSM = towerNum/(AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols);
1057     col = (towerNum/AliEMCALGeoParams::fgkEMCALRows) % (AliEMCALGeoParams::fgkEMCALCols);
1058     row = towerNum % (AliEMCALGeoParams::fgkEMCALRows);
1059     
1060     //DecodeTowerNum(towerNum, &SM, &col, &row);
1061     // then we calculate what the global 2D coord are, based on which SM 
1062     // we are in
1063     iSMSector = iSM / 2;
1064     iSMSide = iSM % 2;
1065     
1066     if (iSMSide == 1) { // C side, shown to the right
1067       col2d = col + AliEMCALGeoParams::fgkEMCALCols;
1068     }
1069     else { // A side, shown to the left 
1070       col2d = col; 
1071     }
1072     
1073     row2d = row + iSMSector * AliEMCALGeoParams::fgkEMCALRows;
1074     
1075     histo->SetBinContent(col2d+1, row2d+1, binContent);
1076   }
1077
1078 //____________________________________________________________________________ 
1079 void AliEMCALQADataMakerRec::GetTruChannelPosition( Int_t &globRow, Int_t &globColumn, Int_t module, Int_t ddl, Int_t branch, Int_t column )
1080 {
1081   Int_t mrow;
1082   Int_t mcol;
1083   Int_t trow;
1084   Int_t tcol;
1085   Int_t drow;
1086   Int_t rcu;
1087   // RCU 0 or 1
1088   rcu = ddl % 2;
1089
1090   // 12 rows of 2x2s in a module (3 TRUs by 4 rows)
1091   mrow = (module/2) * 12;
1092   // 24 columns per module, odd module numbers increased by 24
1093   mcol = (module%2) * 24;
1094
1095   // position within TRU coordinates
1096   tcol = column / 4;
1097   trow = column % 4;
1098
1099   //.combine
1100   if( module%2 == 0 ){   // A side
1101     // mirror rows
1102     trow = 3 - trow;
1103
1104     // TRU in module row addition
1105     drow = (rcu*branch+rcu) * 4;
1106
1107   }
1108   else{   // C side
1109     // mirror columns
1110     tcol = 23 - tcol;
1111
1112     // TRU in module row addition
1113     drow = (2 - (rcu*branch+rcu)) * 4;
1114   }
1115
1116   // output global row/collumn position (0,0 = SMA0, phi = 0, |eta| = max)
1117   globRow = mrow + drow + trow;
1118   globColumn = mcol + tcol;
1119         return;
1120
1121 }
1122 //____________________________________________________________________________ 
1123 void AliEMCALQADataMakerRec::MakeRawsSTU(AliRawReader* rawReader){
1124
1125         AliEMCALTriggerSTURawStream* inSTU = new AliEMCALTriggerSTURawStream(rawReader);
1126         
1127         rawReader->Reset();
1128         rawReader->Select("EMCAL", 44);
1129
1130         
1131         //L1 segmentation
1132         Int_t sizeL1gsubr = 1;
1133         Int_t sizeL1gpatch = 2; 
1134         Int_t sizeL1jsubr = 4; 
1135
1136      Int_t EMCALtrig[AliEMCALGeoParams::fgkEMCALSTUCols][AliEMCALGeoParams::fgkEMCALSTURows];
1137
1138                 memset(EMCALtrig, 0, sizeof(int) * AliEMCALGeoParams::fgkEMCALSTUCols * AliEMCALGeoParams::fgkEMCALSTURows);
1139                 
1140
1141         
1142                 
1143      if (inSTU->ReadPayLoad()) 
1144         {
1145                         
1146           //Fw version (use in case of change in L1 jet 
1147           Int_t fw = inSTU->GetFwVersion();
1148           Int_t sizeL1jpatch = 2+(fw >> 16);
1149
1150           //To check link
1151         
1152           Int_t mask = inSTU->GetFrameReceived();
1153
1154
1155           for (int i = 0; i < 32; i++)
1156             {
1157               if ((mask >> i) & 0x1) FillRawsData(kSTUTRU, i);
1158             }
1159
1160
1161           //V0 signal in STU
1162           Int_t V0Sig = inSTU->GetV0A()+inSTU->GetV0C();
1163                         
1164           //FastOR amplitude receive from TRU
1165           for (Int_t i = 0; i < 32; i++)
1166             {
1167               UInt_t adc[96];
1168               for (Int_t j = 0; j < 96; j++) adc[j] = 0;
1169                                 
1170               inSTU->GetADC(i, adc);
1171                                 
1172               Int_t iTRU = fGeom->GetTRUIndexFromSTUIndex(i);
1173                                 
1174               for (Int_t j = 0; j < 96; j++)
1175                 {
1176                   Int_t idx;
1177                   fGeom->GetAbsFastORIndexFromTRU(iTRU, j, idx);
1178                                 
1179                   Int_t px, py;
1180                   fGeom->GetPositionInEMCALFromAbsFastORIndex(idx, px, py);
1181                                         
1182                   EMCALtrig[px][py] = adc[j];
1183                 }
1184             }
1185                         
1186           //L1 Gamma patches
1187           Int_t iTRU_STU, x, y;
1188           for (Int_t i = 0; i < inSTU->GetNL1GammaPatch(); i++)
1189             {
1190               if (inSTU->GetL1GammaPatch(i, iTRU_STU, x, y)) // col (0..23), row (0..3)
1191                 {
1192                   Int_t iTRU;
1193                   iTRU = fGeom->GetTRUIndexFromSTUIndex(iTRU_STU);
1194                                         
1195                   Int_t etaG = 23-x, phiG = y + 4 * int(iTRU/2); //position in EMCal
1196                                         
1197                   if (iTRU%2) etaG += 24; //C-side
1198                                         
1199                   etaG = etaG - sizeL1gsubr * sizeL1gpatch + 1;
1200                                 
1201                   //Position of patch L1G (bottom-left FastOR of the patch)
1202                         FillRawsData(kGL1, etaG, phiG);
1203                                         
1204                   //loop to sum amplitude of FOR in the gamma patch
1205                   Int_t L1G_PatchAmp = 0;
1206                   for (Int_t L1Gx = 0; L1Gx < sizeL1gpatch; L1Gx ++)
1207                     {
1208                       for (Int_t L1Gy = 0; L1Gy < sizeL1gpatch; L1Gy ++)
1209                         {
1210                           if (etaG+L1Gx < 48 && phiG+L1Gy < 64) L1G_PatchAmp += EMCALtrig[etaG+L1Gx][phiG+L1Gy];
1211                           //cout << EMCALtrig[etaG+L1Gx][phiG+L1Gy] << endl;
1212                         }
1213                     }
1214                         
1215                   //if (L1G_PatchAmp > 500) cout << "L1G amp =" << L1G_PatchAmp << endl;
1216                         FillRawsData(kGL1V0, V0Sig, L1G_PatchAmp);
1217                                         
1218                 }
1219             }
1220                 
1221                         
1222           //L1 Jet patches
1223           for (Int_t i = 0; i < inSTU->GetNL1JetPatch(); i++)
1224             {
1225               if (inSTU->GetL1JetPatch(i, x, y)) // col (0,15), row (0,11)
1226                 {
1227                                 
1228                   Int_t etaJ = sizeL1jsubr * (11-y-sizeL1jpatch + 1);
1229                   Int_t phiJ = sizeL1jsubr * (15-x-sizeL1jpatch + 1);
1230                                         
1231                   //position of patch L1J (FOR bottom-left)
1232                         FillRawsData(kJL1, etaJ, phiJ);
1233                                         
1234                   //loop the sum aplitude of FOR in the jet patch
1235                   Int_t L1J_PatchAmp = 0;
1236                   for (Int_t L1Jx = 0; L1Jx < sizeL1jpatch*4; L1Jx ++)
1237                     {
1238                       for (Int_t L1Jy = 0; L1Jy < sizeL1jpatch*4; L1Jy ++)
1239                         {
1240                           if (etaJ+L1Jx < 48 && phiJ+L1Jy < 64) L1J_PatchAmp += EMCALtrig[etaJ+L1Jx][phiJ+L1Jy];
1241                         }
1242                     }
1243                 
1244                   //cout << "L1J amp =" << L1J_PatchAmp << endl;
1245                         FillRawsData(kJL1V0, V0Sig, L1J_PatchAmp);
1246                                         
1247                 }
1248             }
1249
1250         }
1251                 
1252      //Fill FOR amplitude histo
1253      for (Int_t i = 0; i < 48; i++)
1254         {
1255           for (Int_t j = 0; j < 60; j++)
1256             {
1257               if (EMCALtrig[i][j] != 0) FillRawsData(kAmpL1, i, j, EMCALtrig[i][j]);
1258             }
1259         }
1260                 
1261         delete inSTU;
1262         return;
1263 }
1264
1265