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 const AliPHOSGeometry * phosgeom = gime->PHOSGeometry() ;
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,x,z);
263 AliPHOSSDigitizer* sd = dynamic_cast<AliPHOSSDigitizer*>(gime->SDigitizer());
264 Float_t e = sd->Calibrate(sdigit->GetAmp()) ;
265 nsdig[relid[0]-1]++ ;
267 if(relid[1]==0) //EMC
268 emcSdigits->Fill(x,z,e) ;
270 cpvSdigits->Fill(x,z,e) ;
275 message = "Number of EMC + CPV SDigits per module: \n" ;
276 message += "%d %d %d %d %d\n";
277 AliInfo(Form(message.Data(), nsdig[0], nsdig[1], nsdig[2], nsdig[3], nsdig[4] )) ;
281 AliPHOSDigit * digit ;
282 const TClonesArray * digits = gime->Digits();
284 for(iDigit = 0; iDigit < digits->GetEntriesFast(); iDigit++)
286 digit = (AliPHOSDigit *) digits->At(iDigit) ;
288 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
290 phosgeom->RelPosInModule(relid,x,z) ;
291 AliPHOSSDigitizer* sd = dynamic_cast<AliPHOSSDigitizer*>(gime->SDigitizer());
292 Float_t e = sd->Calibrate(digit->GetAmp()) ;
294 if(relid[1]==0) //EMC
295 emcDigits->Fill(x,z,e) ;
297 cpvDigits->Fill(x,z,e) ;
306 TObjArray * emcrp = gime->EmcRecPoints() ;
308 for(irecp = 0; irecp < emcrp->GetEntriesFast() ; irecp ++){
309 AliPHOSEmcRecPoint * emc = (AliPHOSEmcRecPoint *) emcrp->At(irecp) ;
310 if(emc->GetPHOSMod()==Nmod){
311 emc->GetLocalPosition(pos) ;
312 emcRecPoints->Fill(pos.X(),pos.Z(),emc->GetEnergy());
317 TObjArray * cpvrp = gime->CpvRecPoints() ;
319 for(irecp = 0; irecp < cpvrp->GetEntriesFast() ; irecp ++){
320 AliPHOSRecPoint * cpv = (AliPHOSCpvRecPoint *) cpvrp->At(irecp) ;
321 if(cpv->GetPHOSMod()==Nmod){
322 cpv->GetLocalPosition(pos) ;
323 cpvRecPoints->Fill(pos.X(),pos.Z(),cpv->GetEnergy());
329 AliPHOSRecParticle * recParticle ;
331 TClonesArray * rp = gime->RecParticles() ;
332 TClonesArray * ts = gime->TrackSegments() ;
333 if(rp && ts && emcrp) {
334 for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast() ; iRecParticle++ )
336 recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
337 Int_t moduleNumberRec ;
338 Double_t recX, recZ ;
339 phosgeom->ImpactOnEmc(vtx,recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
340 if(moduleNumberRec == Nmod){
342 Double_t minDistance = 5. ;
343 Int_t closestPrimary = -1 ;
345 //extract list of primaries: it is stored at EMC RecPoints
346 Int_t emcIndex = ((AliPHOSTrackSegment *) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
347 Int_t numberofprimaries ;
348 Int_t * listofprimaries = ((AliPHOSRecPoint*) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries) ;
350 const TParticle * primary ;
351 Double_t distance = minDistance ;
353 for ( index = 0 ; index < numberofprimaries ; index++){
354 primary = fRunLoader->Stack()->Particle(listofprimaries[index]) ;
356 Double_t primX, primZ ;
357 phosgeom->ImpactOnEmc(vtx,primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
358 if(moduleNumberRec == moduleNumber)
359 distance = TMath::Sqrt((recX-primX)*(recX-primX)+(recZ-primZ)*(recZ-primZ) ) ;
360 if(minDistance > distance)
362 minDistance = distance ;
363 closestPrimary = listofprimaries[index] ;
367 if(closestPrimary >=0 ){
369 Int_t primaryType = fRunLoader->Stack()->Particle(closestPrimary)->GetPdgCode() ;
372 recPhot->Fill(recZ,recX,recParticle->Energy()) ;
374 // if(primaryType==-2112)
375 // recNbar->Fill(recZ,recX,recParticle->Energy()) ;
382 //Plot made histograms
383 emcSdigits->Draw("box") ;
384 emcDigits->SetLineColor(5) ;
385 emcDigits->Draw("boxsame") ;
386 emcRecPoints->SetLineColor(2) ;
387 emcRecPoints->Draw("boxsame") ;
388 cpvSdigits->SetLineColor(1) ;
389 cpvSdigits->Draw("boxsame") ;
392 //____________________________________________________________________________
393 void AliPHOSAnalyze::Ls(){
394 //lists branches and titles of PHOS-related branches of TreeR, TreeD, TreeS
396 if (fRunLoader == 0x0)
398 AliError(Form("Error Loading session"));
402 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
405 AliError(Form("Could not obtain the Loader object !"));
411 TObjArray * branches;
413 if (gime->TreeS() == 0x0)
415 if (gime->LoadSDigits("READ"))
417 AliError(Form("Problems with loading summable digits"));
421 branches = gime->TreeS()->GetListOfBranches() ;
424 message = "TreeS:\n" ;
425 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
426 TBranch * branch=(TBranch *) branches->At(ibranch) ;
427 if(strstr(branch->GetName(),"PHOS") ){
429 message += branch->GetName() ;
431 message += branch->GetTitle() ;
435 if (gime->TreeD() == 0x0)
437 if (gime->LoadDigits("READ"))
439 AliError(Form("Problems with loading digits"));
444 branches = gime->TreeD()->GetListOfBranches() ;
446 message += "TreeD:\n" ;
447 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
448 TBranch * branch=(TBranch *) branches->At(ibranch) ;
449 if(strstr(branch->GetName(),"PHOS") ) {
451 message += branch->GetName() ;
453 message += branch->GetTitle() ;
458 if (gime->TreeR() == 0x0)
460 if (gime->LoadRecPoints("READ"))
462 AliError(Form("Problems with loading rec points"));
467 branches = gime->TreeR()->GetListOfBranches() ;
469 message += "TreeR: \n" ;
470 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
471 TBranch * branch=(TBranch *) branches->At(ibranch) ;
472 if(strstr(branch->GetName(),"PHOS") ) {
474 message += branch->GetName() ;
476 message += branch->GetTitle() ;
480 AliInfo(Form(message.Data())) ;
482 //____________________________________________________________________________
483 void AliPHOSAnalyze::InvariantMass()
485 // Calculates Real and Mixed invariant mass distributions
486 if (fRunLoader == 0x0)
488 AliError(Form("Error Loading session"));
492 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
495 AliError(Form("Could not obtain the Loader object !"));
499 gime->LoadRecParticles("READ");
501 Int_t nMixedEvents = 4 ; //# of events used for calculation of 'mixed' distribution
505 TFile * mfile = new TFile("invmass.root","update");
507 //========== Reading /Booking Histograms
509 hRealEM = (TH2D*) mfile->Get("hRealEM") ;
511 hRealEM = new TH2D("hRealEM", "Real for EM particles", 250,0.,1.,40,0.,4.) ;
512 TH2D * hRealPhot = 0 ;
514 hRealPhot = (TH2D*)mfile->Get("hRealPhot");
516 hRealPhot = new TH2D("hRealPhot", "Real for kPhoton particles", 250,0.,1.,40,0.,4.) ;
518 TH2D * hMixedEM = 0 ;
519 hMixedEM = (TH2D*) mfile->Get("hMixedEM") ;
521 hMixedEM = new TH2D("hMixedEM", "Mixed for EM particles", 250,0.,1.,40,0.,4.) ;
523 TH2D * hMixedPhot = 0 ;
524 hMixedPhot = (TH2D*) mfile->Get("hMixedPhot") ;
526 hMixedPhot = new TH2D("hMixedPhot","Mixed for kPhoton particles",250,0.,1.,40,0.,4.) ;
529 //reading event and copyng it to TConesArray of all photons
531 TClonesArray * allRecParticleList = new TClonesArray("AliPHOSRecParticle", 1000) ;
532 Int_t * nRecParticles = new Int_t[nMixedEvents] ; // to mark boundaries of each event in the total list
533 for(Int_t index = 0; index < nMixedEvents; index ++)
534 nRecParticles[index] = 0 ;
535 Int_t iRecPhot = 0 ; // number of EM particles in total list
537 //scan over all events
541 if (fRunLoader->TreeE() == 0x0) fRunLoader->LoadHeader();
544 Int_t maxevent = (Int_t)fRunLoader->TreeE()->GetEntries();
545 // for(event = 0; event < gime->MaxEvent(); event++ ){
549 for(event = 0; event < maxevent; event++ ){
550 fRunLoader->GetEvent(event); //will read only TreeR
552 //copy EM RecParticles to the "total" list
553 const AliPHOSRecParticle * recParticle ;
555 TClonesArray * rp = gime->RecParticles() ;
557 AliError(Form("Can't find RecParticles")) ;
561 for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast(); iRecParticle++ )
563 recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
564 if((recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)||
565 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW))
566 new( (*allRecParticleList)[iRecPhot++] ) AliPHOSRecParticle(*recParticle) ;
569 Int_t mevent = event%nMixedEvents ; //event number in the "mixed" cicle
570 nRecParticles[mevent] = iRecPhot-1 ;
572 //check, if it is time to calculate invariant mass?
573 Int_t maxevent = (Int_t)gAlice->TreeE()->GetEntries() ;
574 if((mevent == 0) && (event +1 == maxevent)){
576 // if((mevent == 0) && (event +1 == gime->MaxEvent())){
578 //calculate invariant mass:
580 Int_t nCurEvent = 0 ;
582 for(irp1 = 0; irp1 < allRecParticleList->GetEntries()-1; irp1++){
583 AliPHOSRecParticle * rp1 = (AliPHOSRecParticle *)allRecParticleList->At(irp1) ;
585 for(irp2 = irp1+1; irp2 < allRecParticleList->GetEntries(); irp2++){
586 AliPHOSRecParticle * rp2 = (AliPHOSRecParticle *)allRecParticleList->At(irp2) ;
589 invMass = (rp1->Energy()+rp2->Energy())*(rp1->Energy()+rp2->Energy())-
590 (rp1->Px()+rp2->Px())*(rp1->Px()+rp2->Px())-
591 (rp1->Py()+rp2->Py())*(rp1->Py()+rp2->Py())-
592 (rp1->Pz()+rp2->Pz())*(rp1->Pz()+rp2->Pz()) ;
595 invMass = TMath::Sqrt(invMass);
598 pt = TMath::Sqrt((rp1->Px()+rp2->Px() )*( rp1->Px()+rp2->Px() ) +
599 (rp1->Py()+rp2->Py() )*( rp1->Py()+rp2->Py() ) );
601 if(irp1 > nRecParticles[nCurEvent])
604 if(irp2 <= nRecParticles[nCurEvent]){ //'Real' event
605 hRealEM->Fill(invMass,pt);
606 if((rp1->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)&&
607 (rp2->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST) )
608 hRealPhot->Fill(invMass,pt);
611 hMixedEM->Fill(invMass,pt);
612 if((rp1->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)&&
613 (rp2->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST) )
614 hMixedPhot->Fill(invMass,pt);
617 } //loop over second rp
618 }//loop over first rp
620 //Make some cleanings
621 for(Int_t index = 0; index < nMixedEvents; index ++)
622 nRecParticles[index] = 0 ;
624 allRecParticleList->Clear() ;
628 delete allRecParticleList ;
633 hRealEM->Write(0,kOverwrite) ;
634 hRealPhot->Write(0,kOverwrite) ;
635 hMixedEM->Write(0,kOverwrite) ;
636 hMixedPhot->Write(0,kOverwrite) ;
641 delete nRecParticles;
645 //____________________________________________________________________________
646 void AliPHOSAnalyze::EnergyResolution()
648 //fills two dimentional histo: energy of primary vs. energy of reconstructed
650 TH2F * hAllEnergy = 0 ; //all reconstructed with primary photon
651 TH2F * hPhotEnergy= 0 ; //kGamma with primary photon
652 TH2F * hEMEnergy = 0 ; //electromagnetic with primary photon
654 //opening file and reading histograms if any
655 TFile * efile = new TFile("energy.root","update");
657 hAllEnergy = (TH2F*)efile->Get("hAllEnergy") ;
659 hAllEnergy = new TH2F("hAllEnergy", "Energy of any RP with primary photon",100, 0., 5., 100, 0., 5.);
661 hPhotEnergy =(TH2F*) efile->Get("hPhotEnergy") ;
663 hPhotEnergy = new TH2F("hPhotEnergy", "Energy of kGAMMA with primary photon",100, 0., 5., 100, 0., 5.);
665 hEMEnergy =(TH2F*) efile->Get("hEMEnergy");
667 hEMEnergy = new TH2F("hEMEnergy", "Energy of EM with primary photon", 100, 0., 5., 100, 0., 5.);
670 if (fRunLoader == 0x0)
672 AliError(Form("Error Loading session"));
676 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
679 AliError(Form("Could not obtain the Loader object !"));
684 const AliPHOSGeometry * phosgeom = gime->PHOSGeometry();
687 Int_t maxevent = (Int_t)fRunLoader->TreeE()->GetEntries();
689 fRunLoader->LoadKinematics("READ");
690 gime->LoadTracks("READ");
692 for ( ievent=0; ievent < maxevent ; ievent++){
694 //read the current event
695 fRunLoader->GetEvent(ievent) ;
697 Double_t vtx[3]={0.,0.,0.} ;
699 const AliPHOSRecParticle * recParticle ;
701 TClonesArray * rp = gime->RecParticles() ;
703 AliError(Form("Event %d, Can't find RecParticles ", ievent)) ;
706 TClonesArray * ts = gime->TrackSegments() ;
708 AliError(Form("Event %d, Can't find TrackSegments", ievent)) ;
711 TObjArray * emcrp = gime->EmcRecPoints() ;
713 AliError(Form("Event %d, Can't find EmcRecPoints")) ;
717 for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast() ;iRecParticle++ ){
718 recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
720 //find the closest primary
721 Int_t moduleNumberRec ;
722 Double_t recX, recZ ;
723 phosgeom->ImpactOnEmc(vtx,recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
725 Double_t minDistance = 100. ;
726 Int_t closestPrimary = -1 ;
728 //extract list of primaries: it is stored at EMC RecPoints
729 Int_t emcIndex = ((AliPHOSTrackSegment*) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
730 Int_t numberofprimaries ;
731 Int_t * listofprimaries = ((AliPHOSEmcRecPoint*) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries) ;
734 const TParticle * primary ;
735 Double_t distance = minDistance ;
738 Double_t dZmin = 0. ;
739 for ( index = 0 ; index < numberofprimaries ; index++){
741 primary = fRunLoader->Stack()->Particle(listofprimaries[index]) ;
744 Double_t primX, primZ ;
745 phosgeom->ImpactOnEmc(vtx,primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
746 if(moduleNumberRec == moduleNumber) {
749 distance = TMath::Sqrt(dX*dX + dZ*dZ) ;
750 if(minDistance > distance) {
751 minDistance = distance ;
754 closestPrimary = listofprimaries[index] ;
759 //if found primary, fill histograms
760 if(closestPrimary >=0 ){
761 const TParticle * primary = fRunLoader->Stack()->Particle(closestPrimary) ;
762 if(primary->GetPdgCode() == 22){
763 hAllEnergy->Fill(primary->Energy(), recParticle->Energy()) ;
764 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST){
765 hPhotEnergy->Fill(primary->Energy(), recParticle->Energy() ) ;
766 hEMEnergy->Fill(primary->Energy(), recParticle->Energy() ) ;
769 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW)
770 hEMEnergy->Fill(primary->Energy(), recParticle->Energy() ) ;
776 //write filled histograms
778 hAllEnergy->Write(0,kOverwrite) ;
779 hPhotEnergy->Write(0,kOverwrite) ;
780 hEMEnergy->Write(0,kOverwrite) ;
786 //____________________________________________________________________________
787 void AliPHOSAnalyze::PositionResolution()
789 //fills two dimentional histo: energy vs. primary - reconstructed distance
793 TH2F * hAllPosition = 0; // Position of any RP with primary photon
794 TH2F * hPhotPosition = 0; // Position of kGAMMA with primary photon
795 TH2F * hEMPosition = 0; // Position of EM with primary photon
797 TH1F * hAllPositionX = 0; // X-Position Resolution of photons with photon primary
798 TH1F * hAllPositionZ = 0; // Z-Position Resolution of photons with photon primary
801 //opening file and reading histograms if any
802 TFile * pfile = new TFile("position.root","update");
804 hAllPosition = (TH2F*)pfile->Get("hAllPosition");
805 if(hAllPosition == 0)
806 hAllPosition = new TH2F("hAllPosition",
807 "Position of any RP with primary photon",100, 0., 5., 100, 0., 5.);
808 hPhotPosition= (TH2F*)pfile->Get("hPhotPosition");
809 if(hPhotPosition == 0)
810 hPhotPosition = new TH2F("hPhotPosition",
811 "Position of kGAMMA with primary photon",100, 0., 5., 100, 0., 5.);
812 hEMPosition= (TH2F*)pfile->Get("hEMPosition") ;
814 hEMPosition = new TH2F("hEMPosition",
815 "Position of EM with primary photon", 100, 0., 5., 100, 0., 5.);
816 hAllPositionX = (TH1F*)pfile->Get("hAllPositionX") ;
817 if(hAllPositionX == 0)
818 hAllPositionX = new TH1F("hAllPositionX",
819 "Delta X of any RP with primary photon",100, -2., 2.);
820 hAllPositionZ =(TH1F*) pfile->Get("hAllPositionZ") ;
821 if(hAllPositionZ == 0)
822 hAllPositionZ = new TH1F("hAllPositionZ",
823 "Delta X of any RP with primary photon",100, -2., 2.);
825 if (fRunLoader == 0x0)
827 AliError(Form("Error Loading session"));
831 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
834 AliError(Form("Could not obtain the Loader object !"));
838 if (fRunLoader->TreeE() == 0x0) fRunLoader->LoadHeader();
840 const AliPHOSGeometry * phosgeom = gime->PHOSGeometry() ;
843 Int_t maxevent = (Int_t)fRunLoader->TreeE()->GetEntries() ;
844 for ( ievent=0; ievent < maxevent ; ievent++){
846 //read the current event
847 fRunLoader->GetEvent(ievent) ;
849 //DP:Extract vertex position
850 Double_t vtx[3]={0.,0.,0.} ;
852 TClonesArray * rp = gime->RecParticles() ;
854 AliError(Form("Event %d, Can't find RecParticles", ievent)) ;
857 TClonesArray * ts = gime->TrackSegments() ;
859 AliError(Form("Event %d, Can't find TrackSegments", ievent)) ;
862 TObjArray * emcrp = gime->EmcRecPoints() ;
864 AliError(Form("Event %d, Can't find EmcRecPoints", ievent)) ;
869 const AliPHOSRecParticle * recParticle ;
871 for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast(); iRecParticle++ ){
872 recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
874 //find the closest primary
875 Int_t moduleNumberRec ;
876 Double_t recX, recZ ;
877 phosgeom->ImpactOnEmc(vtx,recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
879 Double_t minDistance = 100. ;
880 Int_t closestPrimary = -1 ;
882 //extract list of primaries: it is stored at EMC RecPoints
883 Int_t emcIndex = ((AliPHOSTrackSegment*) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
884 Int_t numberofprimaries ;
885 Int_t * listofprimaries = ((AliPHOSEmcRecPoint *) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries) ;
888 const TParticle * primary ;
889 Double_t distance = minDistance ;
890 Double_t dX = 1000; // incredible number
891 Double_t dZ = 1000; // for the case if no primary will be found
893 Double_t dZmin = 0. ;
894 for ( index = 0 ; index < numberofprimaries ; index++){
895 primary = fRunLoader->Stack()->Particle(listofprimaries[index]) ;
897 Double_t primX, primZ ;
898 phosgeom->ImpactOnEmc(vtx,primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
899 if(moduleNumberRec == moduleNumber) {
902 distance = TMath::Sqrt(dX*dX + dZ*dZ) ;
903 if(minDistance > distance) {
904 minDistance = distance ;
907 closestPrimary = listofprimaries[index] ;
912 //if found primary, fill histograms
913 if(closestPrimary >=0 ){
914 const TParticle * primary = fRunLoader->Stack()->Particle(closestPrimary) ;
915 if(primary->GetPdgCode() == 22){
916 hAllPosition->Fill(primary->Energy(), minDistance) ;
917 hAllPositionX->Fill(primary->Energy(), dX) ;
918 hAllPositionZ->Fill(primary->Energy(), dZ) ;
919 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST){
920 hPhotPosition->Fill(primary->Energy(), minDistance ) ;
921 hEMPosition->Fill(primary->Energy(), minDistance ) ;
924 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW)
925 hEMPosition->Fill(primary->Energy(), minDistance ) ;
931 //Write output histgrams
933 hAllPosition->Write(0,kOverwrite) ;
934 hAllPositionX->Write(0,kOverwrite) ;
935 hAllPositionZ->Write(0,kOverwrite) ;
936 hPhotPosition->Write(0,kOverwrite) ;
937 hEMPosition->Write(0,kOverwrite) ;
944 //____________________________________________________________________________
945 void AliPHOSAnalyze::Contamination(){
946 // fills spectra of primary photons and several kinds of
947 // reconstructed particles, so that analyzing them one can
948 // estimate conatmination, efficiency of registration etc.
950 //define several general histograms
951 TH1F * hPrimary = 0; //spectrum (P_t distribution) of primary photons
952 TH1F * hAllRP = 0; //spectrum of all RecParticles in PHOS
953 TH1F * hPhot = 0; //spectrum of kGAMMA RecParticles
954 TH1F * hShape = 0; //spectrum of all EM RecParticles
955 TH1F * hVeto = 0; //spectrum of all neutral RecParticles
957 //Now separate histograms in accoradance with primary
959 TH1F * hPhotReg = 0; //Registeres as photon
960 TH1F * hPhotEM = 0; //Registered as EM
963 TH1F * hNReg = 0; //Registeres as photon
964 TH1F * hNEM = 0; //Registered as EM
967 TH1F * hNBarReg = 0; //Registeres as photon
968 TH1F * hNBarEM = 0; //Registered as EM
970 //primary - charged hadron (pBar excluded)
971 TH1F * hChargedReg = 0; //Registeres as photon
972 TH1F * hChargedEM = 0; //Registered as EM
975 TH1F * hPbarReg = 0; //Registeres as photon
976 TH1F * hPbarEM = 0; //Registered as EM
979 //Reading histograms from the file
980 TFile * cfile = new TFile("contamination.root","update") ;
982 //read general histograms
983 hPrimary = (TH1F*) cfile->Get("hPrimary") ;
985 hPrimary= new TH1F("hPrimary", "Primary photon spectrum", 100, 0., 5.);
986 hAllRP = (TH1F*)cfile->Get("hAllRP") ;
988 hAllRP = new TH1F("hAllRP","All Reconstructed particles", 100, 0., 5.);
989 hPhot = (TH1F*)cfile->Get("hPhot") ;
991 hPhot = new TH1F("hPhot","All kGAMMA RecParticles",100, 0., 5.);
992 hShape = (TH1F*) cfile->Get("hShape") ;
994 hShape = new TH1F("hShape","All particles with EM shower",100, 0., 5.);
995 hVeto= (TH1F*)cfile->Get("hVeto") ;
997 hVeto = new TH1F("hVeto", "All uncharged particles", 100, 0., 5.);
1001 hPhotReg = (TH1F*)cfile->Get("hPhotReg");
1003 hPhotReg = new TH1F("hPhotReg","Photon registered as photon",100, 0., 5.);
1004 hPhotEM =(TH1F*)cfile->Get("hPhotEM");
1006 hPhotEM = new TH1F("hPhotEM", "Photon registered as EM", 100, 0., 5.);
1009 hNReg = (TH1F*)cfile->Get("hNReg");
1011 hNReg = new TH1F("hNReg", "N registered as photon", 100, 0., 5.);
1012 hNEM = (TH1F*)cfile->Get("hNEM");
1014 hNEM = new TH1F("hNEM", "N registered as EM", 100, 0., 5.);
1017 hNBarReg =(TH1F*)cfile->Get("hNBarReg");
1019 hNBarReg = new TH1F("hNBarReg", "NBar registered as photon", 100, 0., 5.);
1020 hNBarEM =(TH1F*)cfile->Get("hNBarEM");
1022 hNBarEM = new TH1F("hNBarEM", "NBar registered as EM", 100, 0., 5.);
1024 //primary - charged hadron (pBar excluded)
1025 hChargedReg = (TH1F*)cfile->Get("hChargedReg");
1027 hChargedReg= new TH1F("hChargedReg", "Charged hadron registered as photon",100, 0., 5.);
1028 hChargedEM = (TH1F*)cfile->Get("hChargedEM");
1030 hChargedEM= new TH1F("hChargedEM","Charged registered as EM",100, 0., 5.);
1033 hPbarReg = (TH1F*)cfile->Get("hPbarReg");
1035 hPbarReg= new TH1F("hPbarReg", "pBar registered as photon",100, 0., 5.);
1036 hPbarEM = (TH1F*)cfile->Get("hPbarEM");
1038 hPbarEM= new TH1F("hPbarEM","Pbar registered as EM",100, 0., 5.);
1041 //Now make some initializations
1043 Int_t counter[8][5] ; //# of registered particles
1045 for(i1 = 0; i1<8; i1++)
1046 for(i2 = 0; i2<5; i2++)
1047 counter[i1][i2] = 0 ;
1051 if (fRunLoader == 0x0)
1053 AliError(Form("Error Loading session"));
1057 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
1060 AliError(Form("Could not obtain the Loader object !"));
1064 if (fRunLoader->TreeE() == 0x0) fRunLoader->LoadHeader();
1065 const AliPHOSGeometry * phosgeom = gime->PHOSGeometry() ;
1068 Int_t maxevent = (Int_t)fRunLoader->TreeE()->GetEntries() ;
1069 for ( ievent=0; ievent < maxevent ; ievent++){
1071 fRunLoader->GetEvent(ievent) ;
1073 //DP:Extract vertex position
1074 Double_t vtx[3]={0.,0.,0.} ;
1076 TClonesArray * rp = gime->RecParticles() ;
1078 AliError(Form("Event %d, Can't find RecParticles", ievent)) ;
1081 TClonesArray * ts = gime->TrackSegments() ;
1083 AliError(Form("Event %d, Can't find TrackSegments", ievent)) ;
1086 TObjArray * emcrp = gime->EmcRecPoints() ;
1088 AliError(Form("Event %d, Can't find EmcRecPoints", ievent)) ;
1093 //=========== Make spectrum of the primary photons
1094 const TParticle * primary ;
1096 for( iPrimary = 0 ; iPrimary < fRunLoader->Stack()->GetNprimary() ; iPrimary++){
1097 primary = fRunLoader->Stack()->Particle(iPrimary) ;
1098 Int_t primaryType = primary->GetPdgCode() ;
1099 if( primaryType == 22 ) {
1100 //check, if photons folls onto PHOS
1101 Int_t moduleNumber ;
1102 Double_t primX, primZ ;
1103 phosgeom->ImpactOnEmc(vtx,primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
1105 hPrimary->Fill(primary->Energy()) ;
1111 //========== Now scan over RecParticles
1112 const AliPHOSRecParticle * recParticle ;
1113 Int_t iRecParticle ;
1114 for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast(); iRecParticle++ ){
1115 recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
1116 //fill histo spectrum of all RecParticles
1117 hAllRP->Fill(CorrectedEnergy(recParticle->Energy())) ;
1119 //==========find the closest primary
1120 Int_t moduleNumberRec ;
1121 Double_t recX, recZ ;
1122 phosgeom->ImpactOnEmc(vtx,recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
1124 Double_t minDistance = 100. ;
1125 Int_t closestPrimary = -1 ;
1127 //extract list of primaries: it is stored at EMC RecPoints
1128 Int_t emcIndex = ((AliPHOSTrackSegment *) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
1129 Int_t numberofprimaries ;
1130 Int_t * listofprimaries = ((AliPHOSEmcRecPoint *) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries) ;
1132 const TParticle * primary ;
1133 Double_t distance = minDistance ;
1135 Double_t dXmin = 0.;
1136 Double_t dZmin = 0. ;
1137 for ( index = 0 ; index < numberofprimaries ; index++){
1138 primary = fRunLoader->Stack()->Particle(listofprimaries[index]) ;
1139 Int_t moduleNumber ;
1140 Double_t primX, primZ ;
1141 phosgeom->ImpactOnEmc(vtx,primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
1142 if(moduleNumberRec == moduleNumber) {
1145 distance = TMath::Sqrt(dX*dX + dZ*dZ) ;
1146 if(minDistance > distance) {
1147 minDistance = distance ;
1150 closestPrimary = listofprimaries[index] ;
1155 //===========define the "type" of closest primary
1156 if(closestPrimary >=0 ){
1157 Int_t primaryCode = -1;
1158 const TParticle * primary = fRunLoader->Stack()->Particle(closestPrimary) ;
1159 Int_t primaryType = primary->GetPdgCode() ;
1160 if(primaryType == 22) // photon ?
1163 if(primaryType == 2112) // neutron
1166 if(primaryType == -2112) // Anti neutron
1169 if(primaryType == -2122) //Anti proton
1172 TParticle tempo(*primary) ;
1173 if(tempo.GetPDG()->Charge())
1177 //==========Now look at the type of RecParticle
1178 Float_t energy = CorrectedEnergy(recParticle->Energy()) ;
1179 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST){
1180 hPhot->Fill(energy ) ;
1181 switch(primaryCode){
1183 hPhotReg->Fill(energy ) ;
1186 hNReg->Fill(energy ) ;
1189 hNBarReg->Fill(energy ) ;
1192 hChargedReg->Fill(energy ) ;
1195 hPbarReg->Fill(energy ) ;
1201 if((recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)||
1202 (recParticle->GetType() == AliPHOSFastRecParticle::kCHARGEDEMFAST)||
1203 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW)||
1204 (recParticle->GetType() == AliPHOSFastRecParticle::kCHARGEDEMSLOW) ){ //with EM shower
1205 hShape->Fill(energy ) ;
1206 switch(primaryCode){
1208 hPhotEM->Fill(energy ) ;
1211 hNEM->Fill(energy ) ;
1214 hNBarEM->Fill(energy ) ;
1217 hChargedEM->Fill(energy ) ;
1220 hPbarEM->Fill(energy ) ;
1227 if((recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)||
1228 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALHAFAST) ||
1229 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW) ||
1230 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALHASLOW) ) //nuetral
1231 hVeto->Fill(energy ) ;
1233 //fill number of primaries identified as ...
1234 if(primaryCode >= 0) // Primary code defined
1235 counter[recParticle->GetType()][primaryCode]++ ;
1239 } // no closest primary found
1243 //=================== SaveHistograms
1245 hPrimary->Write(0,kOverwrite);
1246 hAllRP->Write(0,kOverwrite);
1247 hPhot->Write(0,kOverwrite);
1248 hShape->Write(0,kOverwrite);
1249 hVeto->Write(0,kOverwrite);
1250 hPhotReg->Write(0,kOverwrite);
1251 hPhotEM->Write(0,kOverwrite);
1252 hNReg ->Write(0,kOverwrite);
1253 hNEM ->Write(0,kOverwrite);
1254 hNBarReg ->Write(0,kOverwrite);
1255 hNBarEM ->Write(0,kOverwrite);
1256 hChargedReg ->Write(0,kOverwrite);
1257 hChargedEM ->Write(0,kOverwrite);
1258 hPbarReg ->Write(0,kOverwrite);
1259 hPbarEM ->Write(0,kOverwrite);
1261 cfile->Write(0,kOverwrite);
1267 maxevent = (Int_t)gAlice->TreeE()->GetEntries() ;
1270 message = "Resolutions: Analyzed %d event(s)\n" ;
1272 message += " Primary: Photon Neutron Antineutron Charged hadron AntiProton\n" ;
1273 message += "--------------------------------------------------------------------------------\n" ;
1274 message += " kGAMMA: " ;
1275 message += "%d %d %d %d %d\n" ;
1276 message += " kGAMMAHA: " ;
1277 message += "%d %d %d %d %d\n" ;
1278 message += " kNEUTRALEM: " ;
1279 message += "%d %d %d %d %d\n" ;
1280 message += " kNEUTRALHA: " ;
1281 message += "%d %d %d %d %d\n" ;
1282 message += " kABSURDEM: ";
1283 message += "%d %d %d %d %d\n" ;
1284 message += " kABSURDHA: " ;
1285 message += "%d %d %d %d %d\n" ;
1286 message += " kELECTRON: " ;
1287 message += "%d %d %d %d %d\n" ;
1288 message += " kCHARGEDHA: " ;
1289 message += "%d %d %d %d %d\n" ;
1291 message += "--------------------------------------------------------------------------------" ;
1294 Int_t totalInd = 0 ;
1295 for(i1 = 0; i1<8; i1++)
1296 for(i2 = 0; i2<5; i2++)
1297 totalInd+=counter[i1][i2] ;
1298 message += "Indentified particles: %d" ;
1300 AliInfo(Form(message.Data(), maxevent,
1301 counter[2][0], counter[2][1], counter[2][2], counter[2][3], counter[2][4],
1302 counter[3][0], counter[3][1], counter[3][2], counter[3][3], counter[3][4],
1303 counter[0][0], counter[0][1], counter[0][2], counter[0][3], counter[0][4],
1304 counter[1][0], counter[1][1], counter[1][2], counter[1][3], counter[1][4],
1305 counter[4][0], counter[4][1], counter[4][2], counter[4][3], counter[4][4],
1306 counter[5][0], counter[5][1], counter[5][2], counter[5][3], counter[5][4],
1307 counter[6][0], counter[6][1], counter[6][2], counter[6][3], counter[6][4],
1308 counter[7][0], counter[7][1], counter[7][2], counter[7][3], counter[7][4],