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