]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSDA1.cxx
Update to new method names in AliESDv0, load them in visscan_init.C.
[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)
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(const AliPHOSDA1& da) : TNamed(da),
57   fHistoFile(0),fMod(da.fMod)
58 {
59   // Copy constructor.
60
61   fHistoFile = new TFile(da.GetName(),"update");
62
63   char hname[128];
64   TH1F* hist1=0;
65   TH2F* hist2=0;
66
67   for(Int_t iX=0; iX<64; iX++) {
68     for(Int_t iZ=0; iZ<56; iZ++) {
69
70       sprintf(hname,"%d_%d_%d",fMod,iX,iZ);
71       hist1 = (TH1F*)da.fHistoFile->Get(hname);
72       if(hist1) fHgLgRatio[iX][iZ] = hist1;
73       else
74         fHgLgRatio[iX][iZ] = 0;
75
76       for(Int_t iGain=0; iGain<2; iGain++) {
77         sprintf(hname,"%d_%d_%d_%d",fMod,iX,iZ,iGain);
78         hist2 = (TH2F*)da.fHistoFile->Get(hname);
79         if(hist2) fTimeEnergy[iX][iZ][iGain] = hist2;
80         else
81           fTimeEnergy[iX][iZ][iGain] = 0;
82       }
83
84     }
85   }
86   
87 }
88
89 //-------------------------------------------------------------------
90 AliPHOSDA1& AliPHOSDA1::operator= (const AliPHOSDA1& da)
91 {
92   //Assignment operator.
93
94   if(this != &da) {
95
96     TString oldname(fHistoFile->GetName());
97     TString newname(da.fHistoFile->GetName());
98
99     if(oldname != newname) {
100       delete fHistoFile;
101       fHistoFile = new TFile(da.fHistoFile->GetName(),"update");
102     }
103
104     fMod = da.fMod;
105
106     SetName(da.GetName());
107     SetTitle(da.GetTitle());
108
109     for(Int_t iX=0; iX<64; iX++) {
110       for(Int_t iZ=0; iZ<56; iZ++) {
111
112         if(fHgLgRatio[iX][iZ]) delete fHgLgRatio[iX][iZ];
113         fHgLgRatio[iX][iZ] = da.fHgLgRatio[iX][iZ];
114
115         for(Int_t iGain=0; iGain<2; iGain++) {
116           if (fTimeEnergy[iX][iZ][iGain]) delete fTimeEnergy[iX][iZ][iGain];
117           fTimeEnergy[iX][iZ][iGain] = da.fTimeEnergy[iX][iZ][iGain];
118         }
119
120       }
121     }
122     
123   }
124
125   return *this;
126 }
127
128
129 //-------------------------------------------------------------------
130 AliPHOSDA1::~AliPHOSDA1()
131 {
132   // Destructor
133   
134   UpdateHistoFile();
135   if(fHistoFile) delete fHistoFile;
136   
137 }
138
139 //-------------------------------------------------------------------
140 void AliPHOSDA1::FillHistograms(Float_t e[64][56][2], Float_t t[64][56][2]) 
141 {
142   // Fill energy vs time-of-flight and HG/LG ratio histograms in one event.
143   // Energy and TOF data are encoded as e[X][Z][gain] and t[X][Z][gain], 
144   // where X(0..63) and Z(0..55) - crystal position in the module, 
145   // gain=0 - low gain, gain=1 - high gain.
146   // Energy in ADC counts, time in samples.
147   // If no energy or time read for particular channel, 
148   // the correspondent array entry should be filled by zero.
149   // WARNING: this function should be called once per event!
150
151   char hname[128];
152   char htitl[128];
153
154   for(Int_t iX=0; iX<64; iX++) {
155     for (Int_t iZ=0; iZ<56; iZ++) {
156
157       // HG/LG
158       if(e[iX][iZ][0] && e[iX][iZ][1]) {
159         if(fHgLgRatio[iX][iZ]) {
160 //        printf("iX=%d iZ=%d,e[iX][iZ][1]=%.3f,e[iX][iZ][0]=%.3f, t1=%.3f,t0=%.3f\n",
161 //               iX,iZ,e[iX][iZ][1],e[iX][iZ][0], t[iX][iZ][1],t[iX][iZ][0]);
162           fHgLgRatio[iX][iZ]->Fill(e[iX][iZ][1]/e[iX][iZ][0]); 
163         }
164         else
165           {
166             sprintf(hname,"%d_%d_%d",fMod,iX,iZ);
167             sprintf(htitl,"HG/LG ratio for crystal %d_%d_%d",fMod,iX,iZ);
168             fHgLgRatio[iX][iZ] = new TH1F(hname,htitl,400,14.,18.);
169 //          printf("iX=%d iZ=%d,e[iX][iZ][1]=%.3f,e[iX][iZ][0]=%.3f\n",
170 //                 iX,iZ,e[iX][iZ][1],e[iX][iZ][0]);
171             fHgLgRatio[iX][iZ]->Fill(e[iX][iZ][1]/e[iX][iZ][0]);
172           }
173       }
174           
175       // Energy vs TOF
176       for(Int_t iGain=0; iGain<2; iGain++) {
177         if(!e[iX][iZ][iGain]) continue;
178
179         if(fTimeEnergy[iX][iZ][iGain]) 
180           fTimeEnergy[iX][iZ][iGain]->Fill(e[iX][iZ][iGain],t[iX][iZ][iGain]);
181         else {
182           sprintf(hname,"%d_%d_%d_%d",fMod,iX,iZ,iGain);
183           sprintf(htitl,"Energy vs TOF for crystal %d_%d_%d and gain %d",fMod,iX,iZ,iGain);
184           fTimeEnergy[iX][iZ][iGain] = new TH2F(hname,htitl,1024,0.,1024.,100,0.,100);
185           fTimeEnergy[iX][iZ][iGain]->Fill(e[iX][iZ][iGain],t[iX][iZ][iGain]);
186         }
187       }
188
189     }
190   }
191
192 }
193
194 //-------------------------------------------------------------------
195 void AliPHOSDA1::UpdateHistoFile()
196 {
197   // Write histograms to file
198
199   if(!fHistoFile) return;
200   if(!fHistoFile->IsOpen()) return;
201
202   TH1F* hist1=0;
203   TH2F* hist2=0;
204
205   for(Int_t iX=0; iX<64; iX++) {
206     for(Int_t iZ=0; iZ<56; iZ++) {
207
208       hist1 = fHgLgRatio[iX][iZ];
209       if(hist1) hist1->Write(hist1->GetName(),TObject::kWriteDelete);
210
211       for(Int_t iGain=0; iGain<2; iGain++) {
212         hist2 = fTimeEnergy[iX][iZ][iGain];
213         if(hist2) hist2->Write(hist2->GetName(),TObject::kWriteDelete);
214       } 
215
216     }
217   }
218
219 }
220