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.
33 // root [0] AliPHOSAnalyze * a = new AliPHOSAnalyze("galice.root")
34 // // set the file you want to analyse
35 // root [1] a->DrawRecon(1,3)
36 // // plot RecObjects, made in event 1, PHOS module 3
37 // root [2] a->DrawRecon(1,3,"PHOSRP","another PID")
38 // // plot RecObjets made in the event 1, PHOS module 3,
39 // // produced in the another reconstruction pass,
40 // // which produced PHOS RecParticles ("PHOSRP") with
41 // // title "another PID".
42 // root [3] a->InvariantMass()
43 // // Calculates "REAL" and "MIXED" invariant mass
44 // // distributions of kGAMMA and (kGAMMA+kNEUTRALEM)
45 // // and APPENDS this to the file "invmass.root"
46 // root [4] a->PositionResolution()
47 // // calculates two dimentional histos: energy of the primary
48 // // photon vs distance betwin incedence point and reconstructed
49 // // poisition. One can analyse the produced file position.root
50 // // with macro PhotonPosition.C
51 // root [5] a->EnergyResolution()
52 // // calculates two dimentional histos: energy of the primary
53 // // photon vs energy of the reconstructed particle. One can
54 // // analyse the produced file energy.root
55 // // with macro PhotonEnergy.C
56 // root [6] a->Contamination()
57 // // fills spectra of primary photons and several kinds of
58 // // reconstructed particles, so that analyzing them one can
59 // // estimate conatmination, efficiency of registration etc.
61 //*-- Author: Dmitri Peressounko (SUBATECH & RRC Kurchatov Institute)
62 //////////////////////////////////////////////////////////////////////////////
65 // --- ROOT system ---
71 #include "TParticle.h"
72 #include "TClonesArray.h"
76 // --- Standard library ---
80 // --- AliRoot header files ---
83 #include "AliPHOSv1.h"
84 #include "AliPHOSAnalyze.h"
85 #include "AliPHOSDigit.h"
86 #include "AliPHOSSDigitizer.h"
87 #include "AliPHOSTrackSegment.h"
88 #include "AliPHOSRecParticle.h"
89 #include "AliPHOSIndexToObject.h"
90 #include "AliPHOSCpvRecPoint.h"
91 #include "AliPHOSPpsdRecPoint.h"
94 ClassImp(AliPHOSAnalyze)
96 //____________________________________________________________________________
97 AliPHOSAnalyze::AliPHOSAnalyze()
99 // default ctor (useless)
100 fCorrection = 1.2 ; //Value calculated for default parameters of reconstruction
103 //____________________________________________________________________________
104 AliPHOSAnalyze::AliPHOSAnalyze(Text_t * fileName)
106 // ctor: analyze events from root file "name"
107 ffileName = fileName ;
108 fCorrection = 1.05 ; //Value calculated for default parameters of reconstruction
109 fObjGetter = 0 ; // should be instantiated
112 //____________________________________________________________________________
113 AliPHOSAnalyze::AliPHOSAnalyze(const AliPHOSAnalyze & ana)
116 ( (AliPHOSAnalyze &)ana ).Copy(*this) ;
119 //____________________________________________________________________________
120 AliPHOSAnalyze::~AliPHOSAnalyze()
124 if(fPHOS) {delete fPHOS ; fPHOS =0 ;}
127 //____________________________________________________________________________
128 void AliPHOSAnalyze::DrawRecon(Int_t Nevent,Int_t Nmod,const char * branchName,const char* branchTitle){
129 //Draws pimary particles and reconstructed
130 //digits, RecPoints, RecPartices etc
131 //for event Nevent in the module Nmod.
133 TH2F * digitOccupancy = new TH2F("digitOccupancy","EMC digits", 64,-71.,71.,64,-71.,71.);
134 TH2F * sdigitOccupancy = new TH2F("sdigitOccupancy","EMC sdigits", 64,-71.,71.,64,-71.,71.);
135 TH2F * emcOccupancy = new TH2F("emcOccupancy","EMC RecPoints",64,-71.,71.,64,-71.,71.);
136 TH2F * ppsdUp = new TH2F("ppsdUp","PPSD Up digits", 128,-71.,71.,128,-71.,71.) ;
137 TH2F * ppsdUpCl = new TH2F("ppsdUpCl","PPSD Up RecPoints",128,-71.,71.,128,-71.,71.) ;
138 TH2F * ppsdLow = new TH2F("ppsdLow","PPSD Low digits", 128,-71.,71.,128,-71.,71.) ;
139 TH2F * ppsdLowCl = new TH2F("ppsdLowCl","PPSD Low RecPoints",128,-71.,71.,128,-71.,71.) ;
140 TH2F * nbar = new TH2F("nbar","Primary nbar", 64,-71.,71.,64,-71.,71.);
141 TH2F * phot = new TH2F("phot","Primary Photon", 64,-71.,71.,64,-71.,71.);
142 TH2F * charg = new TH2F("charg","Primary charged",64,-71.,71.,64,-71.,71.);
143 TH2F * recPhot = new TH2F("recPhot","RecParticles with primary Photon",64,-71.,71.,64,-71.,71.);
144 TH2F * recNbar = new TH2F("recNbar","RecParticles with primary Nbar", 64,-71.,71.,64,-71.,71.);
146 //========== Create IndexToObjectGetter
147 fObjGetter = AliPHOSIndexToObject::GetInstance(ffileName.Data(),branchName,branchTitle) ;
148 fObjGetter->GetEvent(Nevent);
150 fPHOS = (AliPHOSv1 *)gAlice->GetDetector("PHOS") ;
151 fGeom = AliPHOSGeometry::GetInstance( fPHOS->GetGeometry()->GetName(), fPHOS->GetGeometry()->GetTitle() );
153 //Plot Primary Particles
154 TParticle * primary ;
156 for ( iPrimary = 0 ; iPrimary < fObjGetter->GimeNPrimaries() ; iPrimary++)
158 primary = fObjGetter->GimePrimary(iPrimary) ;
159 Int_t primaryType = primary->GetPdgCode() ;
160 if( (primaryType == 211)||(primaryType == -211)||(primaryType == 2212)||(primaryType == -2212) ) {
162 Double_t primX, primZ ;
163 fGeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
164 if(moduleNumber==Nmod)
165 charg->Fill(primZ,primX,primary->Energy()) ;
167 if( primaryType == 22 ) {
169 Double_t primX, primZ ;
170 fGeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
171 if(moduleNumber==Nmod)
172 phot->Fill(primZ,primX,primary->Energy()) ;
175 if( primaryType == -2112 ) {
177 Double_t primX, primZ ;
178 fGeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
179 if(moduleNumber==Nmod)
180 nbar->Fill(primZ,primX,primary->Energy()) ;
186 AliPHOSDigit * sdigit ;
188 for(iSDigit = 0; iSDigit < fObjGetter->GimeNSDigits(); iSDigit++)
190 sdigit = fObjGetter->GimeSDigit(iSDigit) ;
192 fGeom->AbsToRelNumbering(sdigit->GetId(), relid) ;
194 fGeom->RelPosInModule(relid,x,z) ;
195 Float_t e = fObjGetter->GimeSDigitizer()->Calibrate(sdigit->GetAmp()) ;
197 if(relid[1]==0) //EMC
198 sdigitOccupancy->Fill(x,z,e) ;
199 if((relid[1]>0)&&(relid[1]<17))
200 ppsdUp->Fill(x,z,e) ;
202 ppsdLow->Fill(x,z,e) ;
208 AliPHOSDigit * digit ;
209 for(iDigit = 0; iDigit < fObjGetter->GimeNDigits(); iDigit++)
211 digit = fObjGetter->GimeDigit(iDigit) ;
213 fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
215 fGeom->RelPosInModule(relid,x,z) ;
216 Float_t e = fObjGetter->GimeSDigitizer()->Calibrate(digit->GetAmp()) ;
218 if(relid[1]==0) //EMC
219 digitOccupancy->Fill(x,z,e) ;
220 if((relid[1]>0)&&(relid[1]<17))
221 ppsdUp->Fill(x,z,e) ;
223 ppsdLow->Fill(x,z,e) ;
232 for(irecp = 0; irecp < fObjGetter->GimeNEmcRecPoints() ; irecp ++){
233 AliPHOSEmcRecPoint * emc= fObjGetter->GimeEmcRecPoint(irecp) ;
234 if(emc->GetPHOSMod()==Nmod){
235 emc->GetLocalPosition(pos) ;
236 emcOccupancy->Fill(pos.X(),pos.Z(),emc->GetEnergy());
241 for(irecp = 0; irecp < fObjGetter->GimeNCpvRecPoints() ; irecp ++){
242 AliPHOSRecPoint * cpv = fObjGetter->GimeCpvRecPoint(irecp) ;
243 if((strcmp(cpv->ClassName(),"AliPHOSPpsdRecPoint" )) == 0){ // PPSD Rec Point
244 AliPHOSPpsdRecPoint * ppsd = (AliPHOSPpsdRecPoint*) cpv ;
245 ppsd->GetLocalPosition(pos) ;
246 if(ppsd->GetPHOSMod()==Nmod){
247 ppsd->GetLocalPosition(pos) ;
249 ppsdUpCl->Fill(pos.X(),pos.Z(),ppsd->GetEnergy());
251 ppsdLowCl->Fill(pos.X(),pos.Z(),ppsd->GetEnergy());
255 AliPHOSCpvRecPoint * cpv1 = (AliPHOSCpvRecPoint*) cpv ;
256 if(cpv1->GetPHOSMod()==Nmod){
257 cpv1->GetLocalPosition(pos) ;
258 ppsdUpCl->Fill(pos.X(),pos.Z(),cpv1->GetEnergy());
265 AliPHOSRecParticle * recParticle ;
267 for(iRecParticle = 0; iRecParticle < fObjGetter->GimeNRecParticles() ;iRecParticle++ )
269 recParticle = fObjGetter->GimeRecParticle(iRecParticle) ;
270 Int_t moduleNumberRec ;
271 Double_t recX, recZ ;
272 fGeom->ImpactOnEmc(recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
273 if(moduleNumberRec == Nmod){
275 Double_t minDistance = 5. ;
276 Int_t closestPrimary = -1 ;
279 //extract list of primaries: it is stored at EMC RecPoints
280 Int_t emcIndex = fObjGetter->GimeTrackSegment(recParticle->GetPHOSTSIndex())->GetEmcIndex() ;
281 Int_t numberofprimaries ;
282 Int_t * listofprimaries = fObjGetter->GimeEmcRecPoint(emcIndex)->GetPrimaries(numberofprimaries) ;
284 TParticle * primary ;
285 Double_t distance = minDistance ;
287 for ( index = 0 ; index < numberofprimaries ; index++){
288 primary = fObjGetter->GimePrimary(listofprimaries[index]) ;
290 Double_t primX, primZ ;
291 fGeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
292 if(moduleNumberRec == moduleNumber)
293 distance = TMath::Sqrt((recX-primX)*(recX-primX)+(recZ-primZ)*(recZ-primZ) ) ;
294 if(minDistance > distance)
296 minDistance = distance ;
297 closestPrimary = listofprimaries[index] ;
301 if(closestPrimary >=0 ){
303 Int_t primaryType = fObjGetter->GimePrimary(closestPrimary)->GetPdgCode() ;
306 recPhot->Fill(recZ,recX,recParticle->Energy()) ;
308 if(primaryType==-2112)
309 recNbar->Fill(recZ,recX,recParticle->Energy()) ;
315 //Plot made histograms
316 digitOccupancy->Draw("box") ;
317 sdigitOccupancy->SetLineColor(5) ;
318 sdigitOccupancy->Draw("box") ;
319 emcOccupancy->SetLineColor(2) ;
320 emcOccupancy->Draw("boxsame") ;
321 ppsdUp->SetLineColor(3) ;
322 ppsdUp->Draw("boxsame") ;
323 ppsdLow->SetLineColor(4) ;
324 ppsdLow->Draw("boxsame") ;
325 phot->SetLineColor(8) ;
326 phot->Draw("boxsame") ;
327 nbar->SetLineColor(6) ;
328 nbar->Draw("boxsame") ;
331 //____________________________________________________________________________
332 void AliPHOSAnalyze::Ls(){
333 //lists branches and titles of PHOS-related branches of TreeR, TreeD, TreeS
336 fObjGetter = AliPHOSIndexToObject::GetInstance(ffileName.Data()) ;
339 TObjArray * branches;
341 branches = gAlice->TreeS()->GetListOfBranches() ;
343 cout << "TreeS: " << endl ;
344 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
345 TBranch * branch=(TBranch *) branches->At(ibranch) ;
346 if(strstr(branch->GetName(),"PHOS") )
347 cout << " " << branch->GetName() << " " << branch->GetTitle() << endl ;
351 branches = gAlice->TreeD()->GetListOfBranches() ;
353 cout << "TreeD: " << endl ;
354 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
355 TBranch * branch=(TBranch *) branches->At(ibranch) ;
356 if(strstr(branch->GetName(),"PHOS") )
357 cout << " " << branch->GetName() << " " << branch->GetTitle() << endl ;
362 branches = gAlice->TreeR()->GetListOfBranches() ;
364 cout << "TreeR: " << endl ;
365 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
366 TBranch * branch=(TBranch *) branches->At(ibranch) ;
367 if(strstr(branch->GetName(),"PHOS") )
368 cout << " " << branch->GetName() << " " << branch->GetTitle() << endl ;
374 //____________________________________________________________________________
375 void AliPHOSAnalyze::InvariantMass(const char* branchTitle)
377 // Calculates Real and Mixed invariant mass distributions
378 fObjGetter = AliPHOSIndexToObject::GetInstance(ffileName.Data(),"PHOSRP",branchTitle) ;
380 Int_t nMixedEvents = 4 ; //# of events used for calculation of 'mixed' distribution
383 TFile * mfile = new TFile("invmass.root","update");
385 //========== Reading /Booking Histograms
387 hRealEM = (TH2D*) mfile->Get("hRealEM") ;
389 hRealEM = new TH2D("hRealEM", "Real for EM particles", 250,0.,1.,40,0.,4.) ;
390 TH2D * hRealPhot = 0 ;
392 hRealPhot = (TH2D*)mfile->Get("hRealPhot");
394 hRealPhot = new TH2D("hRealPhot", "Real for kPhoton particles", 250,0.,1.,40,0.,4.) ;
396 TH2D * hMixedEM = 0 ;
397 hMixedEM = (TH2D*) mfile->Get("hMixedEM") ;
399 hMixedEM = new TH2D("hMixedEM", "Mixed for EM particles", 250,0.,1.,40,0.,4.) ;
401 TH2D * hMixedPhot = 0 ;
402 hMixedPhot = (TH2D*) mfile->Get("hMixedPhot") ;
404 hMixedPhot = new TH2D("hMixedPhot","Mixed for kPhoton particles",250,0.,1.,40,0.,4.) ;
407 //reading event and copyng it to TConesArray of all photons
409 TClonesArray * allRecParticleList = new TClonesArray("AliPHOSRecParticle", 1000) ;
410 Int_t * nRecParticles = new Int_t[nMixedEvents] ; // to mark boundaries of each event in the total list
411 for(Int_t index = 0; index < nMixedEvents; index ++)
412 nRecParticles[index] = 0 ;
413 Int_t iRecPhot = 0 ; // number of EM particles in total list
415 //scan over all events
417 for(event = 0; event < fObjGetter->GetMaxEvent(); event++ ){
419 fObjGetter->GetEvent(event);
421 //copy EM RecParticles to the "total" list
422 AliPHOSRecParticle * recParticle ;
424 for(iRecParticle = 0; iRecParticle < fObjGetter->GimeNRecParticles() ;iRecParticle++ )
426 recParticle = fObjGetter->GimeRecParticle(iRecParticle) ;
427 if((recParticle->GetType() == AliPHOSFastRecParticle::kGAMMA)||
428 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEM))
429 new( (*allRecParticleList)[iRecPhot++] ) AliPHOSRecParticle(*recParticle) ;
432 Int_t mevent = event%nMixedEvents ; //event number in the "mixed" cicle
433 nRecParticles[mevent] = iRecPhot-1 ;
435 //check, if it is time to calculate invariant mass?
436 if((mevent == 0) && (event +1 == fObjGetter->GetMaxEvent())){
438 //calculate invariant mass:
440 Int_t nCurEvent = 0 ;
442 for(irp1 = 0; irp1 < allRecParticleList->GetEntries()-1; irp1++){
443 AliPHOSRecParticle * rp1 = (AliPHOSRecParticle *)allRecParticleList->At(irp1) ;
445 for(irp2 = irp1+1; irp2 < allRecParticleList->GetEntries(); irp2++){
446 AliPHOSRecParticle * rp2 = (AliPHOSRecParticle *)allRecParticleList->At(irp2) ;
449 invMass = (rp1->Energy()+rp2->Energy())*(rp1->Energy()+rp2->Energy())-
450 (rp1->Px()+rp2->Px())*(rp1->Px()+rp2->Px())-
451 (rp1->Py()+rp2->Py())*(rp1->Py()+rp2->Py())-
452 (rp1->Pz()+rp2->Pz())*(rp1->Pz()+rp2->Pz()) ;
455 invMass = TMath::Sqrt(invMass);
458 pt = TMath::Sqrt((rp1->Px()+rp2->Px() )*( rp1->Px()+rp2->Px() ) +
459 (rp1->Py()+rp2->Py() )*( rp1->Py()+rp2->Py() ) );
461 if(irp1 > nRecParticles[nCurEvent])
464 if(irp2 <= nRecParticles[nCurEvent]){ //'Real' event
465 hRealEM->Fill(invMass,pt);
466 if((rp1->GetType() == AliPHOSFastRecParticle::kGAMMA)&&
467 (rp2->GetType() == AliPHOSFastRecParticle::kGAMMA) )
468 hRealPhot->Fill(invMass,pt);
471 hMixedEM->Fill(invMass,pt);
472 if((rp1->GetType() == AliPHOSFastRecParticle::kGAMMA)&&
473 (rp2->GetType() == AliPHOSFastRecParticle::kGAMMA) )
474 hMixedPhot->Fill(invMass,pt);
477 } //loop over second rp
478 }//loop over first rp
480 //Make some cleanings
481 for(Int_t index = 0; index < nMixedEvents; index ++)
482 nRecParticles[index] = 0 ;
484 allRecParticleList->Clear() ;
488 delete allRecParticleList ;
493 hRealEM->Write(0,kOverwrite) ;
494 hRealPhot->Write(0,kOverwrite) ;
495 hMixedEM->Write(0,kOverwrite) ;
496 hMixedPhot->Write(0,kOverwrite) ;
501 delete nRecParticles;
505 //____________________________________________________________________________
506 void AliPHOSAnalyze::EnergyResolution(const char * branchTitle)
508 //fills two dimentional histo: energy of primary vs. energy of reconstructed
510 TH2F * hAllEnergy = 0 ; //all reconstructed with primary photon
511 TH2F * hPhotEnergy= 0 ; //kGamma with primary photon
512 TH2F * hEMEnergy = 0 ; //electromagnetic with primary photon
514 //opening file and reading histograms if any
515 TFile * efile = new TFile("energy.root","update");
517 hAllEnergy = (TH2F*)efile->Get("hAllEnergy") ;
519 hAllEnergy = new TH2F("hAllEnergy", "Energy of any RP with primary photon",100, 0., 5., 100, 0., 5.);
521 hPhotEnergy =(TH2F*) efile->Get("hPhotEnergy") ;
523 hPhotEnergy = new TH2F("hPhotEnergy", "Energy of kGAMMA with primary photon",100, 0., 5., 100, 0., 5.);
525 hEMEnergy =(TH2F*) efile->Get("hEMEnergy");
527 hEMEnergy = new TH2F("hEMEnergy", "Energy of EM with primary photon", 100, 0., 5., 100, 0., 5.);
530 fObjGetter = AliPHOSIndexToObject::GetInstance(ffileName.Data(),"PHOSRP",branchTitle) ;
531 fGeom = AliPHOSGeometry::GetInstance(((AliPHOS*)gAlice->GetModule("PHOS"))->GetGeometry()->GetName(),
532 ((AliPHOS*)gAlice->GetModule("PHOS"))->GetGeometry()->GetTitle() );
535 for ( ievent=0; ievent < fObjGetter->GetMaxEvent() ; ievent++){
537 //read the current event
538 fObjGetter->GetEvent(ievent) ;
540 AliPHOSRecParticle * recParticle ;
542 for(iRecParticle = 0; iRecParticle < fObjGetter->GimeNRecParticles() ;iRecParticle++ ){
543 recParticle = fObjGetter->GimeRecParticle(iRecParticle) ;
545 //find the closest primary
546 Int_t moduleNumberRec ;
547 Double_t recX, recZ ;
548 fGeom->ImpactOnEmc(recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
550 Double_t minDistance = 100. ;
551 Int_t closestPrimary = -1 ;
553 //extract list of primaries: it is stored at EMC RecPoints
554 Int_t emcIndex = fObjGetter->GimeTrackSegment(recParticle->GetPHOSTSIndex())->GetEmcIndex() ;
555 Int_t numberofprimaries ;
556 Int_t * listofprimaries = fObjGetter->GimeEmcRecPoint(emcIndex)->GetPrimaries(numberofprimaries) ;
559 TParticle * primary ;
560 Double_t distance = minDistance ;
563 Double_t dZmin = 0. ;
564 for ( index = 0 ; index < numberofprimaries ; index++){
565 primary = fObjGetter->GimePrimary(listofprimaries[index]) ;
567 Double_t primX, primZ ;
568 fGeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
569 if(moduleNumberRec == moduleNumber) {
572 distance = TMath::Sqrt(dX*dX + dZ*dZ) ;
573 if(minDistance > distance) {
574 minDistance = distance ;
577 closestPrimary = listofprimaries[index] ;
582 //if found primary, fill histograms
583 if(closestPrimary >=0 ){
584 TParticle * primary = fObjGetter->GimePrimary(closestPrimary) ;
585 if(primary->GetPdgCode() == 22){
586 hAllEnergy->Fill(primary->Energy(), recParticle->Energy()) ;
587 if(recParticle->GetType() == AliPHOSFastRecParticle::kGAMMA){
588 hPhotEnergy->Fill(primary->Energy(), recParticle->Energy() ) ;
589 hEMEnergy->Fill(primary->Energy(), recParticle->Energy() ) ;
592 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEM)
593 hEMEnergy->Fill(primary->Energy(), recParticle->Energy() ) ;
599 //write filled histograms
601 hAllEnergy->Write(0,kOverwrite) ;
602 hPhotEnergy->Write(0,kOverwrite) ;
603 hEMEnergy->Write(0,kOverwrite) ;
609 //____________________________________________________________________________
610 void AliPHOSAnalyze::PositionResolution(const char * branchTitle)
612 //fills two dimentional histo: energy vs. primary - reconstructed distance
616 TH2F * hAllPosition = 0; // Position of any RP with primary photon
617 TH2F * hPhotPosition = 0; // Position of kGAMMA with primary photon
618 TH2F * hEMPosition = 0; // Position of EM with primary photon
620 TH1F * hAllPositionX = 0; // X-Position Resolution of photons with photon primary
621 TH1F * hAllPositionZ = 0; // Z-Position Resolution of photons with photon primary
624 //opening file and reading histograms if any
625 TFile * pfile = new TFile("position.root","update");
627 hAllPosition = (TH2F*)pfile->Get("hAllPosition");
628 if(hAllPosition == 0)
629 hAllPosition = new TH2F("hAllPosition",
630 "Position of any RP with primary photon",100, 0., 5., 100, 0., 5.);
631 hPhotPosition= (TH2F*)pfile->Get("hPhotPosition");
632 if(hPhotPosition == 0)
633 hPhotPosition = new TH2F("hPhotPosition",
634 "Position of kGAMMA with primary photon",100, 0., 5., 100, 0., 5.);
635 hEMPosition= (TH2F*)pfile->Get("hEMPosition") ;
637 hEMPosition = new TH2F("hEMPosition",
638 "Position of EM with primary photon", 100, 0., 5., 100, 0., 5.);
639 hAllPositionX = (TH1F*)pfile->Get("hAllPositionX") ;
640 if(hAllPositionX == 0)
641 hAllPositionX = new TH1F("hAllPositionX",
642 "Delta X of any RP with primary photon",100, -2., 2.);
643 hAllPositionZ =(TH1F*) pfile->Get("hAllPositionZ") ;
644 if(hAllPositionZ == 0)
645 hAllPositionZ = new TH1F("hAllPositionZ",
646 "Delta X of any RP with primary photon",100, -2., 2.);
649 fObjGetter = AliPHOSIndexToObject::GetInstance(ffileName.Data(),"PHOSRP",branchTitle) ;
650 fGeom = AliPHOSGeometry::GetInstance(((AliPHOS*)gAlice->GetModule("PHOS"))->GetGeometry()->GetName(),
651 ((AliPHOS*)gAlice->GetModule("PHOS"))->GetGeometry()->GetTitle() );
654 for ( ievent=0; ievent < fObjGetter->GetMaxEvent() ; ievent++){
656 //read the current event
657 fObjGetter->GetEvent(ievent) ;
659 AliPHOSRecParticle * recParticle ;
661 for(iRecParticle = 0; iRecParticle < fObjGetter->GimeNRecParticles() ;iRecParticle++ ){
662 recParticle = fObjGetter->GimeRecParticle(iRecParticle) ;
664 //find the closest primary
665 Int_t moduleNumberRec ;
666 Double_t recX, recZ ;
667 fGeom->ImpactOnEmc(recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
669 Double_t minDistance = 100. ;
670 Int_t closestPrimary = -1 ;
672 //extract list of primaries: it is stored at EMC RecPoints
673 Int_t emcIndex = fObjGetter->GimeTrackSegment(recParticle->GetPHOSTSIndex())->GetEmcIndex() ;
674 Int_t numberofprimaries ;
675 Int_t * listofprimaries = fObjGetter->GimeEmcRecPoint(emcIndex)->GetPrimaries(numberofprimaries) ;
678 TParticle * primary ;
679 Double_t distance = minDistance ;
680 Double_t dX = 1000; // incredible number
681 Double_t dZ = 1000; // for the case if no primary will be found
683 Double_t dZmin = 0. ;
684 for ( index = 0 ; index < numberofprimaries ; index++){
685 primary = fObjGetter->GimePrimary(listofprimaries[index]) ;
687 Double_t primX, primZ ;
688 fGeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
689 if(moduleNumberRec == moduleNumber) {
692 distance = TMath::Sqrt(dX*dX + dZ*dZ) ;
693 if(minDistance > distance) {
694 minDistance = distance ;
697 closestPrimary = listofprimaries[index] ;
702 //if found primary, fill histograms
703 if(closestPrimary >=0 ){
704 TParticle * primary = fObjGetter->GimePrimary(closestPrimary) ;
705 if(primary->GetPdgCode() == 22){
706 hAllPosition->Fill(primary->Energy(), minDistance) ;
707 hAllPositionX->Fill(primary->Energy(), dX) ;
708 hAllPositionZ->Fill(primary->Energy(), dZ) ;
709 if(recParticle->GetType() == AliPHOSFastRecParticle::kGAMMA){
710 hPhotPosition->Fill(primary->Energy(), minDistance ) ;
711 hEMPosition->Fill(primary->Energy(), minDistance ) ;
714 if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEM)
715 hEMPosition->Fill(primary->Energy(), minDistance ) ;
721 //Write output histgrams
723 hAllPosition->Write(0,kOverwrite) ;
724 hAllPositionX->Write(0,kOverwrite) ;
725 hAllPositionZ->Write(0,kOverwrite) ;
726 hPhotPosition->Write(0,kOverwrite) ;
727 hEMPosition->Write(0,kOverwrite) ;
734 //____________________________________________________________________________
735 void AliPHOSAnalyze::Contamination(const char* RecPointsTitle){
736 // fills spectra of primary photons and several kinds of
737 // reconstructed particles, so that analyzing them one can
738 // estimate conatmination, efficiency of registration etc.
740 //define several general histograms
741 TH1F * hPrimary = 0; //spectrum (P_t distribution) of primary photons
742 TH1F * hAllRP = 0; //spectrum of all RecParticles in PHOS
743 TH1F * hPhot = 0; //spectrum of kGAMMA RecParticles
744 TH1F * hPPSD = 0; //spectrum of all RecParticles with PPSD signal
745 TH1F * hShape = 0; //spectrum of all EM RecParticles
746 TH1F * hVeto = 0; //spectrum of all neutral RecParticles
748 //Now separate histograms in accoradance with primary
750 TH1F * hPhotReg = 0; //Registeres as photon
751 TH1F * hPhotEM = 0; //Registered as EM
752 TH1F * hPhotPPSD= 0; //Registered as RecParticle with PPSD signal
755 TH1F * hNReg = 0; //Registeres as photon
756 TH1F * hNEM = 0; //Registered as EM
757 TH1F * hNPPSD= 0; //Registered as RecParticle with PPSD signal
760 TH1F * hNBarReg = 0; //Registeres as photon
761 TH1F * hNBarEM = 0; //Registered as EM
762 TH1F * hNBarPPSD= 0; //Registered as RecParticle with PPSD signal
764 //primary - charged hadron (pBar excluded)
765 TH1F * hChargedReg = 0; //Registeres as photon
766 TH1F * hChargedEM = 0; //Registered as EM
767 TH1F * hChargedPPSD= 0; //Registered as RecParticle with PPSD signal
770 TH1F * hPbarReg = 0; //Registeres as photon
771 TH1F * hPbarEM = 0; //Registered as EM
772 TH1F * hPbarPPSD= 0; //Registered as RecParticle with PPSD signal
775 //Reading histograms from the file
776 TFile * cfile = new TFile("contamination.root","update") ;
778 //read general histograms
779 hPrimary = (TH1F*) cfile->Get("hPrimary") ;
781 hPrimary= new TH1F("hPrimary", "Primary photon spectrum", 100, 0., 5.);
782 hAllRP = (TH1F*)cfile->Get("hAllRP") ;
784 hAllRP = new TH1F("hAllRP","All Reconstructed particles", 100, 0., 5.);
785 hPhot = (TH1F*)cfile->Get("hPhot") ;
787 hPhot = new TH1F("hPhot","All kGAMMA RecParticles",100, 0., 5.);
788 hPPSD = (TH1F*)cfile->Get("hPPSD") ;
790 hPPSD = new TH1F("hPPSD", "All PPSD Recparticles", 100, 0., 5.);
791 hShape = (TH1F*) cfile->Get("hShape") ;
793 hShape = new TH1F("hShape","All particles with EM shower",100, 0., 5.);
794 hVeto= (TH1F*)cfile->Get("hVeto") ;
796 hVeto = new TH1F("hVeto", "All uncharged particles", 100, 0., 5.);
800 hPhotReg = (TH1F*)cfile->Get("hPhotReg");
802 hPhotReg = new TH1F("hPhotReg","Photon registered as photon",100, 0., 5.);
803 hPhotEM =(TH1F*)cfile->Get("hPhotEM");
805 hPhotEM = new TH1F("hPhotEM", "Photon registered as EM", 100, 0., 5.);
806 hPhotPPSD= (TH1F*)cfile->Get("hPhotPPSD");
808 hPhotPPSD = new TH1F("hPhotPPSD","Photon registered as PPSD", 100, 0., 5.);
811 hNReg = (TH1F*)cfile->Get("hNReg");
813 hNReg = new TH1F("hNReg", "N registered as photon", 100, 0., 5.);
814 hNEM = (TH1F*)cfile->Get("hNEM");
816 hNEM = new TH1F("hNEM", "N registered as EM", 100, 0., 5.);
817 hNPPSD=(TH1F*)cfile->Get("hNPPSD");
819 hNPPSD = new TH1F("hNPPSD","N registered as PPSD", 100, 0., 5.);
822 hNBarReg =(TH1F*)cfile->Get("hNBarReg");
824 hNBarReg = new TH1F("hNBarReg", "NBar registered as photon", 100, 0., 5.);
825 hNBarEM =(TH1F*)cfile->Get("hNBarEM");
827 hNBarEM = new TH1F("hNBarEM", "NBar registered as EM", 100, 0., 5.);
828 hNBarPPSD=(TH1F*)cfile->Get("hNBarPPSD");
830 hNBarPPSD = new TH1F("hNBarPPSD","NBar registered as PPSD", 100, 0., 5.);
832 //primary - charged hadron (pBar excluded)
833 hChargedReg = (TH1F*)cfile->Get("hChargedReg");
835 hChargedReg= new TH1F("hChargedReg", "Charged hadron registered as photon",100, 0., 5.);
836 hChargedEM = (TH1F*)cfile->Get("hChargedEM");
838 hChargedEM= new TH1F("hChargedEM","Charged registered as EM",100, 0., 5.);
839 hChargedPPSD= (TH1F*)cfile->Get("hChargedPPSD");
841 hChargedPPSD= new TH1F("hChargedPPSD","Charged registered as PPSD",100, 0., 5.);
844 hPbarReg = (TH1F*)cfile->Get("hPbarReg");
846 hPbarReg= new TH1F("hPbarReg", "pBar registered as photon",100, 0., 5.);
847 hPbarEM = (TH1F*)cfile->Get("hPbarEM");
849 hPbarEM= new TH1F("hPbarEM","Pbar registered as EM",100, 0., 5.);
850 hPbarPPSD= (TH1F*)cfile->Get("hPbarPPSD");
852 hPbarPPSD= new TH1F("hPbarPPSD","Pbar as PPSD",100, 0., 5.);
855 //Now make some initializations
857 Int_t counter[8][5] ; //# of registered particles
859 for(i1 = 0; i1<8; i1++)
860 for(i2 = 0; i2<5; i2++)
861 counter[i1][i2] = 0 ;
865 fObjGetter = AliPHOSIndexToObject::GetInstance(ffileName.Data(),"PHOSRP",RecPointsTitle) ;
866 fGeom = AliPHOSGeometry::GetInstance(((AliPHOS*)gAlice->GetModule("PHOS"))->GetGeometry()->GetName(),
867 ((AliPHOS*)gAlice->GetModule("PHOS"))->GetGeometry()->GetTitle() );
870 for ( ievent=0; ievent < fObjGetter->GetMaxEvent() ; ievent++){
872 fObjGetter->GetEvent(ievent) ;
874 //=========== Make spectrum of the primary photons
875 TParticle * primary ;
877 for( iPrimary = 0 ; iPrimary < fObjGetter->GimeNPrimaries() ; iPrimary++){
878 primary = fObjGetter->GimePrimary(iPrimary) ;
879 Int_t primaryType = primary->GetPdgCode() ;
880 if( primaryType == 22 ) {
881 //check, if photons folls onto PHOS
883 Double_t primX, primZ ;
884 fGeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
886 hPrimary->Fill(primary->Energy()) ;
891 //========== Now scan over RecParticles
892 AliPHOSRecParticle * recParticle ;
894 for(iRecParticle = 0; iRecParticle < fObjGetter->GimeNRecParticles(); iRecParticle++ ){
895 recParticle = fObjGetter->GimeRecParticle(iRecParticle) ;
896 //fill histo spectrum of all RecParticles
897 hAllRP->Fill(CorrectedEnergy(recParticle->Energy())) ;
899 //==========find the closest primary
900 Int_t moduleNumberRec ;
901 Double_t recX, recZ ;
902 fGeom->ImpactOnEmc(recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
904 Double_t minDistance = 100. ;
905 Int_t closestPrimary = -1 ;
907 //extract list of primaries: it is stored at EMC RecPoints
908 Int_t emcIndex = fObjGetter->GimeTrackSegment(recParticle->GetPHOSTSIndex())->GetEmcIndex() ;
909 Int_t numberofprimaries ;
910 Int_t * listofprimaries = fObjGetter->GimeEmcRecPoint(emcIndex)->GetPrimaries(numberofprimaries) ;
912 TParticle * primary ;
913 Double_t distance = minDistance ;
916 Double_t dZmin = 0. ;
917 for ( index = 0 ; index < numberofprimaries ; index++){
918 primary = fObjGetter->GimePrimary(listofprimaries[index]) ;
920 Double_t primX, primZ ;
921 fGeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
922 if(moduleNumberRec == moduleNumber) {
925 distance = TMath::Sqrt(dX*dX + dZ*dZ) ;
926 if(minDistance > distance) {
927 minDistance = distance ;
930 closestPrimary = listofprimaries[index] ;
935 //===========define the "type" of closest primary
936 if(closestPrimary >=0 ){
937 Int_t primaryCode = -1;
938 TParticle * primary = fObjGetter->GimePrimary(closestPrimary) ;
939 Int_t primaryType = primary->GetPdgCode() ;
940 if(primaryType == 22) // photon ?
943 if(primaryType == 2112) // neutron
946 if(primaryType == -2112) // Anti neutron
949 if(primaryType == -2122) //Anti proton
952 if(fObjGetter->GimePrimary(closestPrimary)->GetPDG()->Charge())
955 //==========Now look at the type of RecParticle
956 Float_t energy = CorrectedEnergy(recParticle->Energy()) ;
957 if(recParticle->GetType() == AliPHOSFastRecParticle::kGAMMA){
958 hPhot->Fill(energy ) ;
961 hPhotReg->Fill(energy ) ;
964 hNReg->Fill(energy ) ;
967 hNBarReg->Fill(energy ) ;
970 hChargedReg->Fill(energy ) ;
973 hPbarReg->Fill(energy ) ;
979 if((recParticle->GetType() == AliPHOSFastRecParticle::kGAMMA)||
980 (recParticle->GetType() == AliPHOSFastRecParticle::kGAMMAHA)){ //with PPSD signal
981 hPPSD->Fill(energy ) ;
984 hPhotPPSD->Fill(energy ) ;
987 hNPPSD->Fill(energy ) ;
990 hNBarPPSD->Fill(energy ) ;
993 hChargedPPSD->Fill(energy ) ;
996 hPbarPPSD->Fill(energy ) ;
1002 if((recParticle->GetType() == AliPHOSFastRecParticle::kGAMMA)||
1003 (recParticle->GetType() == AliPHOSFastRecParticle::kELECTRON)||
1004 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEM)||
1005 (recParticle->GetType() == AliPHOSFastRecParticle::kABSURDEM) ){ //with EM shower
1006 hShape->Fill(energy ) ;
1007 switch(primaryCode){
1009 hPhotEM->Fill(energy ) ;
1012 hNEM->Fill(energy ) ;
1015 hNBarEM->Fill(energy ) ;
1018 hChargedEM->Fill(energy ) ;
1021 hPbarEM->Fill(energy ) ;
1028 if((recParticle->GetType() == AliPHOSFastRecParticle::kGAMMA)||
1029 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALHA) ||
1030 (recParticle->GetType() == AliPHOSFastRecParticle::kGAMMAHA) ||
1031 (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEM) ) //nuetral
1032 hVeto->Fill(energy ) ;
1034 //fill number of primaries identified as ...
1035 if(primaryCode >= 0) // Primary code defined
1036 counter[recParticle->GetType()][primaryCode]++ ;
1040 } // no closest primary found
1044 //=================== SaveHistograms
1046 hPrimary->Write(0,kOverwrite);
1047 hAllRP->Write(0,kOverwrite);
1048 hPhot->Write(0,kOverwrite);
1049 hPPSD->Write(0,kOverwrite);
1050 hShape->Write(0,kOverwrite);
1051 hVeto->Write(0,kOverwrite);
1052 hPhotReg->Write(0,kOverwrite);
1053 hPhotEM->Write(0,kOverwrite);
1054 hPhotPPSD->Write(0,kOverwrite);
1055 hNReg ->Write(0,kOverwrite);
1056 hNEM ->Write(0,kOverwrite);
1057 hNPPSD->Write(0,kOverwrite);
1058 hNBarReg ->Write(0,kOverwrite);
1059 hNBarEM ->Write(0,kOverwrite);
1060 hNBarPPSD->Write(0,kOverwrite);
1061 hChargedReg ->Write(0,kOverwrite);
1062 hChargedEM ->Write(0,kOverwrite);
1063 hChargedPPSD->Write(0,kOverwrite);
1064 hPbarReg ->Write(0,kOverwrite);
1065 hPbarEM ->Write(0,kOverwrite);
1066 hPbarPPSD->Write(0,kOverwrite);
1068 cfile->Write(0,kOverwrite);
1075 cout << "Resolutions: Analyzed " << fObjGetter->GetMaxEvent() << " event(s)" << endl ;
1078 cout << " Primary: Photon Neutron Antineutron Charged hadron AntiProton" << endl ;
1079 cout << "--------------------------------------------------------------------------------" << endl ;
1081 << setw(8) << counter[2][0] << setw(9) << counter[2][1] << setw(13) << counter[2][2]
1082 << setw(15)<< counter[2][3] << setw(13) << counter[2][4] << endl ;
1083 cout << " kGAMMAHA: "
1084 << setw(8) << counter[3][0] << setw(9) << counter[3][1] << setw(13) << counter[3][2]
1085 << setw(15)<< counter[3][3] << setw(13) << counter[3][4] << endl ;
1086 cout << " kNEUTRALEM: "
1087 << setw(8) << counter[0][0] << setw(9) << counter[0][1] << setw(13) << counter[0][2]
1088 << setw(15)<< counter[0][3] << setw(13) << counter[0][4] << endl ;
1089 cout << " kNEUTRALHA: "
1090 << setw(8) << counter[1][0] << setw(9) << counter[1][1] << setw(13) << counter[1][2]
1091 << setw(15)<< counter[1][3] << setw(13) << counter[1][4] << endl ;
1092 cout << " kABSURDEM: "
1093 << setw(8) << counter[4][0] << setw(9) << counter[4][1] << setw(13) << counter[4][2]
1094 << setw(15)<< counter[4][3] << setw(13) << counter[4][4] << endl ;
1095 cout << " kABSURDHA: "
1096 << setw(8) << counter[5][0] << setw(9) << counter[5][1] << setw(13) << counter[5][2]
1097 << setw(15)<< counter[5][3] << setw(13) << counter[5][4] << endl ;
1098 cout << " kELECTRON: "
1099 << setw(8) << counter[6][0] << setw(9) << counter[6][1] << setw(13) << counter[6][2]
1100 << setw(15)<< counter[6][3] << setw(13) << counter[6][4] << endl ;
1101 cout << " kCHARGEDHA: "
1102 << setw(8) << counter[7][0] << setw(9) << counter[7][1] << setw(13) << counter[7][2]
1103 << setw(15)<< counter[7][3] << setw(13) << counter[7][4] << endl ;
1104 cout << "--------------------------------------------------------------------------------" << endl ;
1106 Int_t totalInd = 0 ;
1107 for(i1 = 0; i1<8; i1++)
1108 for(i2 = 0; i2<5; i2++)
1109 totalInd+=counter[i1][i2] ;
1110 cout << "Indentified particles: " << totalInd << endl ;