1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //________________________________________________________________________
19 // A basic analysis task to analyse photon detected by PHOS
20 // A basic analysis task to analyse photon detected by PHOS
21 // A basic analysis task to analyse photon detected by PHOS
22 // Adapted for AliAnalysisTaskSE and AOD production
26 //////////////////////////////////////////////////////////////////////////////
35 #include "AliAnalysisTaskPHOSExample.h"
36 #include "AliESDEvent.h"
37 #include "AliESDCaloCluster.h"
38 #include "AliAODEvent.h"
39 #include "AliAODPhoton.h"
43 //______________________________________________________________________________
44 AliAnalysisTaskPHOSExample::AliAnalysisTaskPHOSExample() :
54 fhPHOSRecParticles(0),
56 fhPHOSInvariantMass(0),
61 //______________________________________________________________________________
62 AliAnalysisTaskPHOSExample::AliAnalysisTaskPHOSExample(const char *name) :
63 AliAnalysisTaskSE(name),
73 fhPHOSRecParticles(0),
75 fhPHOSInvariantMass(0),
81 DefineOutput(1, TList::Class()) ;
84 //____________________________________________________________________________
85 AliAnalysisTaskPHOSExample::AliAnalysisTaskPHOSExample(const AliAnalysisTaskPHOSExample& ap) :
86 AliAnalysisTaskSE(ap.GetName()),
88 fAODPhotons(ap.fAODPhotons),
89 fPhotonsInPhos(ap.fPhotonsInPhos),
90 fPhotonId(ap.fPhotonId),
91 fOutputList(ap.fOutputList),
92 fhPHOSPos(ap.fhPHOSPos),
94 fhPHOSEnergy(ap.fhPHOSEnergy),
95 fhPHOSDigits(ap.fhPHOSDigits),
96 fhPHOSRecParticles(ap.fhPHOSRecParticles),
97 fhPHOSPhotons(ap.fhPHOSPhotons),
98 fhPHOSInvariantMass(ap.fhPHOSInvariantMass),
99 fhPHOSDigitsEvent(ap.fhPHOSDigitsEvent)
104 //_____________________________________________________________________________
105 AliAnalysisTaskPHOSExample& AliAnalysisTaskPHOSExample::operator = (const AliAnalysisTaskPHOSExample& ap)
107 // assignment operator
109 this->~AliAnalysisTaskPHOSExample();
110 new(this) AliAnalysisTaskPHOSExample(ap);
113 fAODPhotons = ap.fAODPhotons;
114 fPhotonsInPhos = ap.fPhotonsInPhos;
115 fPhotonId = ap.fPhotonId;
116 fOutputList = ap.fOutputList;
117 fhPHOSPos = ap.fhPHOSPos;
119 fhPHOSEnergy = ap.fhPHOSEnergy;
120 fhPHOSDigits = ap.fhPHOSDigits;
121 fhPHOSRecParticles = ap.fhPHOSRecParticles;
122 fhPHOSPhotons = ap.fhPHOSPhotons;
123 fhPHOSInvariantMass = ap.fhPHOSInvariantMass;
124 fhPHOSDigitsEvent = ap.fhPHOSDigitsEvent;
129 //______________________________________________________________________________
130 AliAnalysisTaskPHOSExample::~AliAnalysisTaskPHOSExample()
134 fOutputList->Clear() ;
140 //________________________________________________________________________
141 void AliAnalysisTaskPHOSExample::UserCreateOutputObjects()
143 // Create the outputs containers
145 fAODPhotons = new TClonesArray("AliAODPhoton", 0);
146 fAODPhotons->SetName("Photons");
147 AddAODBranch("TClonesArray", &fAODPhotons);
151 fhPHOSPos = new TNtuple("PHOSPos" , "Position in PHOS" , "x:y:z");
152 fhPHOS = new TNtuple("PHOS" , "PHOS" , "event:digits:clusters:photons");
153 fhPHOSEnergy = new TH1D("PHOSEnergy" , "PHOSEnergy" , 100, 0., 100. ) ;
154 fhPHOSDigits = new TH1I("PHOSDigitsCluster" , "PHOSDigits" , 20 , 0 , 20 ) ;
155 fhPHOSRecParticles = new TH1D("PHOSRecParticles" , "PHOSRecParticles" , 20 , 0., 20. ) ;
156 fhPHOSPhotons = new TH1I("PHOSPhotons" , "PHOSPhotons" , 20 , 0 , 20 ) ;
157 fhPHOSInvariantMass = new TH1D("PHOSInvariantMass" , "PHOSInvariantMass" , 400, 0., 400.) ;
158 fhPHOSDigitsEvent = new TH1I("PHOSDigitsEvent" , "PHOSDigitsEvent" , 30 , 0 , 30 ) ;
160 // create output container
162 fOutputList = new TList() ;
163 fOutputList->SetName(GetName()) ;
165 fOutputList->AddAt(fhPHOSPos, 0) ;
166 fOutputList->AddAt(fhPHOS, 1) ;
167 fOutputList->AddAt(fhPHOSEnergy, 2) ;
168 fOutputList->AddAt(fhPHOSDigits, 3) ;
169 fOutputList->AddAt(fhPHOSRecParticles, 4) ;
170 fOutputList->AddAt(fhPHOSPhotons, 5) ;
171 fOutputList->AddAt(fhPHOSInvariantMass, 6) ;
172 fOutputList->AddAt(fhPHOSDigitsEvent, 7) ;
177 //______________________________________________________________________________
178 void AliAnalysisTaskPHOSExample::UserExec(Option_t *)
180 // Processing of one event
182 // if ( !((Entry()-1)%100) )
183 AliInfo(Form(" Processing event # %lld", Entry())) ;
184 AliESDEvent* esd = (AliESDEvent*)InputEvent();
186 //************************ PHOS *************************************
187 TRefArray * caloClustersArr = new TRefArray();
188 esd->GetPHOSClusters(caloClustersArr);
190 const Int_t kNumberOfPhosClusters = caloClustersArr->GetEntries() ;
192 TVector3 ** phosVector = new TVector3*[kNumberOfPhosClusters] ;
193 Float_t * phosPhotonsEnergy = new Float_t[kNumberOfPhosClusters] ;
195 Int_t numberOfDigitsInPhos = 0 ;
196 Double_t v[3] ; //vertex ;
197 esd->GetVertex()->GetXYZ(v) ;
200 // loop over the PHOS Cluster
201 for(phosCluster = 0 ; phosCluster < kNumberOfPhosClusters ; phosCluster++) {
202 AliESDCaloCluster * caloCluster = (AliESDCaloCluster *) caloClustersArr->At(phosCluster) ;
204 //AliESDCaloCluster * caloCluster = fESD->GetCaloCluster(phosCluster) ;
207 caloCluster->GetPosition( pos ) ;
208 fhPHOSEnergy->Fill( caloCluster->E() ) ;
209 fhPHOSPos->Fill( pos[0], pos[1], pos[2] ) ;
210 fhPHOSDigits->Fill(Entry(), caloCluster->GetNCells() ) ;
211 numberOfDigitsInPhos += caloCluster->GetNCells() ;
212 const Double_t * pid = caloCluster->GetPID() ;
213 if(pid[AliPID::kPhoton] > GetPhotonId() ) {
214 phosVector[fPhotonsInPhos] = new TVector3(pos[0],pos[1],pos[2]) ;
215 phosPhotonsEnergy[fPhotonsInPhos]=caloCluster->E() ;
216 TLorentzVector momentum ;
217 caloCluster->GetMomentum(momentum, v);
218 new ((*fAODPhotons)[fPhotonsInPhos++]) AliAODPhoton (momentum);
223 fhPHOSRecParticles->Fill(kNumberOfPhosClusters);
224 fhPHOSPhotons->Fill(fPhotonsInPhos);
225 fhPHOSDigitsEvent->Fill(numberOfDigitsInPhos);
226 fhPHOS->Fill(Entry(), numberOfDigitsInPhos, kNumberOfPhosClusters, fPhotonsInPhos) ;
229 if (fPhotonsInPhos > 1 ) {
230 Int_t phosPhoton1, phosPhoton2 ;
231 for(phosPhoton1 = 0 ; phosPhoton1 < fPhotonsInPhos ; phosPhoton1++) {
232 for(phosPhoton2 = phosPhoton1 + 1 ; phosPhoton2 < fPhotonsInPhos ; phosPhoton2++) {
233 Float_t tempMass = TMath::Sqrt( 2 * phosPhotonsEnergy[phosPhoton1] * phosPhotonsEnergy[phosPhoton2] *
234 ( 1 - TMath::Cos(phosVector[phosPhoton1]->Angle(*phosVector[phosPhoton2])) )
236 fhPHOSInvariantMass->Fill(tempMass*1000.);
241 PostData(1, fOutputList);
243 delete [] phosVector ;
244 delete [] phosPhotonsEnergy ;
249 //______________________________________________________________________________
250 void AliAnalysisTaskPHOSExample::Init()
252 // Intialisation of parameters
253 AliInfo("Doing initialisation") ;
257 //______________________________________________________________________________
258 void AliAnalysisTaskPHOSExample::Terminate(Option_t *)
260 // Processing when the event loop is ended
262 // Bool_t problem = kFALSE ;
263 AliInfo(Form(" *** %s Report:", GetName())) ;
264 printf(" PHOSEnergy Mean : %5.3f , RMS : %5.3f \n", fhPHOSEnergy->GetMean(), fhPHOSEnergy->GetRMS() ) ;
265 printf(" PHOSDigits Mean : %5.3f , RMS : %5.3f \n", fhPHOSDigits->GetMean(), fhPHOSDigits->GetRMS() ) ;
266 printf(" PHOSRecParticles Mean : %5.3f , RMS : %5.3f \n", fhPHOSRecParticles->GetMean(), fhPHOSRecParticles->GetRMS() ) ;
267 printf(" PHOSPhotons Mean : %5.3f , RMS : %5.3f \n", fhPHOSPhotons->GetMean(), fhPHOSPhotons->GetRMS() ) ;
268 printf(" PHOSInvariantMass Mean : %5.3f , RMS : %5.3f \n", fhPHOSInvariantMass->GetMean(), fhPHOSInvariantMass->GetRMS() ) ;
269 printf(" PHOSDigitsEvent Mean : %5.3f , RMS : %5.3f \n", fhPHOSDigitsEvent->GetMean(), fhPHOSDigitsEvent->GetRMS() ) ;
271 TCanvas * cPHOS = new TCanvas("cPHOS", "PHOS ESD Test", 400, 10, 600, 700) ;
275 if ( fhPHOSEnergy->GetMaximum() > 0. )
277 fhPHOSEnergy->SetAxisRange(0, 25.);
278 fhPHOSEnergy->SetLineColor(2);
279 fhPHOSEnergy->Draw();
282 fhPHOSDigits->SetAxisRange(0,25.);
283 fhPHOSDigits->SetLineColor(2);
284 fhPHOSDigits->Draw();
287 if ( fhPHOSRecParticles->GetMaximum() > 0. )
289 fhPHOSRecParticles->SetAxisRange(0, 25.);
290 fhPHOSRecParticles->SetLineColor(2);
291 fhPHOSRecParticles->Draw();
294 if ( fhPHOSPhotons->GetMaximum() > 0. )
296 fhPHOSPhotons->SetAxisRange(0,25.);
297 fhPHOSPhotons->SetLineColor(2);
298 fhPHOSPhotons->Draw();
301 fhPHOSInvariantMass->SetLineColor(2);
302 fhPHOSInvariantMass->Draw();
305 if ( fhPHOSDigitsEvent->GetMaximum() > 0. )
307 fhPHOSDigitsEvent->SetAxisRange(0,40.);
308 fhPHOSDigitsEvent->SetLineColor(2);
309 fhPHOSDigitsEvent->Draw();
311 cPHOS->Print("PHOS.eps");
314 snprintf(line,1024, ".!tar -zcf %s.tar.gz *.eps", GetName()) ;
315 gROOT->ProcessLine(line);
316 snprintf(line,1024, ".!rm -fR *.eps");
317 gROOT->ProcessLine(line);
319 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
321 // char * report = 0x0 ;
323 // sprintf(report,"Problems found, please check!!!");
325 // sprintf(report,"OK");
327 // AliInfo(Form("*** %s Summary Report: %s \n",GetName(), report)) ;