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