]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALSuperModule.cxx
Initialization of some data members (Christian)
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALSuperModule.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 /* 
17 $Log$ 
18 */
19
20 //_________________________________________________________________________
21 // Top EMCAL folder which will keep all information about EMCAL itself,
22 // super Modules (SM), modules, towers, set of hists and so on.
23 //  TObjectSet -> TFolder; Sep 6, 2007
24 //
25 //*-- Author: Aleksei Pavlinov (WSU, Detroit, USA) 
26
27 #include "AliEMCALSuperModule.h"
28 #include "AliEMCALFolder.h"
29 #include "AliEMCALCell.h"
30 #include "AliEMCALHistoUtilities.h"
31
32 #include <TROOT.h>
33 #include <TStyle.h>
34 #include <TList.h>
35 #include <TH1.h>
36 #include <TF1.h>
37 #include <TCanvas.h>
38 #include <TLegend.h>
39 #include <TLegendEntry.h>
40
41 typedef  AliEMCALHistoUtilities u;
42
43 ClassImp(AliEMCALSuperModule)
44
45 AliEMCALSuperModule::AliEMCALSuperModule() : TFolder()
46 , fParent(0),fLh(0),fSMNumber(0)
47 {
48 }
49
50 AliEMCALSuperModule::AliEMCALSuperModule(const Int_t m, const char* title) : 
51 TFolder(Form("SM%2.2i",m), title)
52 , fParent(0),fLh(0),fSMNumber(m)
53
54
55
56 AliEMCALSuperModule::~AliEMCALSuperModule()
57 {
58   // dtor
59 }
60
61 void AliEMCALSuperModule::Init()
62 {
63   if(GetHists()==0) {
64     fLh = BookHists();
65     Add(fLh);
66   }
67 }
68
69 void  AliEMCALSuperModule::AddCellToEtaRow(AliEMCALCell *cell, const Int_t etaRow)
70 {
71   static TFolder *set;
72   set = dynamic_cast<TFolder*>(FindObject(Form("ETA%2.2i",etaRow))); 
73   if(set==0) {
74     set = new  TFolder(Form("ETA%2.2i",etaRow),"eta row");
75     Add(set);
76   }
77   set->Add(cell);
78 }
79
80 void AliEMCALSuperModule::FitForAllCells()
81 {
82   Int_t ncells=0;
83   for(int eta=0; eta<48; eta++) { // eta row
84     TFolder *setEta = dynamic_cast<TFolder*>(FindObject(Form("ETA%2.2i",eta)));
85     if(setEta) {
86       TList* l = (TList*)setEta->GetListOfFolders();
87       printf(" eta %i : %s : cells %i \n", eta, setEta->GetName(), l->GetSize());
88       for(int phi=0; phi<l->GetSize(); phi++) { // cycle on cells (phi directions)
89         AliEMCALCell* cell = dynamic_cast<AliEMCALCell*>(l->At(phi));
90         if(cell == 0) continue; 
91
92         cell->FitEffMassHist();
93
94         u::FillH1(GetHists(), 1, cell->GetCcIn()*1.e+3);
95         u::FillH1(GetHists(), 2, cell->GetCcOut()*1.e+3);
96
97         TF1 *f = cell->GetFunction();
98         u::FillH1(GetHists(), 3, f->GetParameter(1));
99         u::FillH1(GetHists(), 4, f->GetParameter(2));
100         u::FillH1(GetHists(), 5, f->GetChisquare()/f->GetNDF());
101         u::FillH1(GetHists(), 6, f->GetParameter(0));
102
103         ncells++;
104       }
105     }
106   }
107   printf(" <I> AliEMCALSuperModule::FitForAllCells() : ncells %i with fit \n", ncells);
108 }
109
110 void AliEMCALSuperModule::FitEffMassHist()
111 {
112   TH1* h = (TH1*)GetHists()->At(0);
113   AliEMCALCell::FitHist(h, GetName());
114 }
115
116
117 void AliEMCALSuperModule::PrintInfo()
118 {
119   printf(" Super Module :   %s    :   %i \n", GetName(), fSMNumber);
120   printf(" # of active cells                %i \n", GetNumberOfCells());
121   TH1* h = (TH1*)GetHists()->At(0);
122   printf(" # h->Integral() of eff.mass hist %i \n", int(h->Integral()));
123 }
124
125 void AliEMCALSuperModule::DrawCC(int iopt)
126 {
127   TCanvas *c=0; 
128   if(iopt==1) c = new TCanvas("ccInOut","ccInOut");
129  
130   gStyle->SetOptStat(0);
131
132   TH1 *hCCIn  = (TH1*)GetHists()->At(1);
133   TH1 *hCCOut = (TH1*)GetHists()->At(2);
134
135   if(hCCIn == 0)              return;
136   if(hCCIn->GetEntries()<10.) return;
137
138   hCCIn->SetStats(kFALSE);
139   hCCOut->SetStats(kFALSE);
140   hCCOut->SetTitle("CC in and out; Iter 1; Jul 26; All Statistics");
141   hCCOut->SetXTitle("cc in MeV");
142   hCCOut->SetYTitle("  N  ");
143
144   u::DrawHist(hCCOut,2);
145   hCCOut->SetAxisRange(10., 24.);
146   u::DrawHist(hCCIn,2, kRed, "same");
147
148   TLegend *leg = new TLegend(0.5,0.36, 0.97,0.80);
149   TLegendEntry *leIn = leg->AddEntry(hCCIn, Form("input cc : %6.2f #pm %6.2f", hCCIn->GetMean(),hCCIn->GetRMS()), "L");
150   leIn->SetTextColor(hCCIn->GetLineColor());
151
152   if(hCCOut->GetEntries()>10.) 
153   leg->AddEntry(hCCOut, Form("output cc : %6.2f #pm %6.2f", hCCOut->GetMean(),hCCOut->GetRMS()), "L");
154
155   leg->Draw();
156
157   if(c) c->Update();
158 }
159
160 Int_t AliEMCALSuperModule::GetNumberOfCells()
161 {
162   Int_t ncells=0;
163   TList* l = (TList*)GetListOfFolders();
164   for(int eta=0; eta<l->GetSize(); eta++) { // cycle on eta row
165     TFolder *setEta = dynamic_cast<TFolder*>(l->At(eta));
166     if(setEta==0) continue;
167
168     TList* le = (TList*)setEta->GetListOfFolders();
169     ncells += le->GetSize();
170   }
171   return ncells;
172 }
173
174 TList* AliEMCALSuperModule::BookHists()
175 {
176   gROOT->cd();
177   TH1::AddDirectory(1);
178
179   AliEMCALFolder* EMCAL = (AliEMCALFolder*)GetParent(); 
180   Int_t it = EMCAL->GetIterationNumber();
181
182   new TH1F("00_EffMass",  "effective mass of #gamma,#gamma(m_{#pi^{0}}=134.98 MeV) ", 250,0.0,0.5);
183   new TH1F("01_CCInput",  "input CC dist.(MEV) ", 200, 5., 25.);
184   new TH1F("02_CCOutput", "output CC dist.(MEV) ", 200, 5., 25.);
185   new TH1F("03_MPI0", "mass of #pi_{0} dist. ", 170, 0.05, 0.22);
186   new TH1F("04_RESPI0", "resolution at #pi_{0} dist. ", 50, 0.0, 0.05);
187   new TH1F("05_XI2/NDF", "#chi^{2} / ndf", 50, 0.0, 5.0);
188   new TH1F("06_NPI0", "number of #pi_{0}", 150, 0.0, 1500.);
189
190   TList *l = AliEMCALHistoUtilities::MoveHistsToList(Form("HistsOfSM_%2.2i",fSMNumber), kFALSE);
191   AliEMCALHistoUtilities::AddToNameAndTitleToList(l, Form("_%2.2i_It%i",fSMNumber, it), 
192                                                   Form(" SM %2.2i, Iter %i",fSMNumber, it));
193
194   TH1::AddDirectory(0);
195   return l;
196 }