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 //////////////////////////////////////////////////////////////////////////////
37 #include "AliAnalysisTaskPHOSExample.h"
38 #include "AliAnalysisManager.h"
39 #include "AliESDEvent.h"
40 #include "AliESDCaloCluster.h"
41 #include "AliAODEvent.h"
42 #include "AliAODPhoton.h"
44 #include "AliESDVertex.h"
46 //______________________________________________________________________________
47 AliAnalysisTaskPHOSExample::AliAnalysisTaskPHOSExample() :
57 fhPHOSRecParticles(0),
59 fhPHOSInvariantMass(0),
64 //______________________________________________________________________________
65 AliAnalysisTaskPHOSExample::AliAnalysisTaskPHOSExample(const char *name) :
66 AliAnalysisTaskSE(name),
76 fhPHOSRecParticles(0),
78 fhPHOSInvariantMass(0),
84 DefineOutput(1, TList::Class()) ;
87 //____________________________________________________________________________
88 AliAnalysisTaskPHOSExample::AliAnalysisTaskPHOSExample(const AliAnalysisTaskPHOSExample& ap) :
89 AliAnalysisTaskSE(ap.GetName()),
91 fAODPhotons(ap.fAODPhotons),
92 fPhotonsInPhos(ap.fPhotonsInPhos),
93 fPhotonId(ap.fPhotonId),
94 fOutputList(ap.fOutputList),
95 fhPHOSPos(ap.fhPHOSPos),
97 fhPHOSEnergy(ap.fhPHOSEnergy),
98 fhPHOSDigits(ap.fhPHOSDigits),
99 fhPHOSRecParticles(ap.fhPHOSRecParticles),
100 fhPHOSPhotons(ap.fhPHOSPhotons),
101 fhPHOSInvariantMass(ap.fhPHOSInvariantMass),
102 fhPHOSDigitsEvent(ap.fhPHOSDigitsEvent)
107 //_____________________________________________________________________________
108 AliAnalysisTaskPHOSExample& AliAnalysisTaskPHOSExample::operator = (const AliAnalysisTaskPHOSExample& ap)
110 // assignment operator
112 this->~AliAnalysisTaskPHOSExample();
113 new(this) AliAnalysisTaskPHOSExample(ap);
116 fAODPhotons = ap.fAODPhotons;
117 fPhotonsInPhos = ap.fPhotonsInPhos;
118 fPhotonId = ap.fPhotonId;
119 fOutputList = ap.fOutputList;
120 fhPHOSPos = ap.fhPHOSPos;
122 fhPHOSEnergy = ap.fhPHOSEnergy;
123 fhPHOSDigits = ap.fhPHOSDigits;
124 fhPHOSRecParticles = ap.fhPHOSRecParticles;
125 fhPHOSPhotons = ap.fhPHOSPhotons;
126 fhPHOSInvariantMass = ap.fhPHOSInvariantMass;
127 fhPHOSDigitsEvent = ap.fhPHOSDigitsEvent;
132 //______________________________________________________________________________
133 AliAnalysisTaskPHOSExample::~AliAnalysisTaskPHOSExample()
137 fOutputList->Clear() ;
143 //________________________________________________________________________
144 void AliAnalysisTaskPHOSExample::UserCreateOutputObjects()
146 // Create the outputs containers
148 fAODPhotons = new TClonesArray("AliAODPhoton", 0);
149 fAODPhotons->SetName("Photons");
150 AddAODBranch("TClonesArray", &fAODPhotons);
154 fhPHOSPos = new TNtuple("PHOSPos" , "Position in PHOS" , "x:y:z");
155 fhPHOS = new TNtuple("PHOS" , "PHOS" , "event:digits:clusters:photons");
156 fhPHOSEnergy = new TH1D("PHOSEnergy" , "PHOSEnergy" , 100, 0., 100. ) ;
157 fhPHOSDigits = new TH1I("PHOSDigitsCluster" , "PHOSDigits" , 20 , 0 , 20 ) ;
158 fhPHOSRecParticles = new TH1D("PHOSRecParticles" , "PHOSRecParticles" , 20 , 0., 20. ) ;
159 fhPHOSPhotons = new TH1I("PHOSPhotons" , "PHOSPhotons" , 20 , 0 , 20 ) ;
160 fhPHOSInvariantMass = new TH1D("PHOSInvariantMass" , "PHOSInvariantMass" , 400, 0., 400.) ;
161 fhPHOSDigitsEvent = new TH1I("PHOSDigitsEvent" , "PHOSDigitsEvent" , 30 , 0 , 30 ) ;
163 // create output container
165 fOutputList = new TList() ;
166 fOutputList->SetName(GetName()) ;
168 fOutputList->AddAt(fhPHOSPos, 0) ;
169 fOutputList->AddAt(fhPHOS, 1) ;
170 fOutputList->AddAt(fhPHOSEnergy, 2) ;
171 fOutputList->AddAt(fhPHOSDigits, 3) ;
172 fOutputList->AddAt(fhPHOSRecParticles, 4) ;
173 fOutputList->AddAt(fhPHOSPhotons, 5) ;
174 fOutputList->AddAt(fhPHOSInvariantMass, 6) ;
175 fOutputList->AddAt(fhPHOSDigitsEvent, 7) ;
180 //______________________________________________________________________________
181 void AliAnalysisTaskPHOSExample::UserExec(Option_t *)
183 // Processing of one event
185 // if ( !((Entry()-1)%100) )
186 AliInfo(Form(" Processing event # %lld", Entry())) ;
187 AliESDEvent* esd = (AliESDEvent*)InputEvent();
189 //************************ PHOS *************************************
190 TRefArray * caloClustersArr = new TRefArray();
191 esd->GetPHOSClusters(caloClustersArr);
193 const Int_t kNumberOfPhosClusters = caloClustersArr->GetEntries() ;
195 TVector3 ** phosVector = new TVector3*[kNumberOfPhosClusters] ;
196 Float_t * phosPhotonsEnergy = new Float_t[kNumberOfPhosClusters] ;
198 Int_t numberOfDigitsInPhos = 0 ;
199 Double_t v[3] ; //vertex ;
200 esd->GetVertex()->GetXYZ(v) ;
203 // loop over the PHOS Cluster
204 for(phosCluster = 0 ; phosCluster < kNumberOfPhosClusters ; phosCluster++) {
205 AliESDCaloCluster * caloCluster = (AliESDCaloCluster *) caloClustersArr->At(phosCluster) ;
207 //AliESDCaloCluster * caloCluster = fESD->GetCaloCluster(phosCluster) ;
210 caloCluster->GetPosition( pos ) ;
211 fhPHOSEnergy->Fill( caloCluster->E() ) ;
212 fhPHOSPos->Fill( pos[0], pos[1], pos[2] ) ;
213 fhPHOSDigits->Fill(Entry(), caloCluster->GetNCells() ) ;
214 numberOfDigitsInPhos += caloCluster->GetNCells() ;
215 Double_t * pid = caloCluster->GetPid() ;
216 if(pid[AliPID::kPhoton] > GetPhotonId() ) {
217 phosVector[fPhotonsInPhos] = new TVector3(pos[0],pos[1],pos[2]) ;
218 phosPhotonsEnergy[fPhotonsInPhos]=caloCluster->E() ;
219 TLorentzVector momentum ;
220 caloCluster->GetMomentum(momentum, v);
221 new ((*fAODPhotons)[fPhotonsInPhos++]) AliAODPhoton (momentum);
226 fhPHOSRecParticles->Fill(kNumberOfPhosClusters);
227 fhPHOSPhotons->Fill(fPhotonsInPhos);
228 fhPHOSDigitsEvent->Fill(numberOfDigitsInPhos);
229 fhPHOS->Fill(Entry(), numberOfDigitsInPhos, kNumberOfPhosClusters, fPhotonsInPhos) ;
232 if (fPhotonsInPhos > 1 ) {
233 Int_t phosPhoton1, phosPhoton2 ;
234 for(phosPhoton1 = 0 ; phosPhoton1 < fPhotonsInPhos ; phosPhoton1++) {
235 for(phosPhoton2 = phosPhoton1 + 1 ; phosPhoton2 < fPhotonsInPhos ; phosPhoton2++) {
236 Float_t tempMass = TMath::Sqrt( 2 * phosPhotonsEnergy[phosPhoton1] * phosPhotonsEnergy[phosPhoton2] *
237 ( 1 - TMath::Cos(phosVector[phosPhoton1]->Angle(*phosVector[phosPhoton2])) )
239 fhPHOSInvariantMass->Fill(tempMass*1000.);
244 PostData(1, fOutputList);
246 delete [] phosVector ;
247 delete [] phosPhotonsEnergy ;
252 //______________________________________________________________________________
253 void AliAnalysisTaskPHOSExample::Init()
255 // Intialisation of parameters
256 AliInfo("Doing initialisation") ;
260 //______________________________________________________________________________
261 void AliAnalysisTaskPHOSExample::Terminate(Option_t *)
263 // Processing when the event loop is ended
265 Bool_t problem = kFALSE ;
266 AliInfo(Form(" *** %s Report:", GetName())) ;
267 printf(" PHOSEnergy Mean : %5.3f , RMS : %5.3f \n", fhPHOSEnergy->GetMean(), fhPHOSEnergy->GetRMS() ) ;
268 printf(" PHOSDigits Mean : %5.3f , RMS : %5.3f \n", fhPHOSDigits->GetMean(), fhPHOSDigits->GetRMS() ) ;
269 printf(" PHOSRecParticles Mean : %5.3f , RMS : %5.3f \n", fhPHOSRecParticles->GetMean(), fhPHOSRecParticles->GetRMS() ) ;
270 printf(" PHOSPhotons Mean : %5.3f , RMS : %5.3f \n", fhPHOSPhotons->GetMean(), fhPHOSPhotons->GetRMS() ) ;
271 printf(" PHOSInvariantMass Mean : %5.3f , RMS : %5.3f \n", fhPHOSInvariantMass->GetMean(), fhPHOSInvariantMass->GetRMS() ) ;
272 printf(" PHOSDigitsEvent Mean : %5.3f , RMS : %5.3f \n", fhPHOSDigitsEvent->GetMean(), fhPHOSDigitsEvent->GetRMS() ) ;
274 TCanvas * cPHOS = new TCanvas("cPHOS", "PHOS ESD Test", 400, 10, 600, 700) ;
278 if ( fhPHOSEnergy->GetMaximum() > 0. )
280 fhPHOSEnergy->SetAxisRange(0, 25.);
281 fhPHOSEnergy->SetLineColor(2);
282 fhPHOSEnergy->Draw();
285 fhPHOSDigits->SetAxisRange(0,25.);
286 fhPHOSDigits->SetLineColor(2);
287 fhPHOSDigits->Draw();
290 if ( fhPHOSRecParticles->GetMaximum() > 0. )
292 fhPHOSRecParticles->SetAxisRange(0, 25.);
293 fhPHOSRecParticles->SetLineColor(2);
294 fhPHOSRecParticles->Draw();
297 if ( fhPHOSPhotons->GetMaximum() > 0. )
299 fhPHOSPhotons->SetAxisRange(0,25.);
300 fhPHOSPhotons->SetLineColor(2);
301 fhPHOSPhotons->Draw();
304 fhPHOSInvariantMass->SetLineColor(2);
305 fhPHOSInvariantMass->Draw();
308 if ( fhPHOSDigitsEvent->GetMaximum() > 0. )
310 fhPHOSDigitsEvent->SetAxisRange(0,40.);
311 fhPHOSDigitsEvent->SetLineColor(2);
312 fhPHOSDigitsEvent->Draw();
314 cPHOS->Print("PHOS.eps");
317 sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ;
318 gROOT->ProcessLine(line);
319 sprintf(line, ".!rm -fR *.eps");
320 gROOT->ProcessLine(line);
322 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
326 report="Problems found, please check!!!";
330 AliInfo(Form("*** %s Summary Report: %s \n",GetName(), report)) ;