]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALFolder.cxx
added pi0 calibration, linearity, shower profile
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALFolder.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, 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 /* $Id$ */
17
18 //_________________________________________________________________________
19 // Top EMCAL folder which will keep all information about EMCAL itself,
20 // super Modules (SM), modules, towers, set of hists and so on.
21 //
22 //*-- Author: Aleksei Pavlinov (WSU, Detroit, USA) 
23
24 #include "AliEMCALFolder.h"
25 #include "AliEMCALHistoUtilities.h"
26 #include "AliEMCALGeometry.h"
27 #include "AliEMCALSuperModule.h"
28 #include "AliEMCALCell.h"
29 #include "AliESDCaloCluster.h"
30
31 #include "AliRun.h"
32
33 #include "AliEMCALCalibData.h"
34 #include "AliCDBMetaData.h"
35 #include "AliCDBId.h"
36 #include "AliCDBEntry.h"
37 #include "AliCDBManager.h"
38 #include "AliCDBStorage.h"
39
40 #include "AliEMCALCalibCoefs.h"
41 #include "AliEMCALDigit.h"
42 #include "AliEMCALRecPoint.h"
43
44 #include "AliEMCALPi0SelectionParam.h"
45
46 #include <cassert>
47
48 #include <TROOT.h>
49 #include <TStyle.h>
50 #include <TList.h>
51 #include <TH1.h>
52 #include <TF1.h>
53 #include <TFile.h>
54 #include <TCanvas.h>
55 #include <TClonesArray.h>
56 #include <TKey.h>
57 #include <TNtuple.h>
58 #include <TLegend.h>
59 #include <TLegendEntry.h>
60 #include <TLine.h>
61
62   const TString AliEMCALFolder::fgkBaseFolderName("EMCAL");
63   const TString AliEMCALFolder::fgkCCFirstName("CCFIRST");
64   const TString AliEMCALFolder::fgkCCinName("CCIN");
65   const TString AliEMCALFolder::fgkCCoutName("CCOUT");
66   const TString AliEMCALFolder::fgkDirOfRootFiles("$HOME/ALICE/SHISHKEBAB/RF/CALIB/JUL16/");
67
68 typedef  AliEMCALHistoUtilities u;
69
70 ClassImp(AliEMCALFolder)
71
72 //AliEMCALGeometry* AliEMCALFolder::fGeometry = 0;
73
74   TList *lObj = 0;
75
76 AliEMCALFolder::AliEMCALFolder() : 
77   TFolder(), 
78   fCounter(0), fGeometry(0), fNumOfCell(0), fLhists(0), fLofCells(0),fPi0SelPar(0),fCalibData(0),
79   fCellNtuple(0)
80 {
81 }
82
83 AliEMCALFolder::AliEMCALFolder(const char* name, const char* title, Bool_t putToBrowser) : 
84   TFolder(name,title),
85   fCounter(-1), fGeometry(0), fNumOfCell(0), fLhists(0), fLofCells(0),fPi0SelPar(0),fCalibData(0),
86   fCellNtuple(0)
87 {
88   Init(putToBrowser);
89 }
90
91 AliEMCALFolder::AliEMCALFolder(const Int_t it, const char* title, Bool_t putToBrowser) : 
92   TFolder(Form("%s_%2.2i", AliEMCALFolder::fgkBaseFolderName.Data(),it),title),
93   fCounter(it), fGeometry(0), fNumOfCell(0), fLhists(0), fLofCells(0),fPi0SelPar(0),fCalibData(0),
94   fCellNtuple(0)
95 {
96   Init(putToBrowser);
97 }
98
99 AliEMCALFolder::~AliEMCALFolder()
100 {
101   // dtor
102 }
103
104 void AliEMCALFolder::Init(Bool_t putToBrowser)
105 {
106   lObj = new TList;
107   lObj->SetName("Objects"); // name is good ?
108   this->Add((TObject*)lObj);
109   //this->AddObject((TObject*)lObj, kTRUE);
110   // Get default geometry - "SHISH_77_TRD1_2X2_FINAL_110DEG"; May 29, 2007
111   fGeometry = AliEMCALGeometry::GetInstance(); // should be define before 
112   lObj->Add(fGeometry);
113
114   // Initial cc with decalibration
115   // Jun 13, 2007;
116   // Jul 13 - See ~/macros/ALICE/sim.C  for choice of CDB
117   AliEMCALCalibData *calData[1];
118   // Get from defined arrea;
119   // First table is table which used in rec.points finder.
120   Add(AliEMCALCalibCoefs::GetCalibTableFromDb(fgkCCFirstName.Data(),calData));
121   fCalibData = calData[0];
122   lObj->Add(fCalibData);
123   if(GetIterationNumber()<=1) {
124     Add(AliEMCALCalibCoefs::GetCalibTableFromDb(fgkCCinName.Data(), calData));
125   }
126
127   // Selection Parameter
128   fPi0SelPar = AliEMCALPi0SelectionParam::Set1();
129   this->Add(fPi0SelPar);
130   //
131   fLhists = BookHists();
132   lObj->Add(fLhists);
133
134   // dimension should be get from geometry - 4*12*24*11); 
135   fNumOfCell = fGeometry->GetNCells();
136   fLofCells  = new AliEMCALCell*[fNumOfCell];
137   for(int i=0; i<fNumOfCell; i++) fLofCells[i] = 0;
138
139   printf("<I> Create AliEMCALFolder : it %i : name %s\n ", fCounter, GetName());
140
141   if(putToBrowser) gROOT->GetListOfBrowsables()->Add(this); // for testing purpuse
142
143
144 AliEMCALSuperModule* AliEMCALFolder::GetSuperModule(const Int_t nm)
145 {
146   AliEMCALSuperModule *sm = 0;
147
148   TObject *set = FindObject(Form("SM%2.2i",nm));
149   if(set) sm = (AliEMCALSuperModule*)set;
150
151   return sm;
152 }
153
154
155 AliEMCALCell* AliEMCALFolder::GetCell(const Int_t absId)
156 { // May 30, 2007
157   if(absId<0 || absId >= fNumOfCell) return 0;
158   else return fLofCells[absId];
159 }  
160
161 void  AliEMCALFolder::SetCell(AliEMCALCell *cell, const Int_t absId) 
162 {
163   if(absId>=0 && absId < fNumOfCell) {
164     fLofCells[absId] = cell;
165   }
166 }  
167
168 pi0SelectionParam* AliEMCALFolder::GetPi0SelectionParRow(Int_t nrow)
169 {
170   pi0SelectionParam* r=0;
171   if(fPi0SelPar) {
172     r = fPi0SelPar->GetTable(nrow);
173   }
174   return r;
175 }
176
177 void AliEMCALFolder::FillPi0Candidate(const Double_t mgg, AliESDCaloCluster* cl1, AliESDCaloCluster* cl2)
178 {
179   static Int_t absIdMax, nm1, nm2;
180
181   u::FillH1(fLhists, 0, 1.);  // number entries 
182   u::FillH1(fLhists, 1, mgg);
183
184   nm1 = GetSMNumber(cl1);
185   nm2 = GetSMNumber(cl2);
186
187   if(nm1==-1 || nm2==-1) assert(0);
188
189   if(nm1 != nm2) return; // Both cluster should be in the same SM
190
191   AliESDCaloCluster* cl = cl1;
192   if(cl1->E() < cl2->E()) cl = cl2; // Get cluster with highest energy
193
194   const Int_t  nDigits  = cl->GetNumberOfDigits();
195   const Short_t* absId   = cl->GetDigitIndex()->GetArray();
196
197   AliEMCALCalibCoefs *tFirst = GetCCFirst();
198   calibCoef *rFirst=0;
199
200   int    indMax = 0, id=0;
201   id = Int_t(absId[0]);
202   rFirst = tFirst->GetTable(id);
203   double emax   = cl->GetTrueDigitEnergy(indMax, rFirst->cc);
204   if(nDigits > 1) {
205     for(int i=1; i<nDigits; i++) {
206       id = Int_t(absId[i]);
207       rFirst = tFirst->GetTable(id);
208       if(emax < cl->GetTrueDigitEnergy(i, rFirst->cc)) {
209         indMax = i;
210         emax   = cl->GetTrueDigitEnergy(i, rFirst->cc);
211       }
212     }
213   }
214   if(emax/cl->E() > 0.5) { // more than 50% of cluster energy
215     u::FillH1(fLhists, 0, 2.); // number "good" entries 
216     absIdMax = Int_t(absId[indMax]);
217     FillPi0Candidate(mgg, absIdMax, nm1);
218   }
219 }
220
221 void AliEMCALFolder::FillPi0Candidate(const Double_t mgg, Int_t absIdMax, Int_t nm)
222 {
223   //
224   // Jun 08;
225   //
226   static Int_t nSupModMax, nModuleMax, nIphiMax, nIetaMax, iphiCellMax, ietaCellMax;
227   static AliEMCALCell* cell;
228   static TFolder *set; 
229   static AliEMCALSuperModule *sm;
230
231   fGeometry->GetCellIndex(absIdMax, nSupModMax, nModuleMax, nIphiMax, nIetaMax);
232   if(nm != nSupModMax) assert(0);
233
234   fGeometry->GetCellPhiEtaIndexInSModule(nSupModMax, nModuleMax, nIphiMax, nIetaMax, iphiCellMax, ietaCellMax);
235
236   cell = 0;
237   set  = 0;
238   sm   = 0;
239   if(GetCell(absIdMax)==0) {
240     cell = new AliEMCALCell(absIdMax, Form("sm%2.2i:phi%2.2i:eta%2.2i(%4.4i)",nSupModMax,iphiCellMax,ietaCellMax,absIdMax));
241     SetCell(cell, absIdMax);
242     // For browser
243     set  = dynamic_cast<TFolder*>(FindObject(Form("SM%2.2i",nSupModMax)));
244     if(set==0) {
245       sm = new  AliEMCALSuperModule(nSupModMax);
246       Add(sm);
247       sm->SetParent(this);
248       sm->Init();
249     } else {
250       sm = dynamic_cast<AliEMCALSuperModule*>(set);
251     }
252     if(sm) {
253       sm->AddCellToEtaRow(cell, ietaCellMax);
254       cell->SetParent(sm);
255     }
256     // 
257     cell->SetCCfromCCTable(GetCCIn());
258   } else {
259     cell = GetCell(absIdMax);
260     set  = dynamic_cast<TFolder*>(FindObject(Form("SM%2.2i",nm)));
261     if(set) sm = (AliEMCALSuperModule*)set;
262   }
263   if(sm == 0) assert(0);
264   if(nm != sm->GetSMNumber()) assert(0);
265
266   u::FillH1(sm->GetHists(), 0, mgg);
267   cell->FillEffMass(mgg);
268 }
269
270 void AliEMCALFolder::FitAllSMs()
271 { // Jun 14, 2007
272   // Only first SM now - should be changed in the future 
273    AliEMCALSuperModule *sm0 = GetSuperModule(0);
274    sm0->FitForAllCells();
275    // Get input calibration table
276    AliEMCALCalibCoefs *ccIn = GetCCIn();
277    if(ccIn == 0) {
278      printf("<E> no input cc \n"); 
279      return;
280    }
281    // New calibration table
282    AliEMCALCalibCoefs *ccOut = new AliEMCALCalibCoefs(fgkCCoutName.Data(), ccIn->GetNRows());
283    calibCoef *rIn=0, rOut;
284    for(Int_t i=0; i<ccIn->GetNRows(); i++){
285      rIn  = ccIn->GetTable(i);
286      rOut = (*rIn);
287      AliEMCALCell* cell = GetCell(rIn->absId);
288      if(cell && cell->GetSupMod() == 0) { // only first module now
289        rOut.cc = cell->GetCcOut();
290      }
291      ccOut->AddAt(&rOut);
292    }
293    ccOut->Purge();
294    Add(ccOut);
295 }
296
297 AliEMCALCalibCoefs* AliEMCALFolder::GetCCTable(const char* name)
298 {
299   TObject *obj = FindObject(name);
300   if(obj) return (AliEMCALCalibCoefs*)obj;
301   else    return 0;
302 }
303
304 Int_t AliEMCALFolder::GetSMNumber(AliESDCaloCluster* cl)
305 {
306   static Int_t absId, nSupMod, nModule, nIphi, nIeta;
307   nSupMod = -1; // if negative something wrong
308   if(cl) {
309     absId = Int_t(cl->GetDigitIndex()->At(0));
310     fGeometry->GetCellIndex(absId, nSupMod, nModule, nIphi, nIeta);
311   }
312   return nSupMod;
313 }
314
315 // Recalibration staf - Jun 18,2007
316 AliEMCALRecPoint* AliEMCALFolder::GetRecPoint(AliESDCaloCluster *cl, AliEMCALCalibCoefs *tOld,AliEMCALCalibCoefs *tNew,
317 TList *l, Double_t deff, Double_t w0, Double_t phiSlope)
318 {
319   //
320   // Static function;
321   // Get recalibrated rec.point from ESD cluster
322   // If tNew == 0 -> get ideal calibration with  ADCCHANNELEC
323   //
324   static Double_t ADCCHANNELEC = 0.0153;  // Update 24 Apr 2007: 250./16/1024 - width of one ADC channel in GeV
325   //static Float_t ECAW0 = 4.5; // hard coded now - see AliEMCALClusterizerv1::InitParameters()
326   static Double_t ECAW0 = 5.5; // Beter case for simulation 
327   Int_t ampDigi=0, indMax=-1;
328   Double_t eDigiNew=0.0, eDigiMax=0.0; 
329   static TArrayD ed;
330
331   if(w0 > 0.5) ECAW0 = w0;
332   AliEMCALRecPoint *rp=0;
333   calibCoef *rOld=0, *rNew=0;
334   //   printf(" AliEMCALFolder::GetRecPoint() : RECALIBRATION : w0 %f ECAW0 %f \n", w0, ECAW0);
335   if(cl && tOld){
336     // cl->PrintClusterInfo(1);
337     const Int_t ndg   = cl->GetNumberOfDigits();
338     //    const Int_t nprim = cl->GetNumberOfPrimaries();
339
340     const Short_t  prim    = cl->GetLabel();
341     const Short_t* dgAbsId = cl->GetDigitIndex()->GetArray();
342     //    const UShort_t* dgAmp   = cl->GetDigitAmplitude(); // This is energy - bad definition
343
344     rp = new AliEMCALRecPoint(""); // opt=""
345     rp->SetClusterType(AliESDCaloCluster::kClusterv1);
346     AliEMCALDigit* dg=0;
347     TClonesArray digits("AliEMCALDigit", ndg);
348     Int_t absId = 0;
349     ed.Set(ndg); // resize array
350     for(Int_t i=0; i<ndg; i++){
351       // Save just abs id and amplitude of digits which will be used for recalculation of
352       // cluster energy and position
353       absId   = Int_t(dgAbsId[i]);
354       rOld    = tOld->GetTable(absId);
355       ampDigi = cl->GetTrueDigitAmplitude(i, rOld->cc); // True amplitude
356
357       new(digits[i]) AliEMCALDigit(Int_t(prim),0, absId, ampDigi, 0.0, i, 0.0); 
358       dg     =  (AliEMCALDigit*)digits[i];
359
360       if(tNew) rNew = tNew->GetTable(absId);
361       if(rNew) {
362         rNew  = tNew->GetTable(absId);
363         eDigiNew = Double_t(ampDigi) * rNew->cc;     // Recalibrate with new cc
364       } else {
365         eDigiNew = Double_t(ampDigi) * ADCCHANNELEC; // Ideal calibration
366       }
367       //eDigiNew = Double_t(cl->GetTrueDigitEnergy(i)); // Copy from ESD for checking
368       rp->AddDigit(*dg, Float_t(eDigiNew));
369       ed[i] = eDigiNew;
370       if(eDigiMax<eDigiNew) {
371         eDigiMax = eDigiNew;
372         indMax   = i;
373       }
374       u::FillH1(l, 14, eDigiNew);
375       u::FillH1(l, 15, Double_t(absId));
376       //printf("<I> digit %i amp %i rOld->cc %6.5f GeV rNew->cc %6.5f GeV\n", i, ampDigi, rOld->cc, rNew->cc);
377     }
378     //printf("<I> recalibration of digits was done ! \n");
379     //    rp->EvalAll(ECAW0, &digits);
380     if(indMax>=0) rp->SetIndMaxDigit(indMax);
381     if(deff > 0.0) { // for fit
382       rp->EvalLocalPositionFit(deff, ECAW0, phiSlope, &digits); // I need just position
383     } else { // get w0 and deff from parametrisation - Sep 4, 2007 
384       rp->EvalLocalPosition2(&digits, ed);
385     }
386     digits.Delete();
387   }
388   //rp->Print("print");
389   return rp;
390 }
391
392 void  AliEMCALFolder::Save(const char *fn, const char *opt)
393
394   //
395   // Jun 5, 2007; See TFileIter and StFMC.cxx
396   // Jul 16 - added fgkDirOfRootFiles
397   // Sep 7, 2007 - should be changed without TFileIter
398   /*
399   TString FN = fgkDirOfRootFiles;
400   FN += fn;
401   if(FN.Contains(".root")==0) FN += ".root";
402   TFileIter f(FN.Data(),opt,"EMCAL object");
403   UInt_t eventNum  = 0; // just one object
404   UInt_t runNumber = 0; // 0 now, - may statistics on selector
405   f.NextEventPut(this, eventNum, runNumber);
406   printf(" Save %s to file %s\n", GetName(), FN.Data());
407   */
408 }
409
410 AliEMCALFolder* AliEMCALFolder::Read(const char *fn, const char *opt)
411
412   //
413   // Jun 27, 2007
414   // Jul 16 - added fgkDirOfRootFiles
415   //
416   printf("<I> AliEMCALFolder::Read(%s,%s)\n",fn,opt); 
417   AliEMCALFolder* EMCAL = 0;
418   TH1::AddDirectory(0); // this is obligatory
419
420   TString FN = fgkDirOfRootFiles;
421   FN += fn;
422   if(FN.Contains(".root")==0) FN += ".root";
423
424   TFile f(FN.Data(),opt);
425   if(f.IsOpen()) {
426     TList *l = f.GetListOfKeys();
427     printf("<I> The total number of the objects: %d \n File %s\n", l->GetSize(), FN.Data());
428
429     TKey *key = (TKey*)l->At(0);
430     EMCAL = (AliEMCALFolder*)key->ReadObj();
431     f.Close();
432     if(EMCAL) EMCAL->InitAfterRead();
433   }
434   return EMCAL;
435 }
436
437
438 void  AliEMCALFolder::InitAfterRead()
439
440   lObj    = (TList*)FindObject("Objects");
441   if(lObj) {
442     fLhists = (TList*)lObj->FindObject("HistsOfEmcal");
443   }
444 }
445
446 void AliEMCALFolder::DrawQA(const int nsm)
447
448   //
449   // Jun 25, 2007
450   //
451
452   AliEMCALSuperModule* sm = GetSuperModule(nsm);
453   if(sm==0) return;
454   TList *l = sm-> GetHists();
455   Int_t nx=2, ny=2, wh=530, ww=750;
456
457   TCanvas *c = new TCanvas(Form("QA_%i",GetIterationNumber()), Form("QA_%i",GetIterationNumber()),
458                            10, 10, ww, wh);
459   c->Divide(nx,ny);
460
461   int ic=1;
462   c->cd(ic++);
463   TH1 *h1 = (TH1*)l->At(0);
464   sm->FitEffMassHist(); 
465   u::DrawHist(h1,2);
466   h1->SetAxisRange(0.03, 0.28);
467
468   c->cd(ic++);
469   sm->DrawCC(0);  
470   TH1 *hccin = (TH1*)l->At(1);
471   hccin->SetAxisRange(14., 20.);
472
473   gStyle->SetOptStat(1111);
474   c->cd(ic++);
475   TH1 *hmass = (TH1*)l->At(3);
476   u::DrawHist(hmass, 2);
477   hmass->SetAxisRange(0.12, 0.16);
478
479   c->cd(ic++);
480   TH1 *hres = (TH1*)l->At(4);
481   u::DrawHist(hres, 2);
482   hres->SetAxisRange(0.05, 0.120);
483
484   if(ic<nx*ny) {
485     c->cd(ic++);
486     u::DrawHist((TH1*)l->At(5), 2);
487
488     c->cd(ic++);
489     u::DrawHist((TH1*)l->At(6), 2);
490   }
491   c->Update();
492 }
493
494 TList* AliEMCALFolder::BookHists()
495 {
496   gROOT->cd();
497   TH1::AddDirectory(1);  
498
499   new TH1F("00_HStat", "hist of common EMCAL statistics", 100, 0.5, 100.5);
500   new TH1F("01_EffMassAll", "effective mass of #gamma,#gamma(m_{#pi^{0}}=134.9766 MeV) - whole EMCAL", 250,0.0,0.5);
501
502   TList *l = AliEMCALHistoUtilities::MoveHistsToList("HistsOfEmcal", kFALSE);
503
504   TH1::AddDirectory(0);
505   return l;
506 }
507
508 void AliEMCALFolder::CreateCellNtuple()
509 {
510   // Jun 28, 2007
511   if(fCellNtuple) { // Already exist
512     fCellNtuple->Print();
513     return;
514   }
515   // Create ntuple
516   Int_t bsize = int(1.e+5);
517   fCellNtuple = new TNtuple("cells","Cells Ntuple for quick analysis",
518                             "fAbsId:fSupMod:fModule:fPhi:fEta:fPhiCell:fEtaCell:fCcIn:fCcOut", bsize);
519   AliEMCALCell *cell=0;
520   // Fill ntuple
521   AliEMCALSuperModule* sm = GetSuperModule(0);
522   if(sm) printf(" TNtuple was created ! sm0 %s \n", sm->GetName());
523   for(int eta=0; eta<48; eta++) { // eta row
524     TFolder *setEta = dynamic_cast<TFolder*>(sm->FindObject(Form("ETA%2.2i",eta)));
525     if(setEta) {
526       printf(" ***** eta row %s ******\n", setEta->GetName());
527       TList* l = (TList*)setEta->GetListOfFolders();
528       for(int phi=0; phi<l->GetSize(); phi++) { // cycle on cells (phi directions)
529         cell = (AliEMCALCell*)l->At(phi);
530         if(cell) {
531           cell->FillCellNtuple(fCellNtuple);
532           //printf(" fill cell %s : %s \n", cell->GetName(), cell->GetTitle());
533         }
534       }
535     }
536   }
537   fCellNtuple->Print();
538   lObj->Add(fCellNtuple); 
539 }
540
541 void AliEMCALFolder::CreateAndFillAdditionalHists()
542 {
543   CreateCellNtuple();
544   TH1::AddDirectory(0);  
545   fLhists->Add(new TH1F("02_CCoutOnEdge", "cc out on edge of calorimeter (in MeV)", 70, 12., 19.));
546   fLhists->Add(new TH1F("03_CCoutInside", "cc out inside  of calorimeter (in MeV)", 70, 12., 19.));
547   fLhists->Add(new TH1F("04_CCoutOnEdge2", "cc out on edge of calorimeter(2) (in MeV)", 70, 12., 19.));
548   // Fill
549   Float_t* args;
550   for(Int_t i=0; i<(Int_t)fCellNtuple->GetEntries(); i++){
551     fCellNtuple->GetEvent(i);
552     args = fCellNtuple->GetArgs();
553     Int_t phi = (Int_t)args[5];
554     Int_t eta = (Int_t)args[6];
555     Double_t cc = (Double_t)args[8]*1000.;
556     if     ((phi==0||phi==23) || (eta==0||eta==47)) u::FillH1(fLhists, 2, cc);
557     else if((phi==1||phi==22) || (eta==1||eta==46)) u::FillH1(fLhists, 4, cc); // next to edge
558     else                                            u::FillH1(fLhists, 3, cc);
559   }
560   // Drawing
561   Int_t wh=530, ww=750;
562   TCanvas *c = new TCanvas("c_edge","CEDGE", 10, 10, ww, wh);
563
564   gStyle->SetOptStat(1100);
565   gStyle->SetOptFit(111);
566   TH1 *h1 = (TH1*)fLhists->At(3);
567   TF1 *g = u::Gausi("ccInside", 14.7, 16.4, h1);
568   g->SetLineColor(kRed);
569   h1->Fit(g,"Q+","", 14.7, 16.4);
570   u::DrawHist(h1,2);
571   h1->SetTitle("CC distribution after #pi^{0} calibration");
572   h1->SetXTitle("  MeV  ");
573   h1->SetYTitle("  N  ");
574   TLatex *lat1 = u::Lat(Form("rel.width = %4.2f%%", 
575   100.*h1->GetRMS()/ h1->GetMean()), 16.5, 100., 12, 0.045);
576   TLatex *lat2 = u::Lat(Form("rel.width = %4.2f%% (from fit)", 
577                              100.*g->GetParameter(2)/ g->GetParameter(1)), 16.5, 70., 12, 0.045);
578
579   if(0) {
580     TH1 *h2 = (TH1*)fLhists->At(2);
581     u::DrawHist(h2,2,1,"same",2);
582   }
583
584   TH1F *hccFirst = AliEMCALCalibCoefs::GetHistOfCalibTableFromDb("ccTmp");
585   u::DrawHist(hccFirst,2,1,"same",3);
586
587   
588   // Ideal calibration - Jul 18, 2007
589   Double_t ADCCHANNELEC = 0.0153, ccIdeal =  ADCCHANNELEC*1.e+3;
590   Double_t ym = h1->GetMaximum();
591   TLine *l = new TLine(ccIdeal,-ym*0.05, ccIdeal,ym*1.05);
592   l->SetLineColor(kBlue);
593   l->Draw(); 
594
595   TLegend *leg = new TLegend(0.1,0.6, 0.45,0.85);
596   leg->AddEntry(hccFirst, "Initial cc ", "L");
597   leg->AddEntry(h1, "Final cc", "L");
598   TLegendEntry *le=0;
599   if(l) {
600     le  = leg->AddEntry(l, "Ideal calibration", "L");
601     le->SetTextColor(l->GetLineColor());
602   }
603
604   TH1 *hCcEdge  = (TH1*)fLhists->At(2);
605   TH1 *hCcEdge2 = (TH1*)fLhists->At(4);
606   if(1) {
607     u::DrawHist(hCcEdge,2,kGreen,"same",1);
608     le = leg->AddEntry(hCcEdge , "Edge cell", "L");
609     le->SetTextColor(hCcEdge->GetLineColor());
610
611     u::DrawHist(hCcEdge2,2, 28,"same",1); // 28 - like brown
612     le = leg->AddEntry(hCcEdge2 , "Edge cell (2)", "L");
613     le->SetTextColor(hCcEdge2->GetLineColor());
614     u::DrawHist(h1,2,1,"same");
615   }
616
617   leg->Draw();
618
619   c->Update();
620 }
621
622 void AliEMCALFolder::TestSMStruct()
623 {
624   // testing May 22, 2007
625   for(int m=0; m<12; m++) {
626     AliEMCALSuperModule *sm = new AliEMCALSuperModule(m);
627     Add(sm);
628     sm->SetParent(this);
629   }
630 }
631