]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSDA1.cxx
Constructor added
[u/mrichter/AliRoot.git] / PHOS / AliPHOSDA1.cxx
1 #include "AliPHOSDA1.h"
2 #include "TString.h"
3
4 ClassImp(AliPHOSDA1)
5
6 //----------------------------------------------------------------
7 AliPHOSDA1::AliPHOSDA1(int module) : TNamed(),
8   fHistoFile(0),fMod(module),fWriteToFile(kTRUE)
9
10 {
11   // Create DA1 ("Calibration DA") object.
12   // module is the PHOS module number (0..4).
13   // Checks existence of histograms which might have been left
14   // from the previous runs to continue their filling.
15   // Histogram names: module_iX_iZ_gain for TH2 and  module_iX_iZ for TH1.
16   // Root file name: PHOS_ModuleX_Calib.root, where X - module number.
17   
18   char name[128];
19   sprintf(name,"PHOS_Module%d_Calib",fMod);
20   SetName(name);
21
22   SetTitle("Calibration Detector Algorithm");
23
24   char rootname[128];
25   sprintf(rootname,"%s.root",GetName());
26
27   fHistoFile =  new TFile(rootname,"update");
28
29   char hname[128];
30   TH1F* hist1=0;
31   TH2F* hist2=0;
32
33   for(Int_t iX=0; iX<64; iX++) {
34     for(Int_t iZ=0; iZ<56; iZ++) {
35
36       sprintf(hname,"%d_%d_%d",fMod,iX,iZ);
37       hist1 = (TH1F*)fHistoFile->Get(hname);
38       if(hist1) fHgLgRatio[iX][iZ] = hist1;
39       else
40         fHgLgRatio[iX][iZ] = 0;
41
42       for(Int_t iGain=0; iGain<2; iGain++) {
43         sprintf(hname,"%d_%d_%d_%d",fMod,iX,iZ,iGain);
44         hist2 = (TH2F*)fHistoFile->Get(hname);
45         if(hist2) fTimeEnergy[iX][iZ][iGain] = hist2;
46         else
47           fTimeEnergy[iX][iZ][iGain] = 0;
48       }
49
50     }
51   }
52  
53 }
54
55 //-------------------------------------------------------------------
56 AliPHOSDA1::AliPHOSDA1(Int_t module, TH2F oldTimeEnergy[64][56][2]) :
57   TNamed(),fHistoFile(0),fMod(module),fWriteToFile(kFALSE)
58 {
59   // Constructor. 
60   // oldTimeEnergy is an array of histograms kept from the previous run.
61   // By default the final histograms will not be saved to the root file.
62
63   char name[128];
64   sprintf(name,"PHOS_Module%d_Calib",fMod);
65   SetName(name);
66   
67   SetTitle("Calibration Detector Algorithm");
68   
69   if(oldTimeEnergy) {
70
71     for(Int_t iX=0; iX<64; iX++) {
72       for(Int_t iZ=0; iZ<56; iZ++) {
73         for(Int_t iGain=0; iGain<2; iGain++) {
74
75           if((&(oldTimeEnergy[iX][iZ][iGain]))->GetEntries()) 
76             fTimeEnergy[iX][iZ][iGain] = &(oldTimeEnergy[iX][iZ][iGain]);
77           else
78             fTimeEnergy[iX][iZ][iGain] = 0;
79         }
80       }
81     }
82   
83   }
84
85 }
86
87 //-------------------------------------------------------------------
88 AliPHOSDA1::AliPHOSDA1(const AliPHOSDA1& da) : TNamed(da),
89   fHistoFile(0),fMod(da.fMod),fWriteToFile(da.fWriteToFile)
90 {
91   // Copy constructor.
92
93   fHistoFile = new TFile(da.GetName(),"update");
94
95   char hname[128];
96   TH1F* hist1=0;
97   TH2F* hist2=0;
98
99   for(Int_t iX=0; iX<64; iX++) {
100     for(Int_t iZ=0; iZ<56; iZ++) {
101
102       sprintf(hname,"%d_%d_%d",fMod,iX,iZ);
103       hist1 = (TH1F*)da.fHistoFile->Get(hname);
104       if(hist1) fHgLgRatio[iX][iZ] = hist1;
105       else
106         fHgLgRatio[iX][iZ] = 0;
107
108       for(Int_t iGain=0; iGain<2; iGain++) {
109         sprintf(hname,"%d_%d_%d_%d",fMod,iX,iZ,iGain);
110         hist2 = (TH2F*)da.fHistoFile->Get(hname);
111         if(hist2) fTimeEnergy[iX][iZ][iGain] = hist2;
112         else
113           fTimeEnergy[iX][iZ][iGain] = 0;
114       }
115
116     }
117   }
118   
119 }
120
121 //-------------------------------------------------------------------
122 AliPHOSDA1& AliPHOSDA1::operator= (const AliPHOSDA1& da)
123 {
124   //Assignment operator.
125
126   if(this != &da) {
127
128     TString oldname(fHistoFile->GetName());
129     TString newname(da.fHistoFile->GetName());
130
131     if(oldname != newname) {
132       delete fHistoFile;
133       fHistoFile = new TFile(da.fHistoFile->GetName(),"update");
134     }
135
136     fMod = da.fMod;
137
138     SetName(da.GetName());
139     SetTitle(da.GetTitle());
140
141     for(Int_t iX=0; iX<64; iX++) {
142       for(Int_t iZ=0; iZ<56; iZ++) {
143
144         if(fHgLgRatio[iX][iZ]) delete fHgLgRatio[iX][iZ];
145         fHgLgRatio[iX][iZ] = da.fHgLgRatio[iX][iZ];
146
147         for(Int_t iGain=0; iGain<2; iGain++) {
148           if (fTimeEnergy[iX][iZ][iGain]) delete fTimeEnergy[iX][iZ][iGain];
149           fTimeEnergy[iX][iZ][iGain] = da.fTimeEnergy[iX][iZ][iGain];
150         }
151
152       }
153     }
154     
155   }
156
157   return *this;
158 }
159
160
161 //-------------------------------------------------------------------
162 AliPHOSDA1::~AliPHOSDA1()
163 {
164   // Destructor
165   
166   UpdateHistoFile();
167   if(fHistoFile) delete fHistoFile;
168   
169 }
170
171 //-------------------------------------------------------------------
172 void AliPHOSDA1::FillHistograms(Float_t e[64][56][2], Float_t t[64][56][2]) 
173 {
174   // Fill energy vs time-of-flight and HG/LG ratio histograms in one event.
175   // Energy and TOF data are encoded as e[X][Z][gain] and t[X][Z][gain], 
176   // where X(0..63) and Z(0..55) - crystal position in the module, 
177   // gain=0 - low gain, gain=1 - high gain.
178   // Energy in ADC counts, time in samples.
179   // If no energy or time read for particular channel, 
180   // the correspondent array entry should be filled by zero.
181   // WARNING: this function should be called once per event!
182
183   char hname[128];
184   char htitl[128];
185
186   for(Int_t iX=0; iX<64; iX++) {
187     for (Int_t iZ=0; iZ<56; iZ++) {
188
189       // HG/LG
190       if(e[iX][iZ][0] && e[iX][iZ][1]) {
191         if(fHgLgRatio[iX][iZ]) {
192 //        printf("iX=%d iZ=%d,e[iX][iZ][1]=%.3f,e[iX][iZ][0]=%.3f, t1=%.3f,t0=%.3f\n",
193 //               iX,iZ,e[iX][iZ][1],e[iX][iZ][0], t[iX][iZ][1],t[iX][iZ][0]);
194           fHgLgRatio[iX][iZ]->Fill(e[iX][iZ][1]/e[iX][iZ][0]); 
195         }
196         else
197           {
198             sprintf(hname,"%d_%d_%d",fMod,iX,iZ);
199             sprintf(htitl,"HG/LG ratio for crystal %d_%d_%d",fMod,iX,iZ);
200             fHgLgRatio[iX][iZ] = new TH1F(hname,htitl,400,14.,18.);
201 //          printf("iX=%d iZ=%d,e[iX][iZ][1]=%.3f,e[iX][iZ][0]=%.3f\n",
202 //                 iX,iZ,e[iX][iZ][1],e[iX][iZ][0]);
203             fHgLgRatio[iX][iZ]->Fill(e[iX][iZ][1]/e[iX][iZ][0]);
204           }
205       }
206           
207       // Energy vs TOF
208       for(Int_t iGain=0; iGain<2; iGain++) {
209         if(!e[iX][iZ][iGain]) continue;
210
211         if(fTimeEnergy[iX][iZ][iGain]) 
212           fTimeEnergy[iX][iZ][iGain]->Fill(e[iX][iZ][iGain],t[iX][iZ][iGain]);
213         else {
214           sprintf(hname,"%d_%d_%d_%d",fMod,iX,iZ,iGain);
215           sprintf(htitl,"Energy vs TOF for crystal %d_%d_%d and gain %d",fMod,iX,iZ,iGain);
216           fTimeEnergy[iX][iZ][iGain] = new TH2F(hname,htitl,1024,0.,1024.,100,0.,100);
217           fTimeEnergy[iX][iZ][iGain]->Fill(e[iX][iZ][iGain],t[iX][iZ][iGain]);
218         }
219       }
220
221     }
222   }
223
224 }
225
226 //-------------------------------------------------------------------
227 void AliPHOSDA1::UpdateHistoFile()
228 {
229   // Write histograms to file
230
231   if(!fWriteToFile) return;
232   if(!fHistoFile) return;
233   if(!fHistoFile->IsOpen()) return;
234
235   TH1F* hist1=0;
236   TH2F* hist2=0;
237
238   for(Int_t iX=0; iX<64; iX++) {
239     for(Int_t iZ=0; iZ<56; iZ++) {
240
241       hist1 = fHgLgRatio[iX][iZ];
242       if(hist1) hist1->Write(hist1->GetName(),TObject::kWriteDelete);
243
244       for(Int_t iGain=0; iGain<2; iGain++) {
245         hist2 = fTimeEnergy[iX][iZ][iGain];
246         if(hist2) hist2->Write(hist2->GetName(),TObject::kWriteDelete);
247       } 
248
249     }
250   }
251
252 }
253
254 //-----------------------------------------------------------------------
255 void AliPHOSDA1::SetWriteToFile(Bool_t write)
256 {
257   if(!write) { 
258     fWriteToFile = write;
259     return;
260   }
261
262   if(!fHistoFile) {
263     char rootname[128];
264     sprintf(rootname,"%s.root",GetName());
265     fHistoFile =  new TFile(rootname,"update");
266   } 
267   
268   fWriteToFile = write;
269 }