]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ESDCheck/AliPHOSQATask.cxx
Adding includes now needed by ROOT
[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
20//
21//*-- Yves Schutz
22//////////////////////////////////////////////////////////////////////////////
23
0b28fd57 24#include <TCanvas.h>
1dfe075f 25#include <TChain.h>
0b28fd57 26#include <TFile.h>
1dfe075f 27#include <TH1.h>
28#include <TH1F.h>
29#include <TH1I.h>
1dfe075f 30#include <TLegend.h>
0b28fd57 31#include <TNtuple.h>
32#include <TROOT.h>
1dfe075f 33#include <TVector3.h>
1dfe075f 34
35#include "AliPHOSQATask.h"
36#include "AliESD.h"
37#include "AliLog.h"
38
39//______________________________________________________________________________
40AliPHOSQATask::AliPHOSQATask(const char *name) :
41 AliAnalysisTask(name,""),
42 fChain(0),
43 fESD(0),
44 fhPHOS(0),
45 fhPHOSEnergy(0),
46 fhPHOSDigits(0),
47 fhPHOSRecParticles(0),
48 fhPHOSPhotons(0),
49 fhPHOSInvariantMass(0),
50 fhPHOSDigitsEvent(0)
51{
52 // Constructor.
53 // Input slot #0 works with an Ntuple
54 DefineInput(0, TChain::Class());
55 // Output slot #0 writes into a TH1 container
56 DefineOutput(0, TObjArray::Class()) ;
57}
58
59//______________________________________________________________________________
60AliPHOSQATask::~AliPHOSQATask()
61{
62 // dtor
63 fOutputContainer->Clear() ;
64 delete fOutputContainer ;
65
66 delete fhPHOSPos ;
67 delete fhPHOS ;
68 delete fhPHOSEnergy ;
69 delete fhPHOSDigits ;
70 delete fhPHOSRecParticles ;
71 delete fhPHOSPhotons ;
72 delete fhPHOSInvariantMass ;
73 delete fhPHOSDigitsEvent ;
74}
75
76
77//______________________________________________________________________________
78void AliPHOSQATask::Init(const Option_t*)
79{
80 // Initialisation of branch container and histograms
81
82 AliInfo(Form("*** Initialization of %s", GetName())) ;
83
84 // Get input data
85 fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
86 if (!fChain) {
87 AliError(Form("Input 0 for %s not found\n", GetName()));
88 return ;
89 }
90
91 if (!fESD) {
92 // One should first check if the branch address was taken by some other task
93 char ** address = (char **)GetBranchAddress(0, "ESD") ;
94 if (address)
95 fESD = (AliESD *)(*address) ;
96 if (!fESD)
97 fChain->SetBranchAddress("ESD", &fESD) ;
98 }
99 // The output objects will be written to
100 TDirectory * cdir = gDirectory ;
101 // Open a file for output #0
102 char outputName[1024] ;
103 sprintf(outputName, "%s.root", GetName() ) ;
104 OpenFile(0, outputName , "RECREATE") ;
105 if (cdir)
106 cdir->cd() ;
107
108 // create histograms
109
110 fhPHOSPos = new TNtuple("PHOSPos" , "Position in PHOS" , "x:y:z");
111 fhPHOS = new TNtuple("PHOS" , "PHOS" , "event:digits:clusters:photons");
112 fhPHOSEnergy = new TH1D("PHOSEnergy" , "PHOSEnergy" , 1000, 0., 10. ) ;
113 fhPHOSDigits = new TH1I("PHOSDigitsCluster" , "PHOSDigits" , 20 , 0 , 20 ) ;
114 fhPHOSRecParticles = new TH1D("PHOSRecParticles" , "PHOSRecParticles" , 20 , 0., 20. ) ;
115 fhPHOSPhotons = new TH1I("PHOSPhotons" , "PHOSPhotons" , 20 , 0 , 20 ) ;
116 fhPHOSInvariantMass = new TH1D("PHOSInvariantMass" , "PHOSInvariantMass" , 400, 0., 400.) ;
117 fhPHOSDigitsEvent = new TH1I("PHOSDigitsEvent" , "PHOSDigitsEvent" , 30 , 0 , 30 ) ;
118
119 // create output container
120
121 fOutputContainer = new TObjArray(8) ;
122 fOutputContainer->SetName(GetName()) ;
123
124 fOutputContainer->AddAt(fhPHOSPos, 0) ;
125 fOutputContainer->AddAt(fhPHOS, 1) ;
126 fOutputContainer->AddAt(fhPHOSEnergy, 2) ;
127 fOutputContainer->AddAt(fhPHOSDigits, 3) ;
128 fOutputContainer->AddAt(fhPHOSRecParticles, 4) ;
129 fOutputContainer->AddAt(fhPHOSPhotons, 5) ;
130 fOutputContainer->AddAt(fhPHOSInvariantMass, 6) ;
131 fOutputContainer->AddAt(fhPHOSDigitsEvent, 7) ;
132}
133
134//______________________________________________________________________________
135void AliPHOSQATask::Exec(Option_t *)
136{
137 // Processing of one event
138
139 Long64_t entry = fChain->GetReadEntry() ;
140
141 if (!fESD) {
142 AliError("fESD is not connected to the input!") ;
143 return ;
144 }
145
146 if ( !((entry-1)%100) )
147 AliInfo(Form("%s ----> Processing event # %lld", (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ;
148
149 //************************ PHOS *************************************
150
151 Int_t firstPhosCluster = fESD->GetFirstPHOSCluster() ;
152 const Int_t numberOfPhosClusters = fESD->GetNumberOfPHOSClusters() ;
153
154 TVector3 ** phosVector = new TVector3*[numberOfPhosClusters] ;
155 Float_t * phosPhotonsEnergy = new Float_t[numberOfPhosClusters] ;
156 Int_t phosCluster ;
157 Int_t numberOfPhotonsInPhos = 0 ;
158 Int_t numberOfDigitsInPhos = 0 ;
159
160 // loop over the PHOS Cluster
161 for(phosCluster = firstPhosCluster ; phosCluster < firstPhosCluster + numberOfPhosClusters ; phosCluster++) {
162 AliESDCaloCluster * caloCluster = fESD->GetCaloCluster(phosCluster) ;
163 if (caloCluster) {
164 Float_t pos[3] ;
165 caloCluster->GetGlobalPosition( pos ) ;
166 fhPHOSEnergy->Fill( caloCluster->GetClusterEnergy() ) ;
167 fhPHOSPos->Fill( pos[0], pos[1], pos[2] ) ;
168 fhPHOSDigits->Fill(entry, caloCluster->GetNumberOfDigits() ) ;
169 numberOfDigitsInPhos += caloCluster->GetNumberOfDigits() ;
170 Float_t * pid = caloCluster->GetPid() ;
171 if(pid[AliPID::kPhoton] > 0.9) {
172 phosVector[numberOfPhotonsInPhos] = new TVector3(pos[0],pos[1],pos[2]) ;
173 phosPhotonsEnergy[numberOfPhotonsInPhos]=caloCluster->GetClusterEnergy() ;
174 numberOfPhotonsInPhos++;
175 }
176 }
177 } //PHOS clusters
178
179 fhPHOSRecParticles->Fill(numberOfPhosClusters);
180 fhPHOSPhotons->Fill(numberOfPhotonsInPhos);
181 fhPHOSDigitsEvent->Fill(numberOfDigitsInPhos);
182 fhPHOS->Fill(entry, numberOfDigitsInPhos, numberOfPhosClusters, numberOfPhotonsInPhos) ;
183
184 // invariant Mass
185 if (numberOfPhotonsInPhos > 1 ) {
186 Int_t phosPhoton1, phosPhoton2 ;
187 for(phosPhoton1 = 0 ; phosPhoton1 < numberOfPhotonsInPhos ; phosPhoton1++) {
188 for(phosPhoton2 = phosPhoton1 + 1 ; phosPhoton2 < numberOfPhotonsInPhos ; phosPhoton2++) {
189 Float_t tempMass = TMath::Sqrt( 2 * phosPhotonsEnergy[phosPhoton1] * phosPhotonsEnergy[phosPhoton2] *
190 ( 1 - TMath::Cos(phosVector[phosPhoton1]->Angle(*phosVector[phosPhoton2])) )
191 );
192 fhPHOSInvariantMass->Fill(tempMass*1000.);
193 }
194 }
195 }
196
197 PostData(0, fOutputContainer);
198
199 delete [] phosVector ;
200 delete [] phosPhotonsEnergy ;
201
202}
203
204//______________________________________________________________________________
205void AliPHOSQATask::Terminate(Option_t *)
206{
207 // Processing when the event loop is ended
208
209 printf("PHOSEnergy Mean : %5.3f , RMS : %5.3f \n", fhPHOSEnergy->GetMean(), fhPHOSEnergy->GetRMS() ) ;
210 printf("PHOSDigits Mean : %5.3f , RMS : %5.3f \n", fhPHOSDigits->GetMean(), fhPHOSDigits->GetRMS() ) ;
211 printf("PHOSRecParticles Mean : %5.3f , RMS : %5.3f \n", fhPHOSRecParticles->GetMean(), fhPHOSRecParticles->GetRMS() ) ;
212 printf("PHOSPhotons Mean : %5.3f , RMS : %5.3f \n", fhPHOSPhotons->GetMean(), fhPHOSPhotons->GetRMS() ) ;
213 printf("PHOSInvariantMass Mean : %5.3f , RMS : %5.3f \n", fhPHOSInvariantMass->GetMean(), fhPHOSInvariantMass->GetRMS() ) ;
214 printf("PHOSDigitsEvent Mean : %5.3f , RMS : %5.3f \n", fhPHOSDigitsEvent->GetMean(), fhPHOSDigitsEvent->GetRMS() ) ;
215
216 TCanvas * cPHOS = new TCanvas("cPHOS", "PHOS ESD Test", 400, 10, 600, 700) ;
217 cPHOS->Divide(3, 2);
218
219 cPHOS->cd(1) ;
220 gPad->SetLogy();
221 fhPHOSEnergy->SetAxisRange(0, 25.);
222 fhPHOSEnergy->SetLineColor(2);
223 fhPHOSEnergy->Draw();
224
225 cPHOS->cd(2) ;
226 fhPHOSDigits->SetAxisRange(0,25.);
227 fhPHOSDigits->SetLineColor(2);
228 fhPHOSDigits->Draw();
229
230 cPHOS->cd(3) ;
231 gPad->SetLogy();
232 fhPHOSRecParticles->SetAxisRange(0, 25.);
233 fhPHOSRecParticles->SetLineColor(2);
234 fhPHOSRecParticles->Draw();
235
236 cPHOS->cd(4) ;
237 gPad->SetLogy();
238 fhPHOSPhotons->SetAxisRange(0,25.);
239 fhPHOSPhotons->SetLineColor(2);
240 fhPHOSPhotons->Draw();
241
242 cPHOS->cd(5) ;
243 fhPHOSInvariantMass->SetLineColor(2);
244 fhPHOSInvariantMass->Draw();
245
246 cPHOS->cd(6) ;
247 gPad->SetLogy();
248 fhPHOSDigitsEvent->SetAxisRange(0,40.);
249 fhPHOSDigitsEvent->SetLineColor(2);
250 fhPHOSDigitsEvent->Draw();
251
252 cPHOS->Print("PHOS.eps");
253
254 char line[1024] ;
255 sprintf(line, ".!tar -zcvf %s.tar.gz *.eps", GetName()) ;
256 gROOT->ProcessLine(line);
257 sprintf(line, ".!rm -fR *.eps");
258 gROOT->ProcessLine(line);
259
260 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!! \n", GetName())) ;
261}