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()
97 // default ctor (useless)
98 fCorrection = 1.2 ; //Value calculated for default parameters of reconstruction
102 //____________________________________________________________________________
103 AliPHOSAnalyze::AliPHOSAnalyze(Text_t * fileName)
105 // ctor: analyze events from root file "name"
106 ffileName = fileName;
107 fCorrection = 1.05 ; //Value calculated for default parameters of reconstruction
108 fRunLoader = AliRunLoader::Open(fileName,"AliPHOSAnalyze");
109 if (fRunLoader == 0x0)
111 AliError(Form("Error Loading session"));
115 //____________________________________________________________________________
116 AliPHOSAnalyze::AliPHOSAnalyze(const AliPHOSAnalyze & ana)
120 ( (AliPHOSAnalyze &)ana ).Copy(*this) ;
123 //____________________________________________________________________________
124 AliPHOSAnalyze::~AliPHOSAnalyze()
129 //____________________________________________________________________________
130 void AliPHOSAnalyze::DrawRecon(Int_t Nevent,Int_t Nmod){
131 //Draws pimary particles and reconstructed
132 //digits, RecPoints, RecPartices etc
133 //for event Nevent in the module Nmod.
135 //========== Create ObjectLoader
136 if (fRunLoader == 0x0)
138 AliError(Form("Error Loading session"));
142 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
145 AliError(Form("Could not obtain the Loader object !"));
150 if(Nevent >= fRunLoader->GetNumberOfEvents() ) {
151 AliError(Form("There is no event %d only %d events available", Nevent, fRunLoader->GetNumberOfEvents() )) ;
154 const AliPHOSGeometry * phosgeom = gime->PHOSGeometry() ;
155 fRunLoader->GetEvent(Nevent);
157 Int_t nx = phosgeom->GetNPhi() ;
158 Int_t nz = phosgeom->GetNZ() ;
159 Float_t * cri= phosgeom->GetEMCAGeometry()->GetCrystalHalfSize() ;
160 Float_t x = nx*cri[0] ;
161 Float_t z = nz*cri[2] ;
162 Int_t nxCPV = (Int_t) (nx*phosgeom->GetPadSizePhi()/(2.*cri[0])) ;
163 Int_t nzCPV = (Int_t) (nz*phosgeom->GetPadSizeZ()/(2.*cri[2])) ;
165 TH2F * emcDigits = (TH2F*) gROOT->FindObject("emcDigits") ;
167 emcDigits->Delete() ;
168 emcDigits = new TH2F("emcDigits","EMC digits", nx,-x,x,nz,-z,z);
169 TH2F * emcSdigits =(TH2F*) gROOT->FindObject("emcSdigits") ;
171 emcSdigits->Delete() ;
172 emcSdigits = new TH2F("emcSdigits","EMC sdigits", nx,-x,x,nz,-z,z);
173 TH2F * emcRecPoints = (TH2F*)gROOT->FindObject("emcRecPoints") ;
175 emcRecPoints->Delete() ;
176 emcRecPoints = new TH2F("emcRecPoints","EMC RecPoints",nx,-x,x,nz,-z,z);
177 TH2F * cpvSdigits =(TH2F*) gROOT->FindObject("cpvSdigits") ;
179 cpvSdigits->Delete() ;
180 cpvSdigits = new TH2F("cpvSdigits","CPV sdigits", nx,-x,x,nz,-z,z);
181 TH2F * cpvDigits = (TH2F*)gROOT->FindObject("cpvDigits") ;
183 cpvDigits->Delete() ;
184 cpvDigits = new TH2F("cpvDigits","CPV digits", nxCPV,-x,x,nzCPV,-z,z) ;
185 TH2F * cpvRecPoints= (TH2F*)gROOT->FindObject("cpvRecPoints") ;
187 cpvRecPoints->Delete() ;
188 cpvRecPoints = new TH2F("cpvRecPoints","CPV RecPoints", nxCPV,-x,x,nzCPV,-z,z) ;
190 TH2F * phot = (TH2F*)gROOT->FindObject("phot") ;
193 phot = new TH2F("phot","Primary Photon", nx,-x,x,nz,-z,z);
194 TH2F * recPhot = (TH2F*)gROOT->FindObject("recPhot") ;
197 recPhot = new TH2F("recPhot","RecParticles with primary Photon",nx,-x,x,nz,-z,z);
200 //Plot Primary Particles
202 if (fRunLoader->Stack() == 0x0) fRunLoader->LoadKinematics("READ");
205 const TParticle * primary ;
207 for ( iPrimary = 0 ; iPrimary < fRunLoader->Stack()->GetNprimary() ; iPrimary++)
209 primary = fRunLoader->Stack()->Particle(iPrimary);
211 Int_t primaryType = primary->GetPdgCode();
212 // if( (primaryType == 211)||(primaryType == -211)||(primaryType == 2212)||(primaryType == -2212)
213 // ||(primaryType == 11)||(primaryType == -11) ) {
214 // Int_t moduleNumber ;
215 // Double_t primX, primZ ;
216 // phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
217 // if(moduleNumber==Nmod)
218 // charg->Fill(primZ,primX,primary->Energy()) ;
220 if( primaryType == 22 ) {
222 Double_t primX, primZ ;
223 phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
224 if(moduleNumber==Nmod)
225 phot->Fill(primZ,primX,primary->Energy()) ;
228 // if( primaryType == -2112 ) {
229 // Int_t moduleNumber ;
230 // Double_t primX, primZ ;
231 // phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
232 // if(moduleNumber==Nmod)
233 // nbar->Fill(primZ,primX,primary->Energy()) ;
240 AliPHOSDigit * sdigit ;
241 const TClonesArray * sdigits = gime->SDigits() ;
242 Int_t nsdig[5] = {0,0,0,0,0} ;
244 for(iSDigit = 0; iSDigit < sdigits->GetEntriesFast() ; iSDigit++)
246 sdigit = (AliPHOSDigit *) sdigits->At(iSDigit) ;
248 phosgeom->AbsToRelNumbering(sdigit->GetId(), relid) ;
250 phosgeom->RelPosInModule(relid,x,z);
251 AliPHOSSDigitizer* sd = dynamic_cast<AliPHOSSDigitizer*>(gime->SDigitizer());
252 Float_t e = sd->Calibrate(sdigit->GetAmp()) ;
253 nsdig[relid[0]-1]++ ;
255 if(relid[1]==0) //EMC
256 emcSdigits->Fill(x,z,e) ;
258 cpvSdigits->Fill(x,z,e) ;
263 message = "Number of EMC + CPV SDigits per module: \n" ;
264 message += "%d %d %d %d %d\n";
265 AliInfo(Form(message.Data(), nsdig[0], nsdig[1], nsdig[2], nsdig[3], nsdig[4] )) ;
269 AliPHOSDigit * digit ;
270 const TClonesArray * digits = gime->Digits();
272 for(iDigit = 0; iDigit < digits->GetEntriesFast(); iDigit++)
274 digit = (AliPHOSDigit *) digits->At(iDigit) ;
276 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
278 phosgeom->RelPosInModule(relid,x,z) ;
279 AliPHOSSDigitizer* sd = dynamic_cast<AliPHOSSDigitizer*>(gime->SDigitizer());
280 Float_t e = sd->Calibrate(digit->GetAmp()) ;
282 if(relid[1]==0) //EMC
283 emcDigits->Fill(x,z,e) ;
285 cpvDigits->Fill(x,z,e) ;
294 TObjArray * emcrp = gime->EmcRecPoints() ;
296 for(irecp = 0; irecp < emcrp->GetEntriesFast() ; irecp ++){
297 AliPHOSEmcRecPoint * emc = (AliPHOSEmcRecPoint *) emcrp->At(irecp) ;
298 if(emc->GetPHOSMod()==Nmod){
299 emc->GetLocalPosition(pos) ;
300 emcRecPoints->Fill(pos.X(),pos.Z(),emc->GetEnergy());
305 TObjArray * cpvrp = gime->CpvRecPoints() ;
307 for(irecp = 0; irecp < cpvrp->GetEntriesFast() ; irecp ++){
308 AliPHOSRecPoint * cpv = (AliPHOSCpvRecPoint *) cpvrp->At(irecp) ;
309 if(cpv->GetPHOSMod()==Nmod){
310 cpv->GetLocalPosition(pos) ;
311 cpvRecPoints->Fill(pos.X(),pos.Z(),cpv->GetEnergy());
317 AliPHOSRecParticle * recParticle ;
319 TClonesArray * rp = gime->RecParticles() ;
320 TClonesArray * ts = gime->TrackSegments() ;
321 if(rp && ts && emcrp) {
322 for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast() ; iRecParticle++ )
324 recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
325 Int_t moduleNumberRec ;
326 Double_t recX, recZ ;
327 phosgeom->ImpactOnEmc(recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
328 if(moduleNumberRec == Nmod){
330 Double_t minDistance = 5. ;
331 Int_t closestPrimary = -1 ;
333 //extract list of primaries: it is stored at EMC RecPoints
334 Int_t emcIndex = ((AliPHOSTrackSegment *) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
335 Int_t numberofprimaries ;
336 Int_t * listofprimaries = ((AliPHOSRecPoint*) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries) ;
338 const TParticle * primary ;
339 Double_t distance = minDistance ;
341 for ( index = 0 ; index < numberofprimaries ; index++){
342 primary = fRunLoader->Stack()->Particle(listofprimaries[index]) ;
344 Double_t primX, primZ ;
345 phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
346 if(moduleNumberRec == moduleNumber)
347 distance = TMath::Sqrt((recX-primX)*(recX-primX)+(recZ-primZ)*(recZ-primZ) ) ;
348 if(minDistance > distance)
350 minDistance = distance ;
351 closestPrimary = listofprimaries[index] ;
355 if(closestPrimary >=0 ){
357 Int_t primaryType = fRunLoader->Stack()->Particle(closestPrimary)->GetPdgCode() ;
360 recPhot->Fill(recZ,recX,recParticle->Energy()) ;
362 // if(primaryType==-2112)
363 // recNbar->Fill(recZ,recX,recParticle->Energy()) ;
370 //Plot made histograms
371 emcSdigits->Draw("box") ;
372 emcDigits->SetLineColor(5) ;
373 emcDigits->Draw("boxsame") ;
374 emcRecPoints->SetLineColor(2) ;
375 emcRecPoints->Draw("boxsame") ;
376 cpvSdigits->SetLineColor(1) ;
377 cpvSdigits->Draw("boxsame") ;
380 //____________________________________________________________________________
381 void AliPHOSAnalyze::Ls(){
382 //lists branches and titles of PHOS-related branches of TreeR, TreeD, TreeS
384 if (fRunLoader == 0x0)
386 AliError(Form("Error Loading session"));
390 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
393 AliError(Form("Could not obtain the Loader object !"));
399 TObjArray * branches;
401 if (gime->TreeS() == 0x0)
403 if (gime->LoadSDigits("READ"))
405 AliError(Form("Problems with loading summable digits"));
409 branches = gime->TreeS()->GetListOfBranches() ;
412 message = "TreeS:\n" ;
413 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
414 TBranch * branch=(TBranch *) branches->At(ibranch) ;
415 if(strstr(branch->GetName(),"PHOS") ){
417 message += branch->GetName() ;
419 message += branch->GetTitle() ;
423 if (gime->TreeD() == 0x0)
425 if (gime->LoadDigits("READ"))
427 AliError(Form("Problems with loading digits"));
432 branches = gime->TreeD()->GetListOfBranches() ;
434 message += "TreeD:\n" ;
435 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
436 TBranch * branch=(TBranch *) branches->At(ibranch) ;
437 if(strstr(branch->GetName(),"PHOS") ) {
439 message += branch->GetName() ;
441 message += branch->GetTitle() ;
446 if (gime->TreeR() == 0x0)
448 if (gime->LoadRecPoints("READ"))
450 AliError(Form("Problems with loading rec points"));
455 branches = gime->TreeR()->GetListOfBranches() ;
457 message += "TreeR: \n" ;
458 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
459 TBranch * branch=(TBranch *) branches->At(ibranch) ;
460 if(strstr(branch->GetName(),"PHOS") ) {
462 message += branch->GetName() ;
464 message += branch->GetTitle() ;
468 AliInfo(Form(message.Data())) ;
470 //____________________________________________________________________________
471 void AliPHOSAnalyze::InvariantMass()
473 // Calculates Real and Mixed invariant mass distributions
474 if (fRunLoader == 0x0)
476 AliError(Form("Error Loading session"));
480 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
483 AliError(Form("Could not obtain the Loader object !"));
487 gime->LoadRecParticles("READ");
489 Int_t nMixedEvents = 4 ; //# of events used for calculation of 'mixed' distribution
493 TFile * mfile = new TFile("invmass.root","update");
495 //========== Reading /Booking Histograms
497 hRealEM = (TH2D*) mfile->Get("hRealEM") ;
499 hRealEM = new TH2D("hRealEM", "Real for EM particles", 250,0.,1.,40,0.,4.) ;
500 TH2D * hRealPhot = 0 ;
502 hRealPhot = (TH2D*)mfile->Get("hRealPhot");
504 hRealPhot = new TH2D("hRealPhot", "Real for kPhoton particles", 250,0.,1.,40,0.,4.) ;
506 TH2D * hMixedEM = 0 ;
507 hMixedEM = (TH2D*) mfile->Get("hMixedEM") ;
509 hMixedEM = new TH2D("hMixedEM", "Mixed for EM particles", 250,0.,1.,40,0.,4.) ;
511 TH2D * hMixedPhot = 0 ;
512 hMixedPhot = (TH2D*) mfile->Get("hMixedPhot") ;
514 hMixedPhot = new TH2D("hMixedPhot","Mixed for kPhoton particles",250,0.,1.,40,0.,4.) ;
517 //reading event and copyng it to TConesArray of all photons
519 TClonesArray * allRecParticleList = new TClonesArray("AliPHOSRecParticle", 1000) ;
520 Int_t * nRecParticles = new Int_t[nMixedEvents] ; // to mark boundaries of each event in the total list
521 for(Int_t index = 0; index < nMixedEvents; index ++)
522 nRecParticles[index] = 0 ;
523 Int_t iRecPhot = 0 ; // number of EM particles in total list
525 //scan over all events
529 if (fRunLoader->TreeE() == 0x0) fRunLoader->LoadHeader();
532 Int_t maxevent = (Int_t)fRunLoader->TreeE()->GetEntries();
533 // for(event = 0; event < gime->MaxEvent(); event++ ){
537 for(event = 0; event < maxevent; event++ ){
538 fRunLoader->GetEvent(event); //will read only TreeR
540 //copy EM RecParticles to the "total" list
541 const AliPHOSRecParticle * recParticle ;
543 TClonesArray * rp = gime->RecParticles() ;
545 AliError(Form("Can't find RecParticles")) ;
549 for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast(); iRecParticle++ )
551 recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
552 if((recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)||
553 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW))
554 new( (*allRecParticleList)[iRecPhot++] ) AliPHOSRecParticle(*recParticle) ;
557 Int_t mevent = event%nMixedEvents ; //event number in the "mixed" cicle
558 nRecParticles[mevent] = iRecPhot-1 ;
560 //check, if it is time to calculate invariant mass?
561 Int_t maxevent = (Int_t)gAlice->TreeE()->GetEntries() ;
562 if((mevent == 0) && (event +1 == maxevent)){
564 // if((mevent == 0) && (event +1 == gime->MaxEvent())){
566 //calculate invariant mass:
568 Int_t nCurEvent = 0 ;
570 for(irp1 = 0; irp1 < allRecParticleList->GetEntries()-1; irp1++){
571 AliPHOSRecParticle * rp1 = (AliPHOSRecParticle *)allRecParticleList->At(irp1) ;
573 for(irp2 = irp1+1; irp2 < allRecParticleList->GetEntries(); irp2++){
574 AliPHOSRecParticle * rp2 = (AliPHOSRecParticle *)allRecParticleList->At(irp2) ;
577 invMass = (rp1->Energy()+rp2->Energy())*(rp1->Energy()+rp2->Energy())-
578 (rp1->Px()+rp2->Px())*(rp1->Px()+rp2->Px())-
579 (rp1->Py()+rp2->Py())*(rp1->Py()+rp2->Py())-
580 (rp1->Pz()+rp2->Pz())*(rp1->Pz()+rp2->Pz()) ;
583 invMass = TMath::Sqrt(invMass);
586 pt = TMath::Sqrt((rp1->Px()+rp2->Px() )*( rp1->Px()+rp2->Px() ) +
587 (rp1->Py()+rp2->Py() )*( rp1->Py()+rp2->Py() ) );
589 if(irp1 > nRecParticles[nCurEvent])
592 if(irp2 <= nRecParticles[nCurEvent]){ //'Real' event
593 hRealEM->Fill(invMass,pt);
594 if((rp1->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)&&
595 (rp2->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST) )
596 hRealPhot->Fill(invMass,pt);
599 hMixedEM->Fill(invMass,pt);
600 if((rp1->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)&&
601 (rp2->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST) )
602 hMixedPhot->Fill(invMass,pt);
605 } //loop over second rp
606 }//loop over first rp
608 //Make some cleanings
609 for(Int_t index = 0; index < nMixedEvents; index ++)
610 nRecParticles[index] = 0 ;
612 allRecParticleList->Clear() ;
616 delete allRecParticleList ;
621 hRealEM->Write(0,kOverwrite) ;
622 hRealPhot->Write(0,kOverwrite) ;
623 hMixedEM->Write(0,kOverwrite) ;
624 hMixedPhot->Write(0,kOverwrite) ;
629 delete nRecParticles;
633 //____________________________________________________________________________
634 void AliPHOSAnalyze::EnergyResolution()
636 //fills two dimentional histo: energy of primary vs. energy of reconstructed
638 TH2F * hAllEnergy = 0 ; //all reconstructed with primary photon
639 TH2F * hPhotEnergy= 0 ; //kGamma with primary photon
640 TH2F * hEMEnergy = 0 ; //electromagnetic with primary photon
642 //opening file and reading histograms if any
643 TFile * efile = new TFile("energy.root","update");
645 hAllEnergy = (TH2F*)efile->Get("hAllEnergy") ;
647 hAllEnergy = new TH2F("hAllEnergy", "Energy of any RP with primary photon",100, 0., 5., 100, 0., 5.);
649 hPhotEnergy =(TH2F*) efile->Get("hPhotEnergy") ;
651 hPhotEnergy = new TH2F("hPhotEnergy", "Energy of kGAMMA with primary photon",100, 0., 5., 100, 0., 5.);
653 hEMEnergy =(TH2F*) efile->Get("hEMEnergy");
655 hEMEnergy = new TH2F("hEMEnergy", "Energy of EM with primary photon", 100, 0., 5., 100, 0., 5.);
658 if (fRunLoader == 0x0)
660 AliError(Form("Error Loading session"));
664 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
667 AliError(Form("Could not obtain the Loader object !"));
672 const AliPHOSGeometry * phosgeom = gime->PHOSGeometry();
675 Int_t maxevent = (Int_t)fRunLoader->TreeE()->GetEntries();
677 fRunLoader->LoadKinematics("READ");
678 gime->LoadTracks("READ");
680 for ( ievent=0; ievent < maxevent ; ievent++){
682 //read the current event
683 fRunLoader->GetEvent(ievent) ;
685 const AliPHOSRecParticle * recParticle ;
687 TClonesArray * rp = gime->RecParticles() ;
689 AliError(Form("Event %d, Can't find RecParticles ", ievent)) ;
692 TClonesArray * ts = gime->TrackSegments() ;
694 AliError(Form("Event %d, Can't find TrackSegments", ievent)) ;
697 TObjArray * emcrp = gime->EmcRecPoints() ;
699 AliError(Form("Event %d, Can't find EmcRecPoints")) ;
703 for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast() ;iRecParticle++ ){
704 recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
706 //find the closest primary
707 Int_t moduleNumberRec ;
708 Double_t recX, recZ ;
709 phosgeom->ImpactOnEmc(recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
711 Double_t minDistance = 100. ;
712 Int_t closestPrimary = -1 ;
714 //extract list of primaries: it is stored at EMC RecPoints
715 Int_t emcIndex = ((AliPHOSTrackSegment*) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
716 Int_t numberofprimaries ;
717 Int_t * listofprimaries = ((AliPHOSEmcRecPoint*) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries) ;
720 const TParticle * primary ;
721 Double_t distance = minDistance ;
724 Double_t dZmin = 0. ;
725 for ( index = 0 ; index < numberofprimaries ; index++){
727 primary = fRunLoader->Stack()->Particle(listofprimaries[index]) ;
730 Double_t primX, primZ ;
731 phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
732 if(moduleNumberRec == moduleNumber) {
735 distance = TMath::Sqrt(dX*dX + dZ*dZ) ;
736 if(minDistance > distance) {
737 minDistance = distance ;
740 closestPrimary = listofprimaries[index] ;
745 //if found primary, fill histograms
746 if(closestPrimary >=0 ){
747 const TParticle * primary = fRunLoader->Stack()->Particle(closestPrimary) ;
748 if(primary->GetPdgCode() == 22){
749 hAllEnergy->Fill(primary->Energy(), recParticle->Energy()) ;
750 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST){
751 hPhotEnergy->Fill(primary->Energy(), recParticle->Energy() ) ;
752 hEMEnergy->Fill(primary->Energy(), recParticle->Energy() ) ;
755 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW)
756 hEMEnergy->Fill(primary->Energy(), recParticle->Energy() ) ;
762 //write filled histograms
764 hAllEnergy->Write(0,kOverwrite) ;
765 hPhotEnergy->Write(0,kOverwrite) ;
766 hEMEnergy->Write(0,kOverwrite) ;
772 //____________________________________________________________________________
773 void AliPHOSAnalyze::PositionResolution()
775 //fills two dimentional histo: energy vs. primary - reconstructed distance
779 TH2F * hAllPosition = 0; // Position of any RP with primary photon
780 TH2F * hPhotPosition = 0; // Position of kGAMMA with primary photon
781 TH2F * hEMPosition = 0; // Position of EM with primary photon
783 TH1F * hAllPositionX = 0; // X-Position Resolution of photons with photon primary
784 TH1F * hAllPositionZ = 0; // Z-Position Resolution of photons with photon primary
787 //opening file and reading histograms if any
788 TFile * pfile = new TFile("position.root","update");
790 hAllPosition = (TH2F*)pfile->Get("hAllPosition");
791 if(hAllPosition == 0)
792 hAllPosition = new TH2F("hAllPosition",
793 "Position of any RP with primary photon",100, 0., 5., 100, 0., 5.);
794 hPhotPosition= (TH2F*)pfile->Get("hPhotPosition");
795 if(hPhotPosition == 0)
796 hPhotPosition = new TH2F("hPhotPosition",
797 "Position of kGAMMA with primary photon",100, 0., 5., 100, 0., 5.);
798 hEMPosition= (TH2F*)pfile->Get("hEMPosition") ;
800 hEMPosition = new TH2F("hEMPosition",
801 "Position of EM with primary photon", 100, 0., 5., 100, 0., 5.);
802 hAllPositionX = (TH1F*)pfile->Get("hAllPositionX") ;
803 if(hAllPositionX == 0)
804 hAllPositionX = new TH1F("hAllPositionX",
805 "Delta X of any RP with primary photon",100, -2., 2.);
806 hAllPositionZ =(TH1F*) pfile->Get("hAllPositionZ") ;
807 if(hAllPositionZ == 0)
808 hAllPositionZ = new TH1F("hAllPositionZ",
809 "Delta X of any RP with primary photon",100, -2., 2.);
811 if (fRunLoader == 0x0)
813 AliError(Form("Error Loading session"));
817 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
820 AliError(Form("Could not obtain the Loader object !"));
824 if (fRunLoader->TreeE() == 0x0) fRunLoader->LoadHeader();
826 const AliPHOSGeometry * phosgeom = gime->PHOSGeometry() ;
829 Int_t maxevent = (Int_t)fRunLoader->TreeE()->GetEntries() ;
830 for ( ievent=0; ievent < maxevent ; ievent++){
832 //read the current event
833 fRunLoader->GetEvent(ievent) ;
834 TClonesArray * rp = gime->RecParticles() ;
836 AliError(Form("Event %d, Can't find RecParticles", ievent)) ;
839 TClonesArray * ts = gime->TrackSegments() ;
841 AliError(Form("Event %d, Can't find TrackSegments", ievent)) ;
844 TObjArray * emcrp = gime->EmcRecPoints() ;
846 AliError(Form("Event %d, Can't find EmcRecPoints", ievent)) ;
851 const AliPHOSRecParticle * recParticle ;
853 for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast(); iRecParticle++ ){
854 recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
856 //find the closest primary
857 Int_t moduleNumberRec ;
858 Double_t recX, recZ ;
859 phosgeom->ImpactOnEmc(recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
861 Double_t minDistance = 100. ;
862 Int_t closestPrimary = -1 ;
864 //extract list of primaries: it is stored at EMC RecPoints
865 Int_t emcIndex = ((AliPHOSTrackSegment*) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
866 Int_t numberofprimaries ;
867 Int_t * listofprimaries = ((AliPHOSEmcRecPoint *) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries) ;
870 const TParticle * primary ;
871 Double_t distance = minDistance ;
872 Double_t dX = 1000; // incredible number
873 Double_t dZ = 1000; // for the case if no primary will be found
875 Double_t dZmin = 0. ;
876 for ( index = 0 ; index < numberofprimaries ; index++){
877 primary = fRunLoader->Stack()->Particle(listofprimaries[index]) ;
879 Double_t primX, primZ ;
880 phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
881 if(moduleNumberRec == moduleNumber) {
884 distance = TMath::Sqrt(dX*dX + dZ*dZ) ;
885 if(minDistance > distance) {
886 minDistance = distance ;
889 closestPrimary = listofprimaries[index] ;
894 //if found primary, fill histograms
895 if(closestPrimary >=0 ){
896 const TParticle * primary = fRunLoader->Stack()->Particle(closestPrimary) ;
897 if(primary->GetPdgCode() == 22){
898 hAllPosition->Fill(primary->Energy(), minDistance) ;
899 hAllPositionX->Fill(primary->Energy(), dX) ;
900 hAllPositionZ->Fill(primary->Energy(), dZ) ;
901 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST){
902 hPhotPosition->Fill(primary->Energy(), minDistance ) ;
903 hEMPosition->Fill(primary->Energy(), minDistance ) ;
906 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW)
907 hEMPosition->Fill(primary->Energy(), minDistance ) ;
913 //Write output histgrams
915 hAllPosition->Write(0,kOverwrite) ;
916 hAllPositionX->Write(0,kOverwrite) ;
917 hAllPositionZ->Write(0,kOverwrite) ;
918 hPhotPosition->Write(0,kOverwrite) ;
919 hEMPosition->Write(0,kOverwrite) ;
926 //____________________________________________________________________________
927 void AliPHOSAnalyze::Contamination(){
928 // fills spectra of primary photons and several kinds of
929 // reconstructed particles, so that analyzing them one can
930 // estimate conatmination, efficiency of registration etc.
932 //define several general histograms
933 TH1F * hPrimary = 0; //spectrum (P_t distribution) of primary photons
934 TH1F * hAllRP = 0; //spectrum of all RecParticles in PHOS
935 TH1F * hPhot = 0; //spectrum of kGAMMA RecParticles
936 TH1F * hShape = 0; //spectrum of all EM RecParticles
937 TH1F * hVeto = 0; //spectrum of all neutral RecParticles
939 //Now separate histograms in accoradance with primary
941 TH1F * hPhotReg = 0; //Registeres as photon
942 TH1F * hPhotEM = 0; //Registered as EM
945 TH1F * hNReg = 0; //Registeres as photon
946 TH1F * hNEM = 0; //Registered as EM
949 TH1F * hNBarReg = 0; //Registeres as photon
950 TH1F * hNBarEM = 0; //Registered as EM
952 //primary - charged hadron (pBar excluded)
953 TH1F * hChargedReg = 0; //Registeres as photon
954 TH1F * hChargedEM = 0; //Registered as EM
957 TH1F * hPbarReg = 0; //Registeres as photon
958 TH1F * hPbarEM = 0; //Registered as EM
961 //Reading histograms from the file
962 TFile * cfile = new TFile("contamination.root","update") ;
964 //read general histograms
965 hPrimary = (TH1F*) cfile->Get("hPrimary") ;
967 hPrimary= new TH1F("hPrimary", "Primary photon spectrum", 100, 0., 5.);
968 hAllRP = (TH1F*)cfile->Get("hAllRP") ;
970 hAllRP = new TH1F("hAllRP","All Reconstructed particles", 100, 0., 5.);
971 hPhot = (TH1F*)cfile->Get("hPhot") ;
973 hPhot = new TH1F("hPhot","All kGAMMA RecParticles",100, 0., 5.);
974 hShape = (TH1F*) cfile->Get("hShape") ;
976 hShape = new TH1F("hShape","All particles with EM shower",100, 0., 5.);
977 hVeto= (TH1F*)cfile->Get("hVeto") ;
979 hVeto = new TH1F("hVeto", "All uncharged particles", 100, 0., 5.);
983 hPhotReg = (TH1F*)cfile->Get("hPhotReg");
985 hPhotReg = new TH1F("hPhotReg","Photon registered as photon",100, 0., 5.);
986 hPhotEM =(TH1F*)cfile->Get("hPhotEM");
988 hPhotEM = new TH1F("hPhotEM", "Photon registered as EM", 100, 0., 5.);
991 hNReg = (TH1F*)cfile->Get("hNReg");
993 hNReg = new TH1F("hNReg", "N registered as photon", 100, 0., 5.);
994 hNEM = (TH1F*)cfile->Get("hNEM");
996 hNEM = new TH1F("hNEM", "N registered as EM", 100, 0., 5.);
999 hNBarReg =(TH1F*)cfile->Get("hNBarReg");
1001 hNBarReg = new TH1F("hNBarReg", "NBar registered as photon", 100, 0., 5.);
1002 hNBarEM =(TH1F*)cfile->Get("hNBarEM");
1004 hNBarEM = new TH1F("hNBarEM", "NBar registered as EM", 100, 0., 5.);
1006 //primary - charged hadron (pBar excluded)
1007 hChargedReg = (TH1F*)cfile->Get("hChargedReg");
1009 hChargedReg= new TH1F("hChargedReg", "Charged hadron registered as photon",100, 0., 5.);
1010 hChargedEM = (TH1F*)cfile->Get("hChargedEM");
1012 hChargedEM= new TH1F("hChargedEM","Charged registered as EM",100, 0., 5.);
1015 hPbarReg = (TH1F*)cfile->Get("hPbarReg");
1017 hPbarReg= new TH1F("hPbarReg", "pBar registered as photon",100, 0., 5.);
1018 hPbarEM = (TH1F*)cfile->Get("hPbarEM");
1020 hPbarEM= new TH1F("hPbarEM","Pbar registered as EM",100, 0., 5.);
1023 //Now make some initializations
1025 Int_t counter[8][5] ; //# of registered particles
1027 for(i1 = 0; i1<8; i1++)
1028 for(i2 = 0; i2<5; i2++)
1029 counter[i1][i2] = 0 ;
1033 if (fRunLoader == 0x0)
1035 AliError(Form("Error Loading session"));
1039 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(fRunLoader->GetLoader("PHOSLoader"));
1042 AliError(Form("Could not obtain the Loader object !"));
1046 if (fRunLoader->TreeE() == 0x0) fRunLoader->LoadHeader();
1047 const AliPHOSGeometry * phosgeom = gime->PHOSGeometry() ;
1050 Int_t maxevent = (Int_t)fRunLoader->TreeE()->GetEntries() ;
1051 for ( ievent=0; ievent < maxevent ; ievent++){
1053 fRunLoader->GetEvent(ievent) ;
1055 TClonesArray * rp = gime->RecParticles() ;
1057 AliError(Form("Event %d, Can't find RecParticles", ievent)) ;
1060 TClonesArray * ts = gime->TrackSegments() ;
1062 AliError(Form("Event %d, Can't find TrackSegments", ievent)) ;
1065 TObjArray * emcrp = gime->EmcRecPoints() ;
1067 AliError(Form("Event %d, Can't find EmcRecPoints", ievent)) ;
1072 //=========== Make spectrum of the primary photons
1073 const TParticle * primary ;
1075 for( iPrimary = 0 ; iPrimary < fRunLoader->Stack()->GetNprimary() ; iPrimary++){
1076 primary = fRunLoader->Stack()->Particle(iPrimary) ;
1077 Int_t primaryType = primary->GetPdgCode() ;
1078 if( primaryType == 22 ) {
1079 //check, if photons folls onto PHOS
1080 Int_t moduleNumber ;
1081 Double_t primX, primZ ;
1082 phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
1084 hPrimary->Fill(primary->Energy()) ;
1090 //========== Now scan over RecParticles
1091 const AliPHOSRecParticle * recParticle ;
1092 Int_t iRecParticle ;
1093 for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast(); iRecParticle++ ){
1094 recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
1095 //fill histo spectrum of all RecParticles
1096 hAllRP->Fill(CorrectedEnergy(recParticle->Energy())) ;
1098 //==========find the closest primary
1099 Int_t moduleNumberRec ;
1100 Double_t recX, recZ ;
1101 phosgeom->ImpactOnEmc(recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
1103 Double_t minDistance = 100. ;
1104 Int_t closestPrimary = -1 ;
1106 //extract list of primaries: it is stored at EMC RecPoints
1107 Int_t emcIndex = ((AliPHOSTrackSegment *) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
1108 Int_t numberofprimaries ;
1109 Int_t * listofprimaries = ((AliPHOSEmcRecPoint *) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries) ;
1111 const TParticle * primary ;
1112 Double_t distance = minDistance ;
1114 Double_t dXmin = 0.;
1115 Double_t dZmin = 0. ;
1116 for ( index = 0 ; index < numberofprimaries ; index++){
1117 primary = fRunLoader->Stack()->Particle(listofprimaries[index]) ;
1118 Int_t moduleNumber ;
1119 Double_t primX, primZ ;
1120 phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
1121 if(moduleNumberRec == moduleNumber) {
1124 distance = TMath::Sqrt(dX*dX + dZ*dZ) ;
1125 if(minDistance > distance) {
1126 minDistance = distance ;
1129 closestPrimary = listofprimaries[index] ;
1134 //===========define the "type" of closest primary
1135 if(closestPrimary >=0 ){
1136 Int_t primaryCode = -1;
1137 const TParticle * primary = fRunLoader->Stack()->Particle(closestPrimary) ;
1138 Int_t primaryType = primary->GetPdgCode() ;
1139 if(primaryType == 22) // photon ?
1142 if(primaryType == 2112) // neutron
1145 if(primaryType == -2112) // Anti neutron
1148 if(primaryType == -2122) //Anti proton
1151 TParticle tempo(*primary) ;
1152 if(tempo.GetPDG()->Charge())
1156 //==========Now look at the type of RecParticle
1157 Float_t energy = CorrectedEnergy(recParticle->Energy()) ;
1158 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST){
1159 hPhot->Fill(energy ) ;
1160 switch(primaryCode){
1162 hPhotReg->Fill(energy ) ;
1165 hNReg->Fill(energy ) ;
1168 hNBarReg->Fill(energy ) ;
1171 hChargedReg->Fill(energy ) ;
1174 hPbarReg->Fill(energy ) ;
1180 if((recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)||
1181 (recParticle->GetType() == AliPHOSFastRecParticle::kCHARGEDEMFAST)||
1182 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW)||
1183 (recParticle->GetType() == AliPHOSFastRecParticle::kCHARGEDEMSLOW) ){ //with EM shower
1184 hShape->Fill(energy ) ;
1185 switch(primaryCode){
1187 hPhotEM->Fill(energy ) ;
1190 hNEM->Fill(energy ) ;
1193 hNBarEM->Fill(energy ) ;
1196 hChargedEM->Fill(energy ) ;
1199 hPbarEM->Fill(energy ) ;
1206 if((recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)||
1207 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALHAFAST) ||
1208 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW) ||
1209 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALHASLOW) ) //nuetral
1210 hVeto->Fill(energy ) ;
1212 //fill number of primaries identified as ...
1213 if(primaryCode >= 0) // Primary code defined
1214 counter[recParticle->GetType()][primaryCode]++ ;
1218 } // no closest primary found
1222 //=================== SaveHistograms
1224 hPrimary->Write(0,kOverwrite);
1225 hAllRP->Write(0,kOverwrite);
1226 hPhot->Write(0,kOverwrite);
1227 hShape->Write(0,kOverwrite);
1228 hVeto->Write(0,kOverwrite);
1229 hPhotReg->Write(0,kOverwrite);
1230 hPhotEM->Write(0,kOverwrite);
1231 hNReg ->Write(0,kOverwrite);
1232 hNEM ->Write(0,kOverwrite);
1233 hNBarReg ->Write(0,kOverwrite);
1234 hNBarEM ->Write(0,kOverwrite);
1235 hChargedReg ->Write(0,kOverwrite);
1236 hChargedEM ->Write(0,kOverwrite);
1237 hPbarReg ->Write(0,kOverwrite);
1238 hPbarEM ->Write(0,kOverwrite);
1240 cfile->Write(0,kOverwrite);
1246 maxevent = (Int_t)gAlice->TreeE()->GetEntries() ;
1249 message = "Resolutions: Analyzed %d event(s)\n" ;
1251 message += " Primary: Photon Neutron Antineutron Charged hadron AntiProton\n" ;
1252 message += "--------------------------------------------------------------------------------\n" ;
1253 message += " kGAMMA: " ;
1254 message += "%d %d %d %d %d\n" ;
1255 message += " kGAMMAHA: " ;
1256 message += "%d %d %d %d %d\n" ;
1257 message += " kNEUTRALEM: " ;
1258 message += "%d %d %d %d %d\n" ;
1259 message += " kNEUTRALHA: " ;
1260 message += "%d %d %d %d %d\n" ;
1261 message += " kABSURDEM: ";
1262 message += "%d %d %d %d %d\n" ;
1263 message += " kABSURDHA: " ;
1264 message += "%d %d %d %d %d\n" ;
1265 message += " kELECTRON: " ;
1266 message += "%d %d %d %d %d\n" ;
1267 message += " kCHARGEDHA: " ;
1268 message += "%d %d %d %d %d\n" ;
1270 message += "--------------------------------------------------------------------------------" ;
1273 Int_t totalInd = 0 ;
1274 for(i1 = 0; i1<8; i1++)
1275 for(i2 = 0; i2<5; i2++)
1276 totalInd+=counter[i1][i2] ;
1277 message += "Indentified particles: %d" ;
1279 AliInfo(Form(message.Data(), maxevent,
1280 counter[2][0], counter[2][1], counter[2][2], counter[2][3], counter[2][4],
1281 counter[3][0], counter[3][1], counter[3][2], counter[3][3], counter[3][4],
1282 counter[0][0], counter[0][1], counter[0][2], counter[0][3], counter[0][4],
1283 counter[1][0], counter[1][1], counter[1][2], counter[1][3], counter[1][4],
1284 counter[4][0], counter[4][1], counter[4][2], counter[4][3], counter[4][4],
1285 counter[5][0], counter[5][1], counter[5][2], counter[5][3], counter[5][4],
1286 counter[6][0], counter[6][1], counter[6][2], counter[6][3], counter[6][4],
1287 counter[7][0], counter[7][1], counter[7][2], counter[7][3], counter[7][4],