]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ZDC/AliZDCCalibTask.cxx
update from pr task : sjena
[u/mrichter/AliRoot.git] / ZDC / AliZDCCalibTask.cxx
CommitLineData
b81fcf27 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16#include <Riostream.h>
17#include <TROOT.h>
18#include <TSystem.h>
19#include <TInterpreter.h>
20#include <TTree.h>
21#include <TCanvas.h>
22#include <TList.h>
23#include <TH1F.h>
24#include <TFile.h>
25
26#include "AliZDCCalibTask.h"
27#include "AliESDEvent.h"
28#include "AliLog.h"
29
30ClassImp(AliZDCCalibTask)
31
32//_______________________________________________________
33AliZDCCalibTask::AliZDCCalibTask(const char *name):
34 AliAnalysisTaskSE(name),
35 fESD(0),
36 fListOfHistos(0x0),
37 fZNCHisto(0x0),
38 fZPCHisto(0x0),
39 fZNAHisto(0x0),
40 fZPAHisto(0x0)
41{
42 // Constructor
43 DefineOutput(1, TList::Class());
44
45}
46
47//_______________________________________________________
48AliZDCCalibTask::AliZDCCalibTask(const AliZDCCalibTask &calibTask):
49 AliAnalysisTaskSE("AliZDCCalibTask"),
50 fESD(calibTask.fESD),
51 fListOfHistos(calibTask.fListOfHistos),
52 fZNCHisto(calibTask.fZNCHisto),
53 fZPCHisto(calibTask.fZPCHisto),
54 fZNAHisto(calibTask.fZNAHisto),
55 fZPAHisto(calibTask.fZPAHisto)
56{
57 // Copy constructor
58
59}
60
61//______________________________________________________________________________
62AliZDCCalibTask:: ~AliZDCCalibTask()
63{
64 // destructor
65 delete fESD;
66 delete fListOfHistos;
67 delete fZNCHisto;
68 delete fZPCHisto;
69 delete fZNAHisto;
70 delete fZPAHisto;
71}
72
73//______________________________________________________________________________
74AliZDCCalibTask& AliZDCCalibTask::operator=(const AliZDCCalibTask &calibtask)
75{
76 //assignment operator
77
78 fESD = calibtask.fESD;
79 fListOfHistos = calibtask.fListOfHistos;
80 fZNCHisto = calibtask.fZNCHisto;
81 fZPCHisto = calibtask.fZPCHisto;
82 fZNAHisto = calibtask.fZNAHisto;
83 fZPAHisto = calibtask.fZPAHisto;
84
85 return *this;
86}
87
88//--------------------------------------------------------------------------
89void AliZDCCalibTask::BookHistos(){
90
91 //booking histos
92
93 AliInfo(Form("*** Booking Histograms %s", GetName())) ;
94
95 fZNCHisto = new TH1F("fZNCHisto", "Energy in ZNC",100,0.,MAXENVALUE);
96 fZNCHisto->SetXTitle("E (GeV)");
97 fZPCHisto = new TH1F("fZPCHisto", "Energy in ZPC",100,0.,MAXENVALUE);
98 fZPCHisto->SetXTitle("E (GeV)");
99 fZNAHisto = new TH1F("fZNAHisto", "Energy in ZNA",100,0.,MAXENVALUE);
100 fZNAHisto->SetXTitle("E (GeV)");
101 fZPAHisto = new TH1F("fZPAHisto", "Energy in ZPA",100,0.,MAXENVALUE);
102 fZPAHisto->SetXTitle("E (GeV)");
103
104 if(!fListOfHistos) fListOfHistos = new TList();
105 fListOfHistos->SetOwner();
106 fListOfHistos->AddAt(fZNCHisto, 0);
107 fListOfHistos->AddAt(fZPCHisto, 1);
108 fListOfHistos->AddAt(fZNAHisto, 2);
109 fListOfHistos->AddAt(fZPAHisto, 3);
110
111}
112
113//--------------------------------------------------------------------------
114void AliZDCCalibTask::DrawHistos(){
115
116 //drawing output histos
117
118 AliInfo(Form("*** Drawing Histograms %s", GetName())) ;
119
120 if(!gROOT->IsBatch()){
121 fListOfHistos = dynamic_cast<TList*>(GetOutputData(1));
122 fZNCHisto = (TH1F*) fListOfHistos->At(0);
123 fZPCHisto = (TH1F*) fListOfHistos->At(1);
124 fZNAHisto = (TH1F*) fListOfHistos->At(2);
125 fZPAHisto = (TH1F*) fListOfHistos->At(3);
126
127 TCanvas * canvPlot = new TCanvas("canvPlot","ZDC energy",0,0,600,600);
128 canvPlot->Divide(2,2);
129 canvPlot->cd(1);
130 fZNCHisto->SetLineColor(kTeal+1);
131 fZNCHisto->DrawCopy("hist");
132 canvPlot->cd(2);
133 fZPCHisto->SetLineColor(kAzure+1);
134 fZPCHisto->DrawCopy("hist");
135 canvPlot->cd(3);
136 fZNCHisto->SetLineColor(kGreen+1);
137 fZNCHisto->DrawCopy("hist");
138 canvPlot->cd(4);
139 fZPCHisto->SetLineColor(kBlue+1);
140 fZPCHisto->DrawCopy("hist");
141
142 }
143
144}
145
146//_______________________________________________________
147void AliZDCCalibTask::UserCreateOutputObjects()
148{
149 // Create histograms
150
151 if(fDebug>1) printf("ZDCCalibTask::UserCreateOutputObjects() \n");
152
153 //Open histograms
154 OpenFile(1);
155 BookHistos();
156
157}
158
159//_______________________________________________________
160void AliZDCCalibTask::UserExec(Option_t */*option*/)
161{
162 // Execute analysis
163 AliVEvent *event = InputEvent();
164 if (!event) {
165 printf("ERROR: Could not retrieve event\n");
166 return;
167 }
168 AliESDEvent* esd = dynamic_cast<AliESDEvent*>(event);
169 AliESDZDC *esdZDC = esd->GetESDZDC();
170
171 Float_t *znctow = (Float_t *) (esdZDC->GetZN1TowerEnergyLR());
172 Float_t *zpctow = (Float_t *) (esdZDC->GetZP1TowerEnergyLR());
173 Float_t *znatow = (Float_t *) (esdZDC->GetZN2TowerEnergyLR());
174 Float_t *zpatow = (Float_t *) (esdZDC->GetZP2TowerEnergyLR());
175
176 Float_t energyZNC=0., energyZPC=0.;
177 Float_t energyZNA=0., energyZPA=0.;
178 for(Int_t j=0; j<5; j++){
179 energyZNC += znctow[j];
180 energyZPC += zpctow[j];
181 energyZNA += znatow[j];
182 energyZPA += zpatow[j];
183 }
184
185 //printf(" Filling histos with values: %f %f %f %f\n",
186 // energyZNC,energyZPC,energyZNA,energyZPA);
187 //
188 fZNCHisto->Fill(energyZNC);
189 fZPCHisto->Fill(energyZPC);
190 fZNAHisto->Fill(energyZNA);
191 fZPAHisto->Fill(energyZPA);
192
193 PostData(1, fListOfHistos);
194}
195
196//_______________________________________________________
197void AliZDCCalibTask::Terminate(Option_t *)
198{
199 // Draw results
200 TH1::AddDirectory(0);
201
202 fListOfHistos = dynamic_cast<TList*>(GetOutputData(1));
203 fZNCHisto = (TH1F*) fListOfHistos->At(0);
204 fZPCHisto = (TH1F*) fListOfHistos->At(1);
205 fZNAHisto = (TH1F*) fListOfHistos->At(2);
206 fZPAHisto = (TH1F*) fListOfHistos->At(3);
207
208 Int_t binMax[4], nBinsx[4];
209 Float_t yMax[4], meanFitVal[4]={0.,0.,0.,0.};
210 TF1 *fitfun[4];
211
212 /*printf("\n # of entries in histos: %d %d %d %d \n",
213 (Int_t) fZNCHisto->GetEntries(),(Int_t) fZPCHisto->GetEntries(),
214 (Int_t) fZNAHisto->GetEntries(),(Int_t) fZPAHisto->GetEntries());*/
215
216 if(fZNCHisto->GetEntries() != 0){
217 binMax[0] = fZNCHisto->GetMaximumBin();
218 yMax[0] = (fZNCHisto->GetXaxis())->GetXmax();
219 nBinsx[0] = (fZNCHisto->GetXaxis())->GetNbins();
220 fZNCHisto->Fit("gaus","Q","",binMax[0]*yMax[0]/nBinsx[0]*0.7,binMax[0]*yMax[0]/nBinsx[0]*1.25);
221 fitfun[0] = fZNCHisto->GetFunction("gaus");
222 if(fitfun[0]) meanFitVal[0] = (Float_t) (fitfun[0]->GetParameter(1));
223 else printf(" !!! No fit on fZNCHisto !!!\n");
224 }
225 else AliWarning(" !!! ZNC empty histo !!!\n");
226
227 if(fZPCHisto->GetEntries() != 0){
228 binMax[1] = fZPCHisto->GetMaximumBin();
229 yMax[1] = (fZPCHisto->GetXaxis())->GetXmax();
230 nBinsx[1] = (fZPCHisto->GetXaxis())->GetNbins();
231 fZPCHisto->Fit("gaus","Q","",binMax[1]*yMax[1]/nBinsx[1]*0.7,binMax[1]*yMax[1]/nBinsx[1]*1.25);
232 fitfun[1] = fZPCHisto->GetFunction("gaus");
233 if(fitfun[1]) meanFitVal[1] = (Float_t) (fitfun[1]->GetParameter(1));
234 else printf(" !!! No fit on fZPCHisto !!!\n");
235 }
236 else AliWarning(" !!! ZPC empty histo !!!\n");
237
238 if(fZNAHisto->GetEntries() != 0){
239 binMax[2] = fZNAHisto->GetMaximumBin();
240 yMax[2] = (fZNAHisto->GetXaxis())->GetXmax();
241 nBinsx[2] = (fZNAHisto->GetXaxis())->GetNbins();
242 fZNAHisto->Fit("gaus","Q","",binMax[2]*yMax[2]/nBinsx[2]*0.7,binMax[2]*yMax[2]/nBinsx[2]*1.25);
243 fitfun[2] = fZNAHisto->GetFunction("gaus");
244 if(fitfun[2]) meanFitVal[2] = (Float_t) (fitfun[2]->GetParameter(1));
245 else printf(" !!! No fit on fZNAHisto !!!\n");
246 }
247 else AliWarning(" !!! ZNA empty histo !!!\n");
248
249 if(fZPAHisto->GetEntries() != 0){
250 binMax[3] = fZPAHisto->GetMaximumBin();
251 yMax[3] = (fZPAHisto->GetXaxis())->GetXmax();
252 nBinsx[3] = (fZPAHisto->GetXaxis())->GetNbins();
253 fZPAHisto->Fit("gaus","Q","",binMax[3]*yMax[3]/nBinsx[3]*0.7,binMax[3]*yMax[3]/nBinsx[3]*1.25);
254 fitfun[3] = fZPAHisto->GetFunction("gaus");
255 if(fitfun[3]) meanFitVal[3] = (Float_t) (fitfun[3]->GetParameter(1));
256 else printf(" !!! No fit on fZPAHisto !!!\n");
257 }
258 else AliWarning(" !!! ZPA empty histo !!!\n");
259
260 FILE *fileOut = fopen("EMDCalibData.dat","w");
261 for(Int_t j=0; j<6; j++){
262 if(j<4) fprintf(fileOut,"\t%f\n",meanFitVal[j]);
263 else fprintf(fileOut,"\t1.0000\n");
264 }
265 fclose(fileOut);
266
267 DrawHistos();
268
269
270
271}