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 // Algorythm class to analyze PHOS events. In this class we demostrate,
20 // how to handle reconstructed objects with AliPHSOIndexToObject.
21 // As an example we propose sulotions for four most frequently used tasks:
22 // DrawRecon(...) - to draw reconstructed objects in the PHOS plane,
23 // very usefull in the debuging
24 // InvarianMass(...) - to calculate "REAL" and "MIXED" photon pairs
25 // invariant mass distributions
26 // EnergyResoluition(...) -\ Energy and position resolutions of the
27 // PositionResolution(...)-/ reconstructed photons
28 // Contamination(...) - calculates contamination of the photon spectrum and
29 // pobability of reconstruction of several primaries as
30 // kGAMMA,kELECTRON etc.
32 // root [0] AliPHOSAnalyze * a = new AliPHOSAnalyze("galice.root")
33 // // set the file you want to analyse
34 // root [1] a->DrawRecon(1,3)
35 // // plot RecObjects, made in event 1, PHOS module 3
36 // root [2] a->DrawRecon(1,3,"PHOSRP","another PID")
37 // // plot RecObjets made in the event 1, PHOS module 3,
38 // // produced in the another reconstruction pass,
39 // // which produced PHOS RecParticles ("PHOSRP") with
40 // // title "another PID".
41 // root [3] a->InvariantMass()
42 // // Calculates "REAL" and "MIXED" invariant mass
43 // // distributions of kGAMMA and (kGAMMA+kNEUTRALEM)
44 // // and APPENDS this to the file "invmass.root"
45 // root [4] a->PositionResolution()
46 // // calculates two dimentional histos: energy of the primary
47 // // photon vs distance betwin incedence point and reconstructed
48 // // poisition. One can analyse the produced file position.root
49 // // with macro PhotonPosition.C
50 // root [5] a->EnergyResolution()
51 // // calculates two dimentional histos: energy of the primary
52 // // photon vs energy of the reconstructed particle. One can
53 // // analyse the produced file energy.root
54 // // with macro PhotonEnergy.C
55 // root [6] a->Contamination()
56 // // fills spectra of primary photons and several kinds of
57 // // reconstructed particles, so that analyzing them one can
58 // // estimate conatmination, efficiency of registration etc.
60 //*-- Author: Dmitri Peressounko (SUBATECH & RRC Kurchatov Institute)
61 //////////////////////////////////////////////////////////////////////////////
64 // --- ROOT system ---
70 #include "TParticle.h"
71 #include "TClonesArray.h"
75 // --- Standard library ---
77 // --- AliRoot header files ---
81 #include "AliPHOSGeometry.h"
82 #include "AliPHOSAnalyze.h"
83 #include "AliPHOSDigit.h"
84 #include "AliPHOSSDigitizer.h"
85 #include "AliPHOSEmcRecPoint.h"
86 #include "AliPHOSCpvRecPoint.h"
87 #include "AliPHOSTrackSegment.h"
88 #include "AliPHOSRecParticle.h"
89 #include "AliPHOSLoader.h"
92 ClassImp(AliPHOSAnalyze)
94 //____________________________________________________________________________
95 AliPHOSAnalyze::AliPHOSAnalyze():
96 fCorrection(1.2), //Value calculated for default parameters of reconstruction
101 // default ctor (useless)
104 //____________________________________________________________________________
105 AliPHOSAnalyze::AliPHOSAnalyze(Text_t * fileName):
106 fCorrection(1.05), //Value calculated for default parameters of reconstruction
111 // ctor: analyze events from root file "name"
112 fRunLoader = AliRunLoader::Open(fileName,"AliPHOSAnalyze");
113 if (fRunLoader == 0x0)
115 AliError(Form("Error Loading session"));
119 //____________________________________________________________________________
120 AliPHOSAnalyze::AliPHOSAnalyze(const AliPHOSAnalyze & ana):
128 ( (AliPHOSAnalyze &)ana ).Copy(*this) ;
131 //____________________________________________________________________________
132 AliPHOSAnalyze::~AliPHOSAnalyze()
137 //____________________________________________________________________________
138 void AliPHOSAnalyze::DrawRecon(Int_t Nevent,Int_t Nmod){
139 //Draws pimary particles and reconstructed
140 //digits, RecPoints, RecPartices etc
141 //for event Nevent in the module Nmod.
143 //========== Create ObjectLoader
144 if (fRunLoader == 0x0)
146 AliError(Form("Error Loading session"));
150 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
153 AliError(Form("Could not obtain the Loader object !"));
158 if(Nevent >= fRunLoader->GetNumberOfEvents() ) {
159 AliError(Form("There is no event %d only %d events available", Nevent, fRunLoader->GetNumberOfEvents() )) ;
162 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
163 fRunLoader->GetEvent(Nevent);
165 Int_t nx = phosgeom->GetNPhi() ;
166 Int_t nz = phosgeom->GetNZ() ;
167 Float_t * cri= phosgeom->GetEMCAGeometry()->GetCrystalHalfSize() ;
168 Float_t x = nx*cri[0] ;
169 Float_t z = nz*cri[2] ;
170 Int_t nxCPV = (Int_t) (nx*phosgeom->GetPadSizePhi()/(2.*cri[0])) ;
171 Int_t nzCPV = (Int_t) (nz*phosgeom->GetPadSizeZ()/(2.*cri[2])) ;
173 TH2F * emcDigits = (TH2F*) gROOT->FindObject("emcDigits") ;
175 emcDigits->Delete() ;
176 emcDigits = new TH2F("emcDigits","EMC digits", nx,-x,x,nz,-z,z);
177 TH2F * emcSdigits =(TH2F*) gROOT->FindObject("emcSdigits") ;
179 emcSdigits->Delete() ;
180 emcSdigits = new TH2F("emcSdigits","EMC sdigits", nx,-x,x,nz,-z,z);
181 TH2F * emcRecPoints = (TH2F*)gROOT->FindObject("emcRecPoints") ;
183 emcRecPoints->Delete() ;
184 emcRecPoints = new TH2F("emcRecPoints","EMC RecPoints",nx,-x,x,nz,-z,z);
185 TH2F * cpvSdigits =(TH2F*) gROOT->FindObject("cpvSdigits") ;
187 cpvSdigits->Delete() ;
188 cpvSdigits = new TH2F("cpvSdigits","CPV sdigits", nx,-x,x,nz,-z,z);
189 TH2F * cpvDigits = (TH2F*)gROOT->FindObject("cpvDigits") ;
191 cpvDigits->Delete() ;
192 cpvDigits = new TH2F("cpvDigits","CPV digits", nxCPV,-x,x,nzCPV,-z,z) ;
193 TH2F * cpvRecPoints= (TH2F*)gROOT->FindObject("cpvRecPoints") ;
195 cpvRecPoints->Delete() ;
196 cpvRecPoints = new TH2F("cpvRecPoints","CPV RecPoints", nxCPV,-x,x,nzCPV,-z,z) ;
198 TH2F * phot = (TH2F*)gROOT->FindObject("phot") ;
201 phot = new TH2F("phot","Primary Photon", nx,-x,x,nz,-z,z);
202 TH2F * recPhot = (TH2F*)gROOT->FindObject("recPhot") ;
205 recPhot = new TH2F("recPhot","RecParticles with primary Photon",nx,-x,x,nz,-z,z);
208 Double_t vtx[3]={0.,0.,0.} ;
209 //DP: extract vertex either from Generator or from data
212 //Plot Primary Particles
214 if (fRunLoader->Stack() == 0x0) fRunLoader->LoadKinematics("READ");
217 const TParticle * primary ;
219 for ( iPrimary = 0 ; iPrimary < fRunLoader->Stack()->GetNprimary() ; iPrimary++)
221 primary = fRunLoader->Stack()->Particle(iPrimary);
223 Int_t primaryType = primary->GetPdgCode();
224 // if( (primaryType == 211)||(primaryType == -211)||(primaryType == 2212)||(primaryType == -2212)
225 // ||(primaryType == 11)||(primaryType == -11) ) {
226 // Int_t moduleNumber ;
227 // Double_t primX, primZ ;
228 // phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
229 // if(moduleNumber==Nmod)
230 // charg->Fill(primZ,primX,primary->Energy()) ;
232 if( primaryType == 22 ) {
234 Double_t primX, primZ ;
235 phosgeom->ImpactOnEmc(vtx,primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
236 if(moduleNumber==Nmod)
237 phot->Fill(primZ,primX,primary->Energy()) ;
240 // if( primaryType == -2112 ) {
241 // Int_t moduleNumber ;
242 // Double_t primX, primZ ;
243 // phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
244 // if(moduleNumber==Nmod)
245 // nbar->Fill(primZ,primX,primary->Energy()) ;
252 AliPHOSDigit * sdigit ;
253 const TClonesArray * sdigits = gime->SDigits() ;
254 Int_t nsdig[5] = {0,0,0,0,0} ;
256 for(iSDigit = 0; iSDigit < sdigits->GetEntriesFast() ; iSDigit++)
258 sdigit = (AliPHOSDigit *) sdigits->At(iSDigit) ;
260 phosgeom->AbsToRelNumbering(sdigit->GetId(), relid) ;
262 phosgeom->RelPosInModule(relid,xd,zd);
263 Float_t e = sdigit->GetEnergy() ;
264 nsdig[relid[0]-1]++ ;
266 if(relid[1]==0) //EMC
267 emcSdigits->Fill(xd,zd,e) ;
269 cpvSdigits->Fill(xd,zd,e) ;
274 message = "Number of EMC + CPV SDigits per module: \n" ;
275 message += "%d %d %d %d %d\n";
276 AliInfo(Form(message.Data(), nsdig[0], nsdig[1], nsdig[2], nsdig[3], nsdig[4] )) ;
280 AliPHOSDigit * digit ;
281 const TClonesArray * digits = gime->Digits();
283 for(iDigit = 0; iDigit < digits->GetEntriesFast(); iDigit++)
285 digit = (AliPHOSDigit *) digits->At(iDigit) ;
287 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
289 phosgeom->RelPosInModule(relid,xd,zd) ;
290 Float_t e = digit->GetEnergy() ;
292 if(relid[1]==0) //EMC
293 emcDigits->Fill(xd,zd,e) ;
295 cpvDigits->Fill(xd,zd,e) ;
304 TObjArray * emcrp = gime->EmcRecPoints() ;
306 for(irecp = 0; irecp < emcrp->GetEntriesFast() ; irecp ++){
307 AliPHOSEmcRecPoint * emc = (AliPHOSEmcRecPoint *) emcrp->At(irecp) ;
308 if(emc->GetPHOSMod()==Nmod){
309 emc->GetLocalPosition(pos) ;
310 emcRecPoints->Fill(pos.X(),pos.Z(),emc->GetEnergy());
315 TObjArray * cpvrp = gime->CpvRecPoints() ;
317 for(irecp = 0; irecp < cpvrp->GetEntriesFast() ; irecp ++){
318 AliPHOSRecPoint * cpv = (AliPHOSCpvRecPoint *) cpvrp->At(irecp) ;
319 if(cpv->GetPHOSMod()==Nmod){
320 cpv->GetLocalPosition(pos) ;
321 cpvRecPoints->Fill(pos.X(),pos.Z(),cpv->GetEnergy());
327 AliPHOSRecParticle * recParticle ;
329 TClonesArray * rp = gime->RecParticles() ;
330 TClonesArray * ts = gime->TrackSegments() ;
331 if(rp && ts && emcrp) {
332 for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast() ; iRecParticle++ )
334 recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
335 Int_t moduleNumberRec ;
336 Double_t recX, recZ ;
337 phosgeom->ImpactOnEmc(vtx,recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
338 if(moduleNumberRec == Nmod){
340 Double_t minDistance = 5. ;
341 Int_t closestPrimary = -1 ;
343 //extract list of primaries: it is stored at EMC RecPoints
344 Int_t emcIndex = ((AliPHOSTrackSegment *) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
345 Int_t numberofprimaries ;
346 Int_t * listofprimaries = ((AliPHOSRecPoint*) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries) ;
348 const TParticle * primPart ;
349 Double_t distance = minDistance ;
351 for ( index = 0 ; index < numberofprimaries ; index++){
352 primPart = fRunLoader->Stack()->Particle(listofprimaries[index]) ;
354 Double_t primX, primZ ;
355 phosgeom->ImpactOnEmc(vtx,primPart->Theta(), primPart->Phi(), moduleNumber, primX, primZ) ;
356 if(moduleNumberRec == moduleNumber)
357 distance = TMath::Sqrt((recX-primX)*(recX-primX)+(recZ-primZ)*(recZ-primZ) ) ;
358 if(minDistance > distance)
360 minDistance = distance ;
361 closestPrimary = listofprimaries[index] ;
365 if(closestPrimary >=0 ){
367 Int_t primaryType = fRunLoader->Stack()->Particle(closestPrimary)->GetPdgCode() ;
370 recPhot->Fill(recZ,recX,recParticle->Energy()) ;
372 // if(primaryType==-2112)
373 // recNbar->Fill(recZ,recX,recParticle->Energy()) ;
380 //Plot made histograms
381 emcSdigits->Draw("box") ;
382 emcDigits->SetLineColor(5) ;
383 emcDigits->Draw("boxsame") ;
384 emcRecPoints->SetLineColor(2) ;
385 emcRecPoints->Draw("boxsame") ;
386 cpvSdigits->SetLineColor(1) ;
387 cpvSdigits->Draw("boxsame") ;
390 //____________________________________________________________________________
391 void AliPHOSAnalyze::Ls(){
392 //lists branches and titles of PHOS-related branches of TreeR, TreeD, TreeS
394 if (fRunLoader == 0x0)
396 AliError(Form("Error Loading session"));
400 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
403 AliError(Form("Could not obtain the Loader object !"));
409 TObjArray * branches;
411 if (gime->TreeS() == 0x0)
413 if (gime->LoadSDigits("READ"))
415 AliError(Form("Problems with loading summable digits"));
419 branches = gime->TreeS()->GetListOfBranches() ;
422 message = "TreeS:\n" ;
423 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
424 TBranch * branch=(TBranch *) branches->At(ibranch) ;
425 if(strstr(branch->GetName(),"PHOS") ){
427 message += branch->GetName() ;
429 message += branch->GetTitle() ;
433 if (gime->TreeD() == 0x0)
435 if (gime->LoadDigits("READ"))
437 AliError(Form("Problems with loading digits"));
442 branches = gime->TreeD()->GetListOfBranches() ;
444 message += "TreeD:\n" ;
445 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
446 TBranch * branch=(TBranch *) branches->At(ibranch) ;
447 if(strstr(branch->GetName(),"PHOS") ) {
449 message += branch->GetName() ;
451 message += branch->GetTitle() ;
456 if (gime->TreeR() == 0x0)
458 if (gime->LoadRecPoints("READ"))
460 AliError(Form("Problems with loading rec points"));
465 branches = gime->TreeR()->GetListOfBranches() ;
467 message += "TreeR: \n" ;
468 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
469 TBranch * branch=(TBranch *) branches->At(ibranch) ;
470 if(strstr(branch->GetName(),"PHOS") ) {
472 message += branch->GetName() ;
474 message += branch->GetTitle() ;
478 AliInfo(Form(message.Data())) ;
480 //____________________________________________________________________________
481 void AliPHOSAnalyze::InvariantMass()
483 // Calculates Real and Mixed invariant mass distributions
484 if (fRunLoader == 0x0)
486 AliError(Form("Error Loading session"));
490 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
493 AliError(Form("Could not obtain the Loader object !"));
497 gime->LoadRecParticles("READ");
499 Int_t nMixedEvents = 4 ; //# of events used for calculation of 'mixed' distribution
503 TFile * mfile = new TFile("invmass.root","update");
505 //========== Reading /Booking Histograms
507 hRealEM = (TH2D*) mfile->Get("hRealEM") ;
509 hRealEM = new TH2D("hRealEM", "Real for EM particles", 250,0.,1.,40,0.,4.) ;
510 TH2D * hRealPhot = 0 ;
512 hRealPhot = (TH2D*)mfile->Get("hRealPhot");
514 hRealPhot = new TH2D("hRealPhot", "Real for kPhoton particles", 250,0.,1.,40,0.,4.) ;
516 TH2D * hMixedEM = 0 ;
517 hMixedEM = (TH2D*) mfile->Get("hMixedEM") ;
519 hMixedEM = new TH2D("hMixedEM", "Mixed for EM particles", 250,0.,1.,40,0.,4.) ;
521 TH2D * hMixedPhot = 0 ;
522 hMixedPhot = (TH2D*) mfile->Get("hMixedPhot") ;
524 hMixedPhot = new TH2D("hMixedPhot","Mixed for kPhoton particles",250,0.,1.,40,0.,4.) ;
527 //reading event and copyng it to TConesArray of all photons
529 TClonesArray * allRecParticleList = new TClonesArray("AliPHOSRecParticle", 1000) ;
530 Int_t * nRecParticles = new Int_t[nMixedEvents] ; // to mark boundaries of each event in the total list
531 for(Int_t index = 0; index < nMixedEvents; index ++)
532 nRecParticles[index] = 0 ;
533 Int_t iRecPhot = 0 ; // number of EM particles in total list
535 //scan over all events
539 if (fRunLoader->TreeE() == 0x0) fRunLoader->LoadHeader();
542 Int_t maxevent = (Int_t)fRunLoader->TreeE()->GetEntries();
543 // for(event = 0; event < gime->MaxEvent(); event++ ){
547 for(event = 0; event < maxevent; event++ ){
548 fRunLoader->GetEvent(event); //will read only TreeR
550 //copy EM RecParticles to the "total" list
551 const AliPHOSRecParticle * recParticle ;
553 TClonesArray * rp = gime->RecParticles() ;
555 AliError(Form("Can't find RecParticles")) ;
559 for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast(); iRecParticle++ )
561 recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
562 if((recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)||
563 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW))
564 new( (*allRecParticleList)[iRecPhot++] ) AliPHOSRecParticle(*recParticle) ;
567 Int_t mevent = event%nMixedEvents ; //event number in the "mixed" cicle
568 nRecParticles[mevent] = iRecPhot-1 ;
570 //check, if it is time to calculate invariant mass?
571 if((mevent == 0) && (event +1 == maxevent)){
573 // if((mevent == 0) && (event +1 == gime->MaxEvent())){
575 //calculate invariant mass:
577 Int_t nCurEvent = 0 ;
579 for(irp1 = 0; irp1 < allRecParticleList->GetEntries()-1; irp1++){
580 AliPHOSRecParticle * rp1 = (AliPHOSRecParticle *)allRecParticleList->At(irp1) ;
582 for(irp2 = irp1+1; irp2 < allRecParticleList->GetEntries(); irp2++){
583 AliPHOSRecParticle * rp2 = (AliPHOSRecParticle *)allRecParticleList->At(irp2) ;
586 invMass = (rp1->Energy()+rp2->Energy())*(rp1->Energy()+rp2->Energy())-
587 (rp1->Px()+rp2->Px())*(rp1->Px()+rp2->Px())-
588 (rp1->Py()+rp2->Py())*(rp1->Py()+rp2->Py())-
589 (rp1->Pz()+rp2->Pz())*(rp1->Pz()+rp2->Pz()) ;
592 invMass = TMath::Sqrt(invMass);
595 pt = TMath::Sqrt((rp1->Px()+rp2->Px() )*( rp1->Px()+rp2->Px() ) +
596 (rp1->Py()+rp2->Py() )*( rp1->Py()+rp2->Py() ) );
598 if(irp1 > nRecParticles[nCurEvent])
601 if(irp2 <= nRecParticles[nCurEvent]){ //'Real' event
602 hRealEM->Fill(invMass,pt);
603 if((rp1->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)&&
604 (rp2->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST) )
605 hRealPhot->Fill(invMass,pt);
608 hMixedEM->Fill(invMass,pt);
609 if((rp1->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)&&
610 (rp2->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST) )
611 hMixedPhot->Fill(invMass,pt);
614 } //loop over second rp
615 }//loop over first rp
617 //Make some cleanings
618 for(Int_t index = 0; index < nMixedEvents; index ++)
619 nRecParticles[index] = 0 ;
621 allRecParticleList->Clear() ;
625 delete allRecParticleList ;
630 hRealEM->Write(0,kOverwrite) ;
631 hRealPhot->Write(0,kOverwrite) ;
632 hMixedEM->Write(0,kOverwrite) ;
633 hMixedPhot->Write(0,kOverwrite) ;
638 delete nRecParticles;
642 //____________________________________________________________________________
643 void AliPHOSAnalyze::EnergyResolution()
645 //fills two dimentional histo: energy of primary vs. energy of reconstructed
647 TH2F * hAllEnergy = 0 ; //all reconstructed with primary photon
648 TH2F * hPhotEnergy= 0 ; //kGamma with primary photon
649 TH2F * hEMEnergy = 0 ; //electromagnetic with primary photon
651 //opening file and reading histograms if any
652 TFile * efile = new TFile("energy.root","update");
654 hAllEnergy = (TH2F*)efile->Get("hAllEnergy") ;
656 hAllEnergy = new TH2F("hAllEnergy", "Energy of any RP with primary photon",100, 0., 5., 100, 0., 5.);
658 hPhotEnergy =(TH2F*) efile->Get("hPhotEnergy") ;
660 hPhotEnergy = new TH2F("hPhotEnergy", "Energy of kGAMMA with primary photon",100, 0., 5., 100, 0., 5.);
662 hEMEnergy =(TH2F*) efile->Get("hEMEnergy");
664 hEMEnergy = new TH2F("hEMEnergy", "Energy of EM with primary photon", 100, 0., 5., 100, 0., 5.);
667 if (fRunLoader == 0x0)
669 AliError(Form("Error Loading session"));
673 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
676 AliError(Form("Could not obtain the Loader object !"));
681 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
684 Int_t maxevent = (Int_t)fRunLoader->TreeE()->GetEntries();
686 fRunLoader->LoadKinematics("READ");
687 gime->LoadTracks("READ");
689 for ( ievent=0; ievent < maxevent ; ievent++){
691 //read the current event
692 fRunLoader->GetEvent(ievent) ;
694 Double_t vtx[3]={0.,0.,0.} ;
696 const AliPHOSRecParticle * recParticle ;
698 TClonesArray * rp = gime->RecParticles() ;
700 AliError(Form("Event %d, Can't find RecParticles ", ievent)) ;
703 TClonesArray * ts = gime->TrackSegments() ;
705 AliError(Form("Event %d, Can't find TrackSegments", ievent)) ;
708 TObjArray * emcrp = gime->EmcRecPoints() ;
710 AliError(Form("Event %d, Can't find EmcRecPoints")) ;
714 for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast() ;iRecParticle++ ){
715 recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
717 //find the closest primary
718 Int_t moduleNumberRec ;
719 Double_t recX, recZ ;
720 phosgeom->ImpactOnEmc(vtx,recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
722 Double_t minDistance = 100. ;
723 Int_t closestPrimary = -1 ;
725 //extract list of primaries: it is stored at EMC RecPoints
726 Int_t emcIndex = ((AliPHOSTrackSegment*) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
727 Int_t numberofprimaries ;
728 Int_t * listofprimaries = ((AliPHOSEmcRecPoint*) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries) ;
731 const TParticle * primary ;
732 Double_t distance = minDistance ;
735 Double_t dZmin = 0. ;
736 for ( index = 0 ; index < numberofprimaries ; index++){
738 primary = fRunLoader->Stack()->Particle(listofprimaries[index]) ;
741 Double_t primX, primZ ;
742 phosgeom->ImpactOnEmc(vtx,primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
743 if(moduleNumberRec == moduleNumber) {
746 distance = TMath::Sqrt(dX*dX + dZ*dZ) ;
747 if(minDistance > distance) {
748 minDistance = distance ;
751 closestPrimary = listofprimaries[index] ;
756 //if found primary, fill histograms
757 if(closestPrimary >=0 ){
758 primary = fRunLoader->Stack()->Particle(closestPrimary) ;
759 if(primary->GetPdgCode() == 22){
760 hAllEnergy->Fill(primary->Energy(), recParticle->Energy()) ;
761 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST){
762 hPhotEnergy->Fill(primary->Energy(), recParticle->Energy() ) ;
763 hEMEnergy->Fill(primary->Energy(), recParticle->Energy() ) ;
766 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW)
767 hEMEnergy->Fill(primary->Energy(), recParticle->Energy() ) ;
773 //write filled histograms
775 hAllEnergy->Write(0,kOverwrite) ;
776 hPhotEnergy->Write(0,kOverwrite) ;
777 hEMEnergy->Write(0,kOverwrite) ;
783 //____________________________________________________________________________
784 void AliPHOSAnalyze::PositionResolution()
786 //fills two dimentional histo: energy vs. primary - reconstructed distance
790 TH2F * hAllPosition = 0; // Position of any RP with primary photon
791 TH2F * hPhotPosition = 0; // Position of kGAMMA with primary photon
792 TH2F * hEMPosition = 0; // Position of EM with primary photon
794 TH1F * hAllPositionX = 0; // X-Position Resolution of photons with photon primary
795 TH1F * hAllPositionZ = 0; // Z-Position Resolution of photons with photon primary
798 //opening file and reading histograms if any
799 TFile * pfile = new TFile("position.root","update");
801 hAllPosition = (TH2F*)pfile->Get("hAllPosition");
802 if(hAllPosition == 0)
803 hAllPosition = new TH2F("hAllPosition",
804 "Position of any RP with primary photon",100, 0., 5., 100, 0., 5.);
805 hPhotPosition= (TH2F*)pfile->Get("hPhotPosition");
806 if(hPhotPosition == 0)
807 hPhotPosition = new TH2F("hPhotPosition",
808 "Position of kGAMMA with primary photon",100, 0., 5., 100, 0., 5.);
809 hEMPosition= (TH2F*)pfile->Get("hEMPosition") ;
811 hEMPosition = new TH2F("hEMPosition",
812 "Position of EM with primary photon", 100, 0., 5., 100, 0., 5.);
813 hAllPositionX = (TH1F*)pfile->Get("hAllPositionX") ;
814 if(hAllPositionX == 0)
815 hAllPositionX = new TH1F("hAllPositionX",
816 "Delta X of any RP with primary photon",100, -2., 2.);
817 hAllPositionZ =(TH1F*) pfile->Get("hAllPositionZ") ;
818 if(hAllPositionZ == 0)
819 hAllPositionZ = new TH1F("hAllPositionZ",
820 "Delta X of any RP with primary photon",100, -2., 2.);
822 if (fRunLoader == 0x0)
824 AliError(Form("Error Loading session"));
828 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
831 AliError(Form("Could not obtain the Loader object !"));
835 if (fRunLoader->TreeE() == 0x0) fRunLoader->LoadHeader();
837 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
840 Int_t maxevent = (Int_t)fRunLoader->TreeE()->GetEntries() ;
841 for ( ievent=0; ievent < maxevent ; ievent++){
843 //read the current event
844 fRunLoader->GetEvent(ievent) ;
846 //DP:Extract vertex position
847 Double_t vtx[3]={0.,0.,0.} ;
849 TClonesArray * rp = gime->RecParticles() ;
851 AliError(Form("Event %d, Can't find RecParticles", ievent)) ;
854 TClonesArray * ts = gime->TrackSegments() ;
856 AliError(Form("Event %d, Can't find TrackSegments", ievent)) ;
859 TObjArray * emcrp = gime->EmcRecPoints() ;
861 AliError(Form("Event %d, Can't find EmcRecPoints", ievent)) ;
866 const AliPHOSRecParticle * recParticle ;
868 for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast(); iRecParticle++ ){
869 recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
871 //find the closest primary
872 Int_t moduleNumberRec ;
873 Double_t recX, recZ ;
874 phosgeom->ImpactOnEmc(vtx,recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
876 Double_t minDistance = 100. ;
877 Int_t closestPrimary = -1 ;
879 //extract list of primaries: it is stored at EMC RecPoints
880 Int_t emcIndex = ((AliPHOSTrackSegment*) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
881 Int_t numberofprimaries ;
882 Int_t * listofprimaries = ((AliPHOSEmcRecPoint *) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries) ;
885 const TParticle * primary ;
886 Double_t distance = minDistance ;
887 Double_t dX = 1000; // incredible number
888 Double_t dZ = 1000; // for the case if no primary will be found
890 Double_t dZmin = 0. ;
891 for ( index = 0 ; index < numberofprimaries ; index++){
892 primary = fRunLoader->Stack()->Particle(listofprimaries[index]) ;
894 Double_t primX, primZ ;
895 phosgeom->ImpactOnEmc(vtx,primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
896 if(moduleNumberRec == moduleNumber) {
899 distance = TMath::Sqrt(dX*dX + dZ*dZ) ;
900 if(minDistance > distance) {
901 minDistance = distance ;
904 closestPrimary = listofprimaries[index] ;
909 //if found primary, fill histograms
910 if(closestPrimary >=0 ){
911 primary = fRunLoader->Stack()->Particle(closestPrimary) ;
912 if(primary->GetPdgCode() == 22){
913 hAllPosition->Fill(primary->Energy(), minDistance) ;
914 hAllPositionX->Fill(primary->Energy(), dX) ;
915 hAllPositionZ->Fill(primary->Energy(), dZ) ;
916 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST){
917 hPhotPosition->Fill(primary->Energy(), minDistance ) ;
918 hEMPosition->Fill(primary->Energy(), minDistance ) ;
921 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW)
922 hEMPosition->Fill(primary->Energy(), minDistance ) ;
928 //Write output histgrams
930 hAllPosition->Write(0,kOverwrite) ;
931 hAllPositionX->Write(0,kOverwrite) ;
932 hAllPositionZ->Write(0,kOverwrite) ;
933 hPhotPosition->Write(0,kOverwrite) ;
934 hEMPosition->Write(0,kOverwrite) ;
941 //____________________________________________________________________________
942 void AliPHOSAnalyze::Contamination(){
943 // fills spectra of primary photons and several kinds of
944 // reconstructed particles, so that analyzing them one can
945 // estimate conatmination, efficiency of registration etc.
947 //define several general histograms
948 TH1F * hPrimary = 0; //spectrum (P_t distribution) of primary photons
949 TH1F * hAllRP = 0; //spectrum of all RecParticles in PHOS
950 TH1F * hPhot = 0; //spectrum of kGAMMA RecParticles
951 TH1F * hShape = 0; //spectrum of all EM RecParticles
952 TH1F * hVeto = 0; //spectrum of all neutral RecParticles
954 //Now separate histograms in accoradance with primary
956 TH1F * hPhotReg = 0; //Registeres as photon
957 TH1F * hPhotEM = 0; //Registered as EM
960 TH1F * hNReg = 0; //Registeres as photon
961 TH1F * hNEM = 0; //Registered as EM
964 TH1F * hNBarReg = 0; //Registeres as photon
965 TH1F * hNBarEM = 0; //Registered as EM
967 //primary - charged hadron (pBar excluded)
968 TH1F * hChargedReg = 0; //Registeres as photon
969 TH1F * hChargedEM = 0; //Registered as EM
972 TH1F * hPbarReg = 0; //Registeres as photon
973 TH1F * hPbarEM = 0; //Registered as EM
976 //Reading histograms from the file
977 TFile * cfile = new TFile("contamination.root","update") ;
979 //read general histograms
980 hPrimary = (TH1F*) cfile->Get("hPrimary") ;
982 hPrimary= new TH1F("hPrimary", "Primary photon spectrum", 100, 0., 5.);
983 hAllRP = (TH1F*)cfile->Get("hAllRP") ;
985 hAllRP = new TH1F("hAllRP","All Reconstructed particles", 100, 0., 5.);
986 hPhot = (TH1F*)cfile->Get("hPhot") ;
988 hPhot = new TH1F("hPhot","All kGAMMA RecParticles",100, 0., 5.);
989 hShape = (TH1F*) cfile->Get("hShape") ;
991 hShape = new TH1F("hShape","All particles with EM shower",100, 0., 5.);
992 hVeto= (TH1F*)cfile->Get("hVeto") ;
994 hVeto = new TH1F("hVeto", "All uncharged particles", 100, 0., 5.);
998 hPhotReg = (TH1F*)cfile->Get("hPhotReg");
1000 hPhotReg = new TH1F("hPhotReg","Photon registered as photon",100, 0., 5.);
1001 hPhotEM =(TH1F*)cfile->Get("hPhotEM");
1003 hPhotEM = new TH1F("hPhotEM", "Photon registered as EM", 100, 0., 5.);
1006 hNReg = (TH1F*)cfile->Get("hNReg");
1008 hNReg = new TH1F("hNReg", "N registered as photon", 100, 0., 5.);
1009 hNEM = (TH1F*)cfile->Get("hNEM");
1011 hNEM = new TH1F("hNEM", "N registered as EM", 100, 0., 5.);
1014 hNBarReg =(TH1F*)cfile->Get("hNBarReg");
1016 hNBarReg = new TH1F("hNBarReg", "NBar registered as photon", 100, 0., 5.);
1017 hNBarEM =(TH1F*)cfile->Get("hNBarEM");
1019 hNBarEM = new TH1F("hNBarEM", "NBar registered as EM", 100, 0., 5.);
1021 //primary - charged hadron (pBar excluded)
1022 hChargedReg = (TH1F*)cfile->Get("hChargedReg");
1024 hChargedReg= new TH1F("hChargedReg", "Charged hadron registered as photon",100, 0., 5.);
1025 hChargedEM = (TH1F*)cfile->Get("hChargedEM");
1027 hChargedEM= new TH1F("hChargedEM","Charged registered as EM",100, 0., 5.);
1030 hPbarReg = (TH1F*)cfile->Get("hPbarReg");
1032 hPbarReg= new TH1F("hPbarReg", "pBar registered as photon",100, 0., 5.);
1033 hPbarEM = (TH1F*)cfile->Get("hPbarEM");
1035 hPbarEM= new TH1F("hPbarEM","Pbar registered as EM",100, 0., 5.);
1038 //Now make some initializations
1040 Int_t counter[8][5] ; //# of registered particles
1042 for(i1 = 0; i1<8; i1++)
1043 for(i2 = 0; i2<5; i2++)
1044 counter[i1][i2] = 0 ;
1048 if (fRunLoader == 0x0)
1050 AliError(Form("Error Loading session"));
1054 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
1057 AliError(Form("Could not obtain the Loader object !"));
1061 if (fRunLoader->TreeE() == 0x0) fRunLoader->LoadHeader();
1062 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
1065 Int_t maxevent = (Int_t)fRunLoader->TreeE()->GetEntries() ;
1066 for ( ievent=0; ievent < maxevent ; ievent++){
1068 fRunLoader->GetEvent(ievent) ;
1070 //DP:Extract vertex position
1071 Double_t vtx[3]={0.,0.,0.} ;
1073 TClonesArray * rp = gime->RecParticles() ;
1075 AliError(Form("Event %d, Can't find RecParticles", ievent)) ;
1078 TClonesArray * ts = gime->TrackSegments() ;
1080 AliError(Form("Event %d, Can't find TrackSegments", ievent)) ;
1083 TObjArray * emcrp = gime->EmcRecPoints() ;
1085 AliError(Form("Event %d, Can't find EmcRecPoints", ievent)) ;
1090 //=========== Make spectrum of the primary photons
1091 const TParticle * primary ;
1093 for( iPrimary = 0 ; iPrimary < fRunLoader->Stack()->GetNprimary() ; iPrimary++){
1094 primary = fRunLoader->Stack()->Particle(iPrimary) ;
1095 Int_t primaryType = primary->GetPdgCode() ;
1096 if( primaryType == 22 ) {
1097 //check, if photons folls onto PHOS
1098 Int_t moduleNumber ;
1099 Double_t primX, primZ ;
1100 phosgeom->ImpactOnEmc(vtx,primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
1102 hPrimary->Fill(primary->Energy()) ;
1108 //========== Now scan over RecParticles
1109 const AliPHOSRecParticle * recParticle ;
1110 Int_t iRecParticle ;
1111 for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast(); iRecParticle++ ){
1112 recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
1113 //fill histo spectrum of all RecParticles
1114 hAllRP->Fill(CorrectedEnergy(recParticle->Energy())) ;
1116 //==========find the closest primary
1117 Int_t moduleNumberRec ;
1118 Double_t recX, recZ ;
1119 phosgeom->ImpactOnEmc(vtx,recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
1121 Double_t minDistance = 100. ;
1122 Int_t closestPrimary = -1 ;
1124 //extract list of primaries: it is stored at EMC RecPoints
1125 Int_t emcIndex = ((AliPHOSTrackSegment *) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
1126 Int_t numberofprimaries ;
1127 Int_t * listofprimaries = ((AliPHOSEmcRecPoint *) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries) ;
1129 Double_t distance = minDistance ;
1131 Double_t dXmin = 0.;
1132 Double_t dZmin = 0. ;
1133 for ( index = 0 ; index < numberofprimaries ; index++){
1134 primary = fRunLoader->Stack()->Particle(listofprimaries[index]) ;
1135 Int_t moduleNumber ;
1136 Double_t primX, primZ ;
1137 phosgeom->ImpactOnEmc(vtx,primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
1138 if(moduleNumberRec == moduleNumber) {
1141 distance = TMath::Sqrt(dX*dX + dZ*dZ) ;
1142 if(minDistance > distance) {
1143 minDistance = distance ;
1146 closestPrimary = listofprimaries[index] ;
1151 //===========define the "type" of closest primary
1152 if(closestPrimary >=0 ){
1153 Int_t primaryCode = -1;
1154 primary = fRunLoader->Stack()->Particle(closestPrimary) ;
1155 Int_t primaryType = primary->GetPdgCode() ;
1156 if(primaryType == 22) // photon ?
1159 if(primaryType == 2112) // neutron
1162 if(primaryType == -2112) // Anti neutron
1165 if(primaryType == -2122) //Anti proton
1168 TParticle tempo(*primary) ;
1169 if(tempo.GetPDG()->Charge())
1173 //==========Now look at the type of RecParticle
1174 Float_t energy = CorrectedEnergy(recParticle->Energy()) ;
1175 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST){
1176 hPhot->Fill(energy ) ;
1177 switch(primaryCode){
1179 hPhotReg->Fill(energy ) ;
1182 hNReg->Fill(energy ) ;
1185 hNBarReg->Fill(energy ) ;
1188 hChargedReg->Fill(energy ) ;
1191 hPbarReg->Fill(energy ) ;
1197 if((recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)||
1198 (recParticle->GetType() == AliPHOSFastRecParticle::kCHARGEDEMFAST)||
1199 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW)||
1200 (recParticle->GetType() == AliPHOSFastRecParticle::kCHARGEDEMSLOW) ){ //with EM shower
1201 hShape->Fill(energy ) ;
1202 switch(primaryCode){
1204 hPhotEM->Fill(energy ) ;
1207 hNEM->Fill(energy ) ;
1210 hNBarEM->Fill(energy ) ;
1213 hChargedEM->Fill(energy ) ;
1216 hPbarEM->Fill(energy ) ;
1223 if((recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)||
1224 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALHAFAST) ||
1225 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW) ||
1226 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALHASLOW) ) //nuetral
1227 hVeto->Fill(energy ) ;
1229 //fill number of primaries identified as ...
1230 if(primaryCode >= 0) // Primary code defined
1231 counter[recParticle->GetType()][primaryCode]++ ;
1235 } // no closest primary found
1239 //=================== SaveHistograms
1241 hPrimary->Write(0,kOverwrite);
1242 hAllRP->Write(0,kOverwrite);
1243 hPhot->Write(0,kOverwrite);
1244 hShape->Write(0,kOverwrite);
1245 hVeto->Write(0,kOverwrite);
1246 hPhotReg->Write(0,kOverwrite);
1247 hPhotEM->Write(0,kOverwrite);
1248 hNReg ->Write(0,kOverwrite);
1249 hNEM ->Write(0,kOverwrite);
1250 hNBarReg ->Write(0,kOverwrite);
1251 hNBarEM ->Write(0,kOverwrite);
1252 hChargedReg ->Write(0,kOverwrite);
1253 hChargedEM ->Write(0,kOverwrite);
1254 hPbarReg ->Write(0,kOverwrite);
1255 hPbarEM ->Write(0,kOverwrite);
1257 cfile->Write(0,kOverwrite);
1263 maxevent = (Int_t)gAlice->TreeE()->GetEntries() ;
1266 message = "Resolutions: Analyzed %d event(s)\n" ;
1268 message += " Primary: Photon Neutron Antineutron Charged hadron AntiProton\n" ;
1269 message += "--------------------------------------------------------------------------------\n" ;
1270 message += " kGAMMA: " ;
1271 message += "%d %d %d %d %d\n" ;
1272 message += " kGAMMAHA: " ;
1273 message += "%d %d %d %d %d\n" ;
1274 message += " kNEUTRALEM: " ;
1275 message += "%d %d %d %d %d\n" ;
1276 message += " kNEUTRALHA: " ;
1277 message += "%d %d %d %d %d\n" ;
1278 message += " kABSURDEM: ";
1279 message += "%d %d %d %d %d\n" ;
1280 message += " kABSURDHA: " ;
1281 message += "%d %d %d %d %d\n" ;
1282 message += " kELECTRON: " ;
1283 message += "%d %d %d %d %d\n" ;
1284 message += " kCHARGEDHA: " ;
1285 message += "%d %d %d %d %d\n" ;
1287 message += "--------------------------------------------------------------------------------" ;
1290 Int_t totalInd = 0 ;
1291 for(i1 = 0; i1<8; i1++)
1292 for(i2 = 0; i2<5; i2++)
1293 totalInd+=counter[i1][i2] ;
1294 message += "Indentified particles: %d" ;
1296 AliInfo(Form(message.Data(), maxevent,
1297 counter[2][0], counter[2][1], counter[2][2], counter[2][3], counter[2][4],
1298 counter[3][0], counter[3][1], counter[3][2], counter[3][3], counter[3][4],
1299 counter[0][0], counter[0][1], counter[0][2], counter[0][3], counter[0][4],
1300 counter[1][0], counter[1][1], counter[1][2], counter[1][3], counter[1][4],
1301 counter[4][0], counter[4][1], counter[4][2], counter[4][3], counter[4][4],
1302 counter[5][0], counter[5][1], counter[5][2], counter[5][3], counter[5][4],
1303 counter[6][0], counter[6][1], counter[6][2], counter[6][3], counter[6][4],
1304 counter[7][0], counter[7][1], counter[7][2], counter[7][3], counter[7][4],