]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ESDCheck/AliPHOSQATask.cxx
Updating the AliAnalysisFemtoTask to work with
[u/mrichter/AliRoot.git] / ESDCheck / AliPHOSQATask.cxx
CommitLineData
1dfe075f 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 **************************************************************************/
0b28fd57 15
16/* $Id$ */
17
1dfe075f 18//_________________________________________________________________________
19// An analysis task to check the PHOS photon data in simulated data
d40a7452 20// An analysis task to check the PHOS photon data in simulated data
21// An analysis task to check the PHOS photon data in simulated data
22// An analysis task to check the PHOS photon data in simulated data
23// An analysis task to check the PHOS photon data in simulated data
24// An analysis task to check the PHOS photon data in simulated data
1dfe075f 25//
26//*-- Yves Schutz
27//////////////////////////////////////////////////////////////////////////////
28
0b28fd57 29#include <TCanvas.h>
1dfe075f 30#include <TChain.h>
0b28fd57 31#include <TFile.h>
1dfe075f 32#include <TH1.h>
0b28fd57 33#include <TNtuple.h>
34#include <TROOT.h>
1dfe075f 35#include <TVector3.h>
1f588058 36#include <TString.h>
1dfe075f 37
38#include "AliPHOSQATask.h"
39#include "AliESD.h"
40#include "AliLog.h"
41
42//______________________________________________________________________________
43AliPHOSQATask::AliPHOSQATask(const char *name) :
44 AliAnalysisTask(name,""),
45 fChain(0),
46 fESD(0),
47 fhPHOS(0),
48 fhPHOSEnergy(0),
49 fhPHOSDigits(0),
50 fhPHOSRecParticles(0),
51 fhPHOSPhotons(0),
52 fhPHOSInvariantMass(0),
53 fhPHOSDigitsEvent(0)
54{
55 // Constructor.
56 // Input slot #0 works with an Ntuple
57 DefineInput(0, TChain::Class());
58 // Output slot #0 writes into a TH1 container
59 DefineOutput(0, TObjArray::Class()) ;
60}
61
d40a7452 62//____________________________________________________________________________
63AliPHOSQATask::AliPHOSQATask(const AliPHOSQATask& ta) :
64 AliAnalysisTask(ta.GetName(),""),
65 fChain(ta.fChain),
66 fESD(ta.fESD),
67 fOutputContainer(ta.fOutputContainer),
68 fhPHOSPos(ta.fhPHOSPos),
69 fhPHOS(ta.fhPHOS),
70 fhPHOSEnergy(ta.fhPHOSEnergy),
71 fhPHOSDigits(ta.fhPHOSDigits),
72 fhPHOSRecParticles(ta.fhPHOSRecParticles),
73 fhPHOSPhotons(ta.fhPHOSPhotons),
74 fhPHOSInvariantMass(ta.fhPHOSInvariantMass),
75 fhPHOSDigitsEvent(ta.fhPHOSDigitsEvent)
76{
77 // cpy ctor
78}
79
80//_____________________________________________________________________________
81AliPHOSQATask& AliPHOSQATask::operator = (const AliPHOSQATask& ap)
82{
83// assignment operator
84
85 this->~AliPHOSQATask();
86 new(this) AliPHOSQATask(ap);
87 return *this;
88}
89
1dfe075f 90//______________________________________________________________________________
91AliPHOSQATask::~AliPHOSQATask()
92{
93 // dtor
94 fOutputContainer->Clear() ;
95 delete fOutputContainer ;
96
97 delete fhPHOSPos ;
98 delete fhPHOS ;
99 delete fhPHOSEnergy ;
100 delete fhPHOSDigits ;
101 delete fhPHOSRecParticles ;
102 delete fhPHOSPhotons ;
103 delete fhPHOSInvariantMass ;
104 delete fhPHOSDigitsEvent ;
105}
106
107
108//______________________________________________________________________________
c52c2132 109void AliPHOSQATask::ConnectInputData(const Option_t*)
1dfe075f 110{
111 // Initialisation of branch container and histograms
112
113 AliInfo(Form("*** Initialization of %s", GetName())) ;
114
115 // Get input data
116 fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
117 if (!fChain) {
118 AliError(Form("Input 0 for %s not found\n", GetName()));
119 return ;
120 }
121
c52c2132 122 // One should first check if the branch address was taken by some other task
123 char ** address = (char **)GetBranchAddress(0, "ESD");
124 if (address) {
125 fESD = (AliESD*)(*address);
126 } else {
127 fESD = new AliESD();
128 SetBranchAddress(0, "ESD", &fESD);
1dfe075f 129 }
c52c2132 130}
131
132//________________________________________________________________________
133void AliPHOSQATask::CreateOutputObjects()
134{
1dfe075f 135 // create histograms
136
1e20f195 137 OpenFile(0) ;
138
1dfe075f 139 fhPHOSPos = new TNtuple("PHOSPos" , "Position in PHOS" , "x:y:z");
140 fhPHOS = new TNtuple("PHOS" , "PHOS" , "event:digits:clusters:photons");
141 fhPHOSEnergy = new TH1D("PHOSEnergy" , "PHOSEnergy" , 1000, 0., 10. ) ;
142 fhPHOSDigits = new TH1I("PHOSDigitsCluster" , "PHOSDigits" , 20 , 0 , 20 ) ;
143 fhPHOSRecParticles = new TH1D("PHOSRecParticles" , "PHOSRecParticles" , 20 , 0., 20. ) ;
144 fhPHOSPhotons = new TH1I("PHOSPhotons" , "PHOSPhotons" , 20 , 0 , 20 ) ;
145 fhPHOSInvariantMass = new TH1D("PHOSInvariantMass" , "PHOSInvariantMass" , 400, 0., 400.) ;
146 fhPHOSDigitsEvent = new TH1I("PHOSDigitsEvent" , "PHOSDigitsEvent" , 30 , 0 , 30 ) ;
147
148 // create output container
149
150 fOutputContainer = new TObjArray(8) ;
151 fOutputContainer->SetName(GetName()) ;
152
153 fOutputContainer->AddAt(fhPHOSPos, 0) ;
154 fOutputContainer->AddAt(fhPHOS, 1) ;
155 fOutputContainer->AddAt(fhPHOSEnergy, 2) ;
156 fOutputContainer->AddAt(fhPHOSDigits, 3) ;
157 fOutputContainer->AddAt(fhPHOSRecParticles, 4) ;
158 fOutputContainer->AddAt(fhPHOSPhotons, 5) ;
159 fOutputContainer->AddAt(fhPHOSInvariantMass, 6) ;
160 fOutputContainer->AddAt(fhPHOSDigitsEvent, 7) ;
161}
162
163//______________________________________________________________________________
164void AliPHOSQATask::Exec(Option_t *)
165{
166 // Processing of one event
167
168 Long64_t entry = fChain->GetReadEntry() ;
169
170 if (!fESD) {
171 AliError("fESD is not connected to the input!") ;
172 return ;
173 }
174
175 if ( !((entry-1)%100) )
176 AliInfo(Form("%s ----> Processing event # %lld", (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ;
177
178 //************************ PHOS *************************************
179
180 Int_t firstPhosCluster = fESD->GetFirstPHOSCluster() ;
d40a7452 181 const Int_t kNumberOfPhosClusters = fESD->GetNumberOfPHOSClusters() ;
1dfe075f 182
d40a7452 183 TVector3 ** phosVector = new TVector3*[kNumberOfPhosClusters] ;
184 Float_t * phosPhotonsEnergy = new Float_t[kNumberOfPhosClusters] ;
1dfe075f 185 Int_t phosCluster ;
186 Int_t numberOfPhotonsInPhos = 0 ;
187 Int_t numberOfDigitsInPhos = 0 ;
188
189 // loop over the PHOS Cluster
d40a7452 190 for(phosCluster = firstPhosCluster ; phosCluster < firstPhosCluster + kNumberOfPhosClusters ; phosCluster++) {
1dfe075f 191 AliESDCaloCluster * caloCluster = fESD->GetCaloCluster(phosCluster) ;
192 if (caloCluster) {
193 Float_t pos[3] ;
1f588058 194 caloCluster->GetPosition( pos ) ;
195 fhPHOSEnergy->Fill( caloCluster->E() ) ;
1dfe075f 196 fhPHOSPos->Fill( pos[0], pos[1], pos[2] ) ;
197 fhPHOSDigits->Fill(entry, caloCluster->GetNumberOfDigits() ) ;
198 numberOfDigitsInPhos += caloCluster->GetNumberOfDigits() ;
f941f7fa 199 Double_t * pid = caloCluster->GetPid() ;
1dfe075f 200 if(pid[AliPID::kPhoton] > 0.9) {
201 phosVector[numberOfPhotonsInPhos] = new TVector3(pos[0],pos[1],pos[2]) ;
1f588058 202 phosPhotonsEnergy[numberOfPhotonsInPhos]=caloCluster->E() ;
1dfe075f 203 numberOfPhotonsInPhos++;
204 }
205 }
206 } //PHOS clusters
207
d40a7452 208 fhPHOSRecParticles->Fill(kNumberOfPhosClusters);
1dfe075f 209 fhPHOSPhotons->Fill(numberOfPhotonsInPhos);
210 fhPHOSDigitsEvent->Fill(numberOfDigitsInPhos);
d40a7452 211 fhPHOS->Fill(entry, numberOfDigitsInPhos, kNumberOfPhosClusters, numberOfPhotonsInPhos) ;
1dfe075f 212
213 // invariant Mass
214 if (numberOfPhotonsInPhos > 1 ) {
215 Int_t phosPhoton1, phosPhoton2 ;
216 for(phosPhoton1 = 0 ; phosPhoton1 < numberOfPhotonsInPhos ; phosPhoton1++) {
217 for(phosPhoton2 = phosPhoton1 + 1 ; phosPhoton2 < numberOfPhotonsInPhos ; phosPhoton2++) {
218 Float_t tempMass = TMath::Sqrt( 2 * phosPhotonsEnergy[phosPhoton1] * phosPhotonsEnergy[phosPhoton2] *
219 ( 1 - TMath::Cos(phosVector[phosPhoton1]->Angle(*phosVector[phosPhoton2])) )
220 );
221 fhPHOSInvariantMass->Fill(tempMass*1000.);
222 }
223 }
224 }
225
226 PostData(0, fOutputContainer);
227
228 delete [] phosVector ;
229 delete [] phosPhotonsEnergy ;
230
231}
232
233//______________________________________________________________________________
234void AliPHOSQATask::Terminate(Option_t *)
235{
236 // Processing when the event loop is ended
c52c2132 237
238 fOutputContainer = (TObjArray*)GetOutputData(0);
239 fhPHOSPos = (TNtuple*)fOutputContainer->At(0);
240 fhPHOS = (TNtuple*)fOutputContainer->At(1);
241 fhPHOSEnergy = (TH1D*)fOutputContainer->At(2);
242 fhPHOSDigits = (TH1I*)fOutputContainer->At(3);
243 fhPHOSRecParticles = (TH1D*)fOutputContainer->At(4);
244 fhPHOSPhotons = (TH1I*)fOutputContainer->At(5);
245 fhPHOSInvariantMass = (TH1D*)fOutputContainer->At(6);
246 fhPHOSDigitsEvent = (TH1I*)fOutputContainer->At(7);
247
2704006a 248 Bool_t problem = kFALSE ;
84eb42a1 249 AliInfo(Form(" *** %s Report:", GetName())) ;
250 printf(" PHOSEnergy Mean : %5.3f , RMS : %5.3f \n", fhPHOSEnergy->GetMean(), fhPHOSEnergy->GetRMS() ) ;
251 printf(" PHOSDigits Mean : %5.3f , RMS : %5.3f \n", fhPHOSDigits->GetMean(), fhPHOSDigits->GetRMS() ) ;
252 printf(" PHOSRecParticles Mean : %5.3f , RMS : %5.3f \n", fhPHOSRecParticles->GetMean(), fhPHOSRecParticles->GetRMS() ) ;
253 printf(" PHOSPhotons Mean : %5.3f , RMS : %5.3f \n", fhPHOSPhotons->GetMean(), fhPHOSPhotons->GetRMS() ) ;
254 printf(" PHOSInvariantMass Mean : %5.3f , RMS : %5.3f \n", fhPHOSInvariantMass->GetMean(), fhPHOSInvariantMass->GetRMS() ) ;
255 printf(" PHOSDigitsEvent Mean : %5.3f , RMS : %5.3f \n", fhPHOSDigitsEvent->GetMean(), fhPHOSDigitsEvent->GetRMS() ) ;
1dfe075f 256
257 TCanvas * cPHOS = new TCanvas("cPHOS", "PHOS ESD Test", 400, 10, 600, 700) ;
258 cPHOS->Divide(3, 2);
259
260 cPHOS->cd(1) ;
84eb42a1 261 if ( fhPHOSEnergy->GetMaximum() > 0. )
262 gPad->SetLogy();
1dfe075f 263 fhPHOSEnergy->SetAxisRange(0, 25.);
264 fhPHOSEnergy->SetLineColor(2);
265 fhPHOSEnergy->Draw();
266
267 cPHOS->cd(2) ;
268 fhPHOSDigits->SetAxisRange(0,25.);
269 fhPHOSDigits->SetLineColor(2);
270 fhPHOSDigits->Draw();
271
272 cPHOS->cd(3) ;
84eb42a1 273 if ( fhPHOSRecParticles->GetMaximum() > 0. )
274 gPad->SetLogy();
1dfe075f 275 fhPHOSRecParticles->SetAxisRange(0, 25.);
276 fhPHOSRecParticles->SetLineColor(2);
277 fhPHOSRecParticles->Draw();
278
279 cPHOS->cd(4) ;
84eb42a1 280 if ( fhPHOSPhotons->GetMaximum() > 0. )
281 gPad->SetLogy();
1dfe075f 282 fhPHOSPhotons->SetAxisRange(0,25.);
283 fhPHOSPhotons->SetLineColor(2);
284 fhPHOSPhotons->Draw();
285
286 cPHOS->cd(5) ;
287 fhPHOSInvariantMass->SetLineColor(2);
288 fhPHOSInvariantMass->Draw();
289
290 cPHOS->cd(6) ;
84eb42a1 291 if ( fhPHOSDigitsEvent->GetMaximum() > 0. )
292 gPad->SetLogy();
1dfe075f 293 fhPHOSDigitsEvent->SetAxisRange(0,40.);
294 fhPHOSDigitsEvent->SetLineColor(2);
295 fhPHOSDigitsEvent->Draw();
296
297 cPHOS->Print("PHOS.eps");
298
299 char line[1024] ;
84eb42a1 300 sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ;
1dfe075f 301 gROOT->ProcessLine(line);
302 sprintf(line, ".!rm -fR *.eps");
303 gROOT->ProcessLine(line);
304
2704006a 305 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
306
1f588058 307 TString report ;
2704006a 308 if(problem)
309 report="Problems found, please check!!!";
310 else
311 report="OK";
312
1f588058 313 AliInfo(Form("*** %s Summary Report: %s \n",GetName(), report.Data())) ;
1dfe075f 314}