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