]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSDA2.cxx
Minor change for LHC11e
[u/mrichter/AliRoot.git] / PHOS / AliPHOSDA2.cxx
1 #include "AliPHOSDA2.h"
2 #include "TString.h"
3
4 ClassImp(AliPHOSDA2)
5
6 //----------------------------------------------------------------
7 AliPHOSDA2::AliPHOSDA2(int module) : TNamed(),
8  fHistoFile(0),fFiredCells(0),fMod(module)
9
10 {
11   // Create AliPHOSDA2 ("Bad channels finder") object.
12   // module is the PHOS module number (0..4).
13   // Quality histogram names: module_iX_iZ_gain.
14   // Root file name: PHOS_ModuleX_BCM.root, where X - module number.
15   
16   char name[128];
17   TString sname="PHOS_Module%d_BCM";
18   snprintf(name,sname.Length(),sname.Data(),fMod);
19   SetName(name);
20
21   SetTitle("Detector Algorithm to check for PHOS channels quality");
22
23   char rootname[128];
24   TString srootname="%s.root";
25   snprintf(rootname,srootname.Length(),srootname.Data(),GetName());
26
27   fHistoFile =  new TFile(rootname,"recreate"); // new file!
28   
29   for(Int_t iX=0; iX<64; iX++) {
30     for(Int_t iZ=0; iZ<56; iZ++) {
31       for(Int_t iGain=0; iGain<2; iGain++) {
32         fHQuality[iX][iZ][iGain] = 0;
33       }
34     }
35   }
36
37   fMaps[0]=0;
38   fMaps[1]=0;
39
40   fFiredCells = new TH1I("fFiredCells","Number of fired cells per event",100,0,1000);
41
42 }
43
44 //----------------------------------------------------------------
45 AliPHOSDA2::AliPHOSDA2(Int_t module, TObjArray* oldHistos) : TNamed(),
46  fHistoFile(0),fFiredCells(0),fMod(module)
47
48 {
49   // Create AliPHOSDA2 ("Bad channels finder") object.
50   // module is the PHOS module number (0..4).
51   // Quality histogram names: module_iX_iZ_gain.
52   // Read histograms from array oldHistos (if any).
53   // Do not produce an output file!
54
55   char name[128];
56   TString sname="PHOS_Module%d_BCM";
57   snprintf(name,sname.Length(),sname.Data(),fMod);
58   SetName(name);
59
60   SetTitle("Detector Algorithm to check for PHOS channels quality");
61
62   char hname[128];
63   TH1F* hist1=0;
64   TString shname = "%d_%d_%d_%d";
65
66   for(Int_t iX=0; iX<64; iX++) {
67     for(Int_t iZ=0; iZ<56; iZ++) {
68       for(Int_t iGain=0; iGain<2; iGain++) {
69         snprintf(hname,shname.Length(),shname.Data(),fMod,iX,iZ,iGain);
70         if(oldHistos) 
71           hist1 = (TH1F*)oldHistos->FindObject(hname);
72         if(hist1) fHQuality[iX][iZ][iGain] = hist1;
73         else
74           fHQuality[iX][iZ][iGain] = 0;
75       }
76     }
77   }
78   
79   fMaps[0]=0;
80   fMaps[1]=0;
81   
82   fFiredCells = new TH1I("fFiredCells","Number of fired cells per event",100,0,1000);
83   
84 }
85
86 //-------------------------------------------------------------------
87 AliPHOSDA2::AliPHOSDA2(const AliPHOSDA2& da) : TNamed(da),
88   fHistoFile(0),fFiredCells(0),fMod(da.fMod)
89 {
90   // Copy constructor.
91
92   char hname[128];
93   TH1F* hist1=0;
94   TString shname = "%d_%d_%d_%d";
95
96   for(Int_t iX=0; iX<64; iX++) {
97     for(Int_t iZ=0; iZ<56; iZ++) {
98       for(Int_t iGain=0; iGain<2; iGain++) {
99
100         snprintf(hname,shname.Length(),shname.Data(),fMod,iX,iZ,iGain);
101         hist1 = (TH1F*)da.fHistoFile->Get(hname);
102         if(hist1) fHQuality[iX][iZ][iGain] = new TH1F(*hist1);
103         else
104           fHQuality[iX][iZ][iGain] = 0;
105       }
106     }
107   }
108   
109   if(da.fMaps[0]) 
110     fMaps[0] = new TH2F(*da.fMaps[0]);
111   else
112     fMaps[0] = 0;
113
114   if(da.fMaps[1]) 
115     fMaps[1] = new TH2F(*da.fMaps[1]);
116   else
117     fMaps[1] = 0;
118   
119   fHistoFile = new TFile(da.GetName(),"recreate");
120   fFiredCells = new TH1I(*da.fFiredCells);
121   
122 }
123
124 //-------------------------------------------------------------------
125 AliPHOSDA2& AliPHOSDA2::operator= (const AliPHOSDA2& da)
126 {
127   //Assignment operator.
128
129   if(this != &da) {
130
131     TString oldname(fHistoFile->GetName());
132     TString newname(da.fHistoFile->GetName());
133
134     if(oldname != newname) {
135       delete fHistoFile;
136       fHistoFile = new TFile(da.fHistoFile->GetName(),"update");
137     }
138
139     fMod = da.fMod;
140
141     SetName(da.GetName());
142     SetTitle(da.GetTitle());
143
144     for(Int_t iX=0; iX<64; iX++) {
145       for(Int_t iZ=0; iZ<56; iZ++) {
146         for(Int_t iGain=0; iGain<2; iGain++) {
147           if (fHQuality[iX][iZ][iGain]) delete fHQuality[iX][iZ][iGain];
148           fHQuality[iX][iZ][iGain] = da.fHQuality[iX][iZ][iGain];
149         }
150       }
151     }
152
153     if(fMaps[0]) { 
154       delete fMaps[0];
155       fMaps[0] = da.fMaps[0];
156     } 
157
158     if(fMaps[1]) { 
159       delete fMaps[1];
160       fMaps[1] = da.fMaps[1];
161     } 
162     
163     if(fFiredCells) {
164       delete fFiredCells;
165       fFiredCells = da.fFiredCells;
166     }
167     
168   }
169   
170   return *this;
171 }
172
173
174 //-------------------------------------------------------------------
175 AliPHOSDA2::~AliPHOSDA2()
176 {
177   // Destructor
178   
179   UpdateHistoFile();
180   if(fHistoFile) delete fHistoFile;
181   
182 }
183
184 //-------------------------------------------------------------------
185 void AliPHOSDA2::FillQualityHistograms(Float_t quality[64][56][2]) 
186 {
187   // Fills quality histograms.
188   // _By definition_, qood quality is 0<quality<1, 
189   // all outside that is a bad quality.
190   // If no quality value read for particular channel, 
191   // the correspondent array entry should be filled by zero.
192   // WARNING: this function should be called once per event!
193
194   char hname[128];
195   char htitl[128];
196
197   for(Int_t iX=0; iX<64; iX++) {
198     for (Int_t iZ=0; iZ<56; iZ++) {
199
200       for(Int_t iGain=0; iGain<2; iGain++) {
201         if(!quality[iX][iZ][iGain]) continue;
202         
203         if(fHQuality[iX][iZ][iGain]) 
204           fHQuality[iX][iZ][iGain]->Fill(quality[iX][iZ][iGain]);
205         else {
206           snprintf(hname,128,"%d_%d_%d_%d",fMod,iX,iZ,iGain);
207           snprintf(htitl,128,"Quality for crystal %d_%d_%d and gain %d",fMod,iX,iZ,iGain);
208           fHQuality[iX][iZ][iGain] = new TH1F(hname,htitl,100,1.e-6,10.);
209           fHQuality[iX][iZ][iGain]->Fill(quality[iX][iZ][iGain]);
210         }
211       }
212
213     }
214   }
215
216 }
217
218 //-------------------------------------------------------------------
219 void  AliPHOSDA2::FillFiredCellsHistogram(Int_t nCells)
220 {
221   fFiredCells->Fill(nCells);
222 }
223
224 //-------------------------------------------------------------------
225 void AliPHOSDA2::UpdateHistoFile()
226 {
227   // Write histograms to file
228
229   if(!fHistoFile) return;
230   if(!fHistoFile->IsOpen()) return;
231
232   char titl[128];
233
234   if(fMaps[0]) 
235     fMaps[0]->Reset();
236   else {
237     snprintf(titl,128,"Quality map for Low gain");
238     fMaps[0] = new TH2F("gmaplow",  titl, 64,0.,64.,56,0.,56.);
239   }
240
241   if(fMaps[1]) 
242     fMaps[1]->Reset();
243   else {
244     snprintf(titl,128,"Quality map for High gain");
245     fMaps[1] = new TH2F("gmaphigh", titl, 64,0.,64.,56,0.,56.);
246   }
247     
248   TH1F* hist1=0;
249
250   for(Int_t iX=0; iX<64; iX++) {
251     for(Int_t iZ=0; iZ<56; iZ++) {
252
253       for(Int_t iGain=0; iGain<2; iGain++) {
254         hist1 = fHQuality[iX][iZ][iGain];
255         if(hist1) { 
256           hist1->Write(hist1->GetName(),TObject::kWriteDelete);
257           Double_t mean = hist1->GetMean();
258           fMaps[iGain]->SetBinContent(iX+1,iZ+1,mean);
259         }
260       } 
261
262     }
263   }
264
265   fMaps[0]->Write(fMaps[0]->GetName(),TObject::kWriteDelete);
266   fMaps[1]->Write(fMaps[1]->GetName(),TObject::kWriteDelete);
267
268   fFiredCells->Write();
269
270 }
271