]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSAnalyze.cxx
Updated PID with photon/pi0 recognition by M2x at high-pt
[u/mrichter/AliRoot.git] / PHOS / AliPHOSAnalyze.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id$ */
17
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.
31 ////    User Case:
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.
59 //*--
60 //*-- Author: Dmitri Peressounko (SUBATECH & RRC Kurchatov Institute)
61 //////////////////////////////////////////////////////////////////////////////
62
63
64 // --- ROOT system ---
65
66 #include "TFile.h"
67 #include "TH1.h"
68 #include "TH2.h"
69 #include "TH2.h"
70 #include "TParticle.h"
71 #include "TClonesArray.h"
72 #include "TTree.h"
73 #include "TMath.h"
74 #include "TROOT.h"
75 #include "TFolder.h"
76
77 // --- Standard library ---
78
79 // --- AliRoot header files ---
80
81 #include "AliRun.h"
82 #include "AliPHOSv1.h"
83 #include "AliPHOSAnalyze.h"
84 #include "AliPHOSDigit.h"
85 #include "AliPHOSSDigitizer.h"
86 #include "AliPHOSTrackSegment.h"
87 #include "AliPHOSRecParticle.h"
88 #include "AliPHOSCpvRecPoint.h"
89 #include "AliPHOSGetter.h"
90
91
92 ClassImp(AliPHOSAnalyze)
93
94 //____________________________________________________________________________
95   AliPHOSAnalyze::AliPHOSAnalyze()
96 {
97   // default ctor (useless)
98   fCorrection = 1.2 ;  //Value calculated for default parameters of reconstruction  
99 }
100
101 //____________________________________________________________________________
102 AliPHOSAnalyze::AliPHOSAnalyze(Text_t * fileName)
103 {
104   // ctor: analyze events from root file "name"
105   ffileName = fileName ;
106   fCorrection = 1.05 ;  //Value calculated for default parameters of reconstruction   
107 }
108
109 //____________________________________________________________________________
110 AliPHOSAnalyze::AliPHOSAnalyze(const AliPHOSAnalyze & ana)
111 {
112   // copy ctor
113   ( (AliPHOSAnalyze &)ana ).Copy(*this) ;
114 }
115
116 //____________________________________________________________________________
117 AliPHOSAnalyze::~AliPHOSAnalyze()
118 {
119   // dtor
120
121 }
122 //____________________________________________________________________________
123 void AliPHOSAnalyze::DrawRecon(Int_t Nevent,Int_t Nmod,const char * what,const char* branchTitle){
124   //Draws pimary particles and reconstructed 
125   //digits, RecPoints, RecPartices etc 
126   //for event Nevent in the module Nmod.
127
128   //========== Create ObjectGetter
129   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(ffileName.Data(),branchTitle) ;
130   if(Nevent >= gAlice->TreeE()->GetEntries() ) {
131     Error("DrawRecon", "There is no event %d only %d events available", Nevent, gAlice->TreeE()->GetEntries() ) ;
132     return ;
133   }
134   const AliPHOSGeometry * phosgeom = gime->PHOSGeometry() ; 
135   gime->Event(Nevent);
136
137   Int_t nx = phosgeom->GetNPhi() ;
138   Int_t nz = phosgeom->GetNZ() ;
139   Float_t * cri= phosgeom->GetEMCAGeometry()->GetCrystalHalfSize() ;
140   Float_t x = nx*cri[0] ;
141   Float_t z = nz*cri[2] ;
142   Int_t nxCPV = (Int_t) (nx*phosgeom->GetPadSizePhi()/(2.*cri[0])) ;
143   Int_t nzCPV = (Int_t) (nz*phosgeom->GetPadSizeZ()/(2.*cri[2])) ;
144   
145   TH2F * emcDigits = (TH2F*) gROOT->FindObject("emcDigits") ;
146   if(emcDigits)
147     emcDigits->Delete() ;
148   emcDigits = new TH2F("emcDigits","EMC digits",  nx,-x,x,nz,-z,z);
149   TH2F * emcSdigits =(TH2F*) gROOT->FindObject("emcSdigits") ;  
150   if(emcSdigits)
151     emcSdigits->Delete() ;
152   emcSdigits = new TH2F("emcSdigits","EMC sdigits", nx,-x,x,nz,-z,z);
153   TH2F * emcRecPoints = (TH2F*)gROOT->FindObject("emcRecPoints") ; 
154   if(emcRecPoints)
155     emcRecPoints->Delete() ;
156   emcRecPoints = new TH2F("emcRecPoints","EMC RecPoints",nx,-x,x,nz,-z,z);
157   TH2F * cpvSdigits =(TH2F*) gROOT->FindObject("cpvSdigits") ;
158   if(cpvSdigits)
159     cpvSdigits->Delete() ;
160   cpvSdigits = new TH2F("cpvSdigits","CPV sdigits", nx,-x,x,nz,-z,z);
161   TH2F * cpvDigits = (TH2F*)gROOT->FindObject("cpvDigits") ;
162   if(cpvDigits)
163     cpvDigits->Delete() ;
164   cpvDigits = new TH2F("cpvDigits","CPV digits",   nxCPV,-x,x,nzCPV,-z,z) ;
165   TH2F * cpvRecPoints= (TH2F*)gROOT->FindObject("cpvRecPoints") ; 
166   if(cpvRecPoints)
167     cpvRecPoints->Delete() ;
168   cpvRecPoints = new TH2F("cpvRecPoints","CPV RecPoints",    nxCPV,-x,x,nzCPV,-z,z) ;
169
170   TH2F * phot = (TH2F*)gROOT->FindObject("phot") ;
171   if(phot)
172     phot->Delete() ;
173   phot = new TH2F("phot","Primary Photon",  nx,-x,x,nz,-z,z);
174   TH2F * recPhot = (TH2F*)gROOT->FindObject("recPhot") ; 
175   if(recPhot)
176     recPhot->Delete() ;
177   recPhot = new TH2F("recPhot","RecParticles with primary Photon",nx,-x,x,nz,-z,z);
178   
179   if(strstr(what,"P") || strstr(what,"p")){ 
180   //Plot Primary Particles
181   const TParticle * primary ;
182   Int_t iPrimary ;
183   for ( iPrimary = 0 ; iPrimary < gime->NPrimaries() ; iPrimary++)
184     {
185       primary = gime->Primary(iPrimary) ;
186       Int_t primaryType = primary->GetPdgCode() ;
187 //       if( (primaryType == 211)||(primaryType == -211)||(primaryType == 2212)||(primaryType == -2212)
188 //        ||(primaryType == 11)||(primaryType == -11) ) {
189 //         Int_t moduleNumber ;
190 //         Double_t primX, primZ ;
191 //         phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
192 //         if(moduleNumber==Nmod)
193 //           charg->Fill(primZ,primX,primary->Energy()) ;
194 //       }
195       if( primaryType == 22 ) {
196         Int_t moduleNumber ;
197         Double_t primX, primZ ;
198         phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
199         if(moduleNumber==Nmod) 
200           phot->Fill(primZ,primX,primary->Energy()) ;
201       }
202 //       else{
203 //         if( primaryType == -2112 ) {
204 //           Int_t moduleNumber ;
205 //           Double_t primX, primZ ;
206 //           phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
207 //           if(moduleNumber==Nmod)
208 //             nbar->Fill(primZ,primX,primary->Energy()) ;
209 //         }
210 //       }
211     }  
212     phot->Draw("box") ;
213   }
214
215   
216   if(strstr(what,"S") || strstr(what,"s")){
217   Int_t iSDigit ;
218   AliPHOSDigit * sdigit ;
219   const TClonesArray * sdigits = gime->SDigits() ;
220   Int_t nsdig[5] = {0,0,0,0,0} ;
221   if(sdigits){
222     for(iSDigit = 0; iSDigit < sdigits->GetEntriesFast() ; iSDigit++)
223       {
224         sdigit = (AliPHOSDigit *) sdigits->At(iSDigit) ;
225         Int_t relid[4];
226         phosgeom->AbsToRelNumbering(sdigit->GetId(), relid) ;
227         Float_t x,z ;
228         phosgeom->RelPosInModule(relid,x,z) ;
229         Float_t e = gime->SDigitizer()->Calibrate(sdigit->GetAmp()) ;
230         nsdig[relid[0]-1]++ ;
231         if(relid[0]==Nmod){
232           if(relid[1]==0)  //EMC
233             emcSdigits->Fill(x,z,e) ;
234           if( relid[1]!=0 )
235             cpvSdigits->Fill(x,z,e) ;
236         }
237       }
238   }
239   TString message ; 
240   message  = "Number of EMC + CPV SDigits per module: \n" ;
241   message += "%d %d %d %d %d\n"; 
242   Info("DrawRecon", message.Data(), nsdig[0], nsdig[1], nsdig[2], nsdig[3], nsdig[4] ) ;
243   emcSdigits->Draw("box") ;
244   }
245
246   if(strstr(what,"D") || strstr(what,"d")){
247   //Plot digits
248   Int_t iDigit ;
249   AliPHOSDigit * digit ;
250   const TClonesArray * digits = gime->Digits(); 
251   if(digits) {
252     for(iDigit = 0; iDigit < digits->GetEntriesFast(); iDigit++)
253       {
254         digit = (AliPHOSDigit *) digits->At(iDigit) ;
255         Int_t relid[4];
256         phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
257         Float_t x,z ;
258         phosgeom->RelPosInModule(relid,x,z) ;
259         Float_t e = gime->SDigitizer()->Calibrate(digit->GetAmp()) ;
260         if(relid[0]==Nmod){
261           if(relid[1]==0)  //EMC
262             emcDigits->Fill(x,z,e) ;
263           if( relid[1]!=0 )
264             cpvDigits->Fill(x,z,e) ;
265         }
266       }
267   }
268   emcDigits->SetLineColor(5) ;
269   emcDigits->Draw("boxsame") ;
270   }
271   
272   if(strstr(what,"R") || strstr(what,"r")){
273   //Plot RecPoints
274   Int_t irecp ;
275   TVector3 pos ;
276   TObjArray * emcrp = gime->EmcRecPoints() ;
277   if(emcrp) {
278     for(irecp = 0; irecp < emcrp->GetEntriesFast() ; irecp ++){
279       AliPHOSEmcRecPoint * emc = (AliPHOSEmcRecPoint *) emcrp->At(irecp) ;
280       if(emc->GetPHOSMod()==Nmod){
281         emc->GetLocalPosition(pos) ;
282         emcRecPoints->Fill(pos.X(),pos.Z(),emc->GetEnergy());
283       }
284     }
285   }
286   
287   TObjArray * cpvrp = gime->CpvRecPoints() ;
288   if(cpvrp) {
289     for(irecp = 0; irecp < cpvrp->GetEntriesFast() ; irecp ++){
290       AliPHOSRecPoint * cpv = (AliPHOSCpvRecPoint *) cpvrp->At(irecp) ;
291       if(cpv->GetPHOSMod()==Nmod){
292         cpv->GetLocalPosition(pos) ;
293         cpvRecPoints->Fill(pos.X(),pos.Z(),cpv->GetEnergy());
294       }
295     }
296   }
297     
298   //Plot RecParticles
299   AliPHOSRecParticle * recParticle ;
300   Int_t iRecParticle ;
301   TClonesArray * rp = gime->RecParticles() ;
302   TClonesArray * ts = gime->TrackSegments() ;
303   if(rp && ts && emcrp) {
304     for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast() ; iRecParticle++ )
305       {
306         recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
307         Int_t moduleNumberRec ;
308         Double_t recX, recZ ;
309         phosgeom->ImpactOnEmc(recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
310         if(moduleNumberRec == Nmod){
311           
312           Double_t minDistance = 5. ;
313           Int_t closestPrimary = -1 ;   
314
315           //extract list of primaries: it is stored at EMC RecPoints
316           Int_t emcIndex = ((AliPHOSTrackSegment *) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
317           Int_t numberofprimaries ;
318           Int_t * listofprimaries  = ((AliPHOSRecPoint*) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries)  ;
319           Int_t index ;
320           const TParticle * primary ;
321           Double_t distance = minDistance ;
322           
323           for ( index = 0 ; index < numberofprimaries ; index++){
324             primary = gime->Primary(listofprimaries[index]) ;
325             Int_t moduleNumber ;
326             Double_t primX, primZ ;
327             phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
328             if(moduleNumberRec == moduleNumber)
329               distance = TMath::Sqrt((recX-primX)*(recX-primX)+(recZ-primZ)*(recZ-primZ) ) ;
330             if(minDistance > distance)
331               {
332                 minDistance = distance ;
333                 closestPrimary = listofprimaries[index] ;
334               }
335           }
336           
337           if(closestPrimary >=0 ){
338             
339             Int_t primaryType = gime->Primary(closestPrimary)->GetPdgCode() ;
340             
341             if(primaryType==22)
342               recPhot->Fill(recZ,recX,recParticle->Energy()) ;
343 //          else
344 //            if(primaryType==-2112)
345 //              recNbar->Fill(recZ,recX,recParticle->Energy()) ; 
346           }
347         }
348       }
349
350   }
351   //Plot made histograms
352   emcRecPoints->SetLineColor(2) ;
353   emcRecPoints->Draw("boxsame") ;
354   cpvSdigits->SetLineColor(1) ;
355   cpvSdigits->Draw("boxsame") ;
356   }  
357
358
359   
360 }
361 //____________________________________________________________________________
362 void AliPHOSAnalyze::Ls(){
363   //lists branches and titles of PHOS-related branches of TreeR, TreeD, TreeS
364   
365   AliPHOSGetter::GetInstance(ffileName.Data()) ;
366
367   Int_t ibranch;
368   TObjArray * branches; 
369   
370   branches = gAlice->TreeS()->GetListOfBranches() ;
371  
372   TString message ; 
373   message  = "TreeS:\n" ;
374   for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
375     TBranch * branch=(TBranch *) branches->At(ibranch) ;
376     if(strstr(branch->GetName(),"PHOS") ){
377       message += "       " ; 
378       message += branch->GetName() ; 
379       message += "     " ; 
380       message += branch->GetTitle() ;
381       message += "\n" ; 
382     }
383   }
384   branches = gAlice->TreeD()->GetListOfBranches() ;
385   
386   message += "TreeD:\n" ;
387   for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
388     TBranch * branch=(TBranch *) branches->At(ibranch) ;
389     if(strstr(branch->GetName(),"PHOS") ) {
390       message += "       "; 
391       message += branch->GetName() ; 
392       message += "     " ; 
393       message += branch->GetTitle() ; 
394       message +="\n" ;
395     }
396   }
397   
398   branches = gAlice->TreeR()->GetListOfBranches() ;
399   
400   message += "TreeR: \n" ;
401   for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
402     TBranch * branch=(TBranch *) branches->At(ibranch) ;
403     if(strstr(branch->GetName(),"PHOS") ) {
404       message += "       " ; 
405       message += branch->GetName() ; 
406       message += "     " ;
407       message += branch->GetTitle() ; 
408       message += "\n" ;
409     }
410   }
411   Info("LS", message.Data()) ;  
412 }
413 //____________________________________________________________________________
414  void AliPHOSAnalyze::InvariantMass(const char* branchTitle)    
415 {
416   // Calculates Real and Mixed invariant mass distributions
417   AliPHOSGetter * gime  = AliPHOSGetter::GetInstance(ffileName.Data(),branchTitle) ;
418
419   Int_t nMixedEvents = 4 ; //# of events used for calculation of 'mixed' distribution 
420
421   //opening file
422   TFile * mfile = new TFile("invmass.root","update");
423   
424   //========== Reading /Booking Histograms
425   TH2D * hRealEM = 0 ;
426   hRealEM = (TH2D*) mfile->Get("hRealEM") ;
427   if(hRealEM == 0) 
428     hRealEM = new TH2D("hRealEM",   "Real for EM particles",      250,0.,1.,40,0.,4.) ;
429   TH2D * hRealPhot = 0 ;
430
431   hRealPhot = (TH2D*)mfile->Get("hRealPhot");
432   if(hRealPhot == 0)
433     hRealPhot = new TH2D("hRealPhot", "Real for kPhoton particles", 250,0.,1.,40,0.,4.) ;
434
435   TH2D * hMixedEM = 0 ;
436   hMixedEM = (TH2D*) mfile->Get("hMixedEM") ;
437   if(hMixedEM == 0)
438     hMixedEM = new TH2D("hMixedEM",  "Mixed for EM particles",     250,0.,1.,40,0.,4.) ;
439
440   TH2D * hMixedPhot = 0 ;
441   hMixedPhot = (TH2D*) mfile->Get("hMixedPhot") ;
442   if(hMixedPhot == 0)
443     hMixedPhot = new TH2D("hMixedPhot","Mixed for kPhoton particles",250,0.,1.,40,0.,4.) ;
444   
445
446   //reading event and copyng it to TConesArray of all photons
447
448   TClonesArray * allRecParticleList  = new TClonesArray("AliPHOSRecParticle", 1000) ;
449   Int_t * nRecParticles = new Int_t[nMixedEvents] ; // to mark boundaries of each event in the total list
450   for(Int_t index = 0; index < nMixedEvents; index ++)
451     nRecParticles[index] = 0 ;
452   Int_t iRecPhot = 0 ;                // number of EM particles in total list
453   
454   //scan over all events
455   Int_t event ;
456   Int_t maxevent = (Int_t)gAlice->TreeE()->GetEntries() ; 
457   //  for(event = 0; event < gime->MaxEvent(); event++  ){
458   for(event = 0; event < maxevent; event++  ){
459     gime->Event(event,"R");  //will read only TreeR 
460     
461     //copy EM RecParticles to the "total" list        
462     const AliPHOSRecParticle * recParticle ;
463     Int_t iRecParticle ;
464     TClonesArray * rp = gime->RecParticles() ;
465     if(!rp){
466       Error("InvariantMass", "Can't find RecParticles") ; 
467       return ;
468     }
469
470     for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast(); iRecParticle++ )
471       {
472         recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
473         if((recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)||
474            (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW))
475           new( (*allRecParticleList)[iRecPhot++] ) AliPHOSRecParticle(*recParticle) ;
476       }
477     
478     Int_t mevent = event%nMixedEvents ; //event number in the "mixed" cicle
479     nRecParticles[mevent] = iRecPhot-1 ;  
480     
481     //check, if it is time to calculate invariant mass?
482     Int_t maxevent = (Int_t)gAlice->TreeE()->GetEntries() ; 
483     if((mevent == 0) && (event +1 == maxevent)){
484       
485       //   if((mevent == 0) && (event +1 == gime->MaxEvent())){
486       
487       //calculate invariant mass:
488       Int_t irp1,irp2 ;
489       Int_t nCurEvent = 0 ;
490       
491       for(irp1 = 0; irp1 < allRecParticleList->GetEntries()-1; irp1++){
492         AliPHOSRecParticle * rp1 = (AliPHOSRecParticle *)allRecParticleList->At(irp1) ;
493         
494         for(irp2 = irp1+1; irp2 < allRecParticleList->GetEntries(); irp2++){
495           AliPHOSRecParticle * rp2 = (AliPHOSRecParticle *)allRecParticleList->At(irp2) ;
496           
497           Double_t invMass ;
498           invMass = (rp1->Energy()+rp2->Energy())*(rp1->Energy()+rp2->Energy())-
499             (rp1->Px()+rp2->Px())*(rp1->Px()+rp2->Px())-
500             (rp1->Py()+rp2->Py())*(rp1->Py()+rp2->Py())-
501             (rp1->Pz()+rp2->Pz())*(rp1->Pz()+rp2->Pz()) ;
502           
503           if(invMass> 0)
504             invMass = TMath::Sqrt(invMass);
505           
506           Double_t pt ; 
507           pt = TMath::Sqrt((rp1->Px()+rp2->Px() )*( rp1->Px()+rp2->Px() ) + 
508                            (rp1->Py()+rp2->Py() )*( rp1->Py()+rp2->Py() ) );
509           
510           if(irp1 > nRecParticles[nCurEvent])
511             nCurEvent++;
512           
513           if(irp2 <= nRecParticles[nCurEvent]){ //'Real' event
514             hRealEM->Fill(invMass,pt);
515             if((rp1->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)&&
516                (rp2->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST) )
517               hRealPhot->Fill(invMass,pt);
518           }
519           else{
520             hMixedEM->Fill(invMass,pt);
521             if((rp1->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)&&
522                (rp2->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST) )
523               hMixedPhot->Fill(invMass,pt);
524           } //real-mixed
525           
526       } //loop over second rp
527       }//loop over first rp
528       
529       //Make some cleanings 
530       for(Int_t index = 0; index < nMixedEvents; index ++)
531         nRecParticles[index] = 0 ;
532       iRecPhot = 0 ;              
533       allRecParticleList->Clear() ;
534       
535     }
536   }
537   delete allRecParticleList ;
538   
539   //writing output
540   mfile->cd();
541   
542   hRealEM->Write(0,kOverwrite) ;
543   hRealPhot->Write(0,kOverwrite) ;
544   hMixedEM->Write(0,kOverwrite) ;
545   hMixedPhot->Write(0,kOverwrite) ;
546   
547   mfile->Write();
548   mfile->Close();
549   delete mfile ;
550   delete nRecParticles;
551
552 }
553
554 //____________________________________________________________________________
555  void AliPHOSAnalyze::EnergyResolution(const char * branchTitle)    
556 {
557   //fills two dimentional histo: energy of primary vs. energy of reconstructed
558
559   TH2F * hAllEnergy = 0 ;  //all reconstructed with primary photon   
560   TH2F * hPhotEnergy= 0 ;  //kGamma  with primary photon   
561   TH2F * hEMEnergy  = 0 ;  //electromagnetic with primary photon   
562
563   //opening file and reading histograms if any
564   TFile * efile = new TFile("energy.root","update");
565
566   hAllEnergy = (TH2F*)efile->Get("hAllEnergy") ;
567   if(hAllEnergy == 0)
568     hAllEnergy  = new TH2F("hAllEnergy",  "Energy of any RP with primary photon",100, 0., 5., 100, 0., 5.);
569
570   hPhotEnergy =(TH2F*) efile->Get("hPhotEnergy") ;
571   if(hPhotEnergy == 0)
572     hPhotEnergy = new TH2F("hPhotEnergy", "Energy of kGAMMA with primary photon",100, 0., 5., 100, 0., 5.);
573
574   hEMEnergy =(TH2F*) efile->Get("hEMEnergy");
575   if(hEMEnergy == 0)
576     hEMEnergy   = new TH2F("hEMEnergy",   "Energy of EM with primary photon",    100, 0., 5., 100, 0., 5.);
577
578
579   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(ffileName.Data(),branchTitle) ;
580   const AliPHOSGeometry * phosgeom = gime->PHOSGeometry() ; 
581
582   Int_t ievent;
583   Int_t maxevent = (Int_t)gAlice->TreeE()->GetEntries() ; 
584   for ( ievent=0; ievent < maxevent ; ievent++){
585
586     //read the current event
587     gime->Event(ievent) ;
588
589     const AliPHOSRecParticle * recParticle ;
590     Int_t iRecParticle ;
591     TClonesArray * rp = gime->RecParticles() ;
592     if(!rp) {
593       Error("EnergyResolution", "Event %d,  Can't find RecParticles ", ievent) ;  
594       return ;
595     }
596     TClonesArray * ts = gime->TrackSegments() ;
597     if(!ts) {
598       Error("EnergyResolution", "Event %d,  Can't find TrackSegments", ievent) ;  
599       return ;
600     }
601     TObjArray * emcrp = gime->EmcRecPoints() ;
602     if(!emcrp){
603       Error("EnergyResolution", "Event %d,  Can't find EmcRecPoints") ; 
604       return ;
605     }
606       
607     for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast() ;iRecParticle++ ){
608       recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
609       
610       //find the closest primary
611       Int_t moduleNumberRec ;
612       Double_t recX, recZ ;
613       phosgeom->ImpactOnEmc(recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
614       
615       Double_t minDistance  = 100. ;
616       Int_t closestPrimary = -1 ;
617       
618       //extract list of primaries: it is stored at EMC RecPoints
619       Int_t emcIndex = ((AliPHOSTrackSegment*) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
620       Int_t numberofprimaries ;
621       Int_t * listofprimaries  = ((AliPHOSEmcRecPoint*) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries)  ;
622       
623       Int_t index ;
624       const TParticle * primary ;
625       Double_t distance = minDistance ;
626       Double_t dX, dZ; 
627       Double_t dXmin = 0.; 
628       Double_t dZmin = 0. ;
629       for ( index = 0 ; index < numberofprimaries ; index++){
630         primary = gime->Primary(listofprimaries[index]) ;
631         Int_t moduleNumber ;
632         Double_t primX, primZ ;
633         phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
634         if(moduleNumberRec == moduleNumber) {
635           dX = recX - primX;
636           dZ = recZ - primZ;
637           distance = TMath::Sqrt(dX*dX + dZ*dZ) ;
638           if(minDistance > distance) {
639             minDistance = distance ;
640             dXmin = dX;
641             dZmin = dZ;
642             closestPrimary = listofprimaries[index] ;
643           }
644         }
645       }
646
647       //if found primary, fill histograms
648       if(closestPrimary >=0 ){
649         const TParticle * primary = gime->Primary(closestPrimary) ;
650         if(primary->GetPdgCode() == 22){
651           hAllEnergy->Fill(primary->Energy(), recParticle->Energy()) ;
652           if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST){
653             hPhotEnergy->Fill(primary->Energy(), recParticle->Energy() ) ; 
654             hEMEnergy->Fill(primary->Energy(), recParticle->Energy() ) ; 
655           }
656           else
657             if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW)
658               hEMEnergy->Fill(primary->Energy(), recParticle->Energy() ) ; 
659         }
660       }
661     }
662   }
663
664   //write filled histograms
665   efile->cd() ;
666   hAllEnergy->Write(0,kOverwrite) ;
667   hPhotEnergy->Write(0,kOverwrite) ;
668   hEMEnergy->Write(0,kOverwrite)  ;
669   //  efile->Write() ;
670   efile->Close() ;
671   delete efile ;
672
673 }
674 //____________________________________________________________________________
675 void AliPHOSAnalyze::PositionResolution(const char * branchTitle)    
676 {
677   //fills two dimentional histo: energy vs. primary - reconstructed distance  
678
679
680
681   TH2F * hAllPosition  = 0;     // Position of any RP with primary photon
682   TH2F * hPhotPosition = 0;    // Position of kGAMMA with primary photon
683   TH2F * hEMPosition   = 0;      // Position of EM with primary photon
684
685   TH1F * hAllPositionX = 0;    // X-Position Resolution of photons with photon primary
686   TH1F * hAllPositionZ = 0;    // Z-Position Resolution of photons with photon primary
687
688
689   //opening file and reading histograms if any
690   TFile * pfile = new TFile("position.root","update");
691
692   hAllPosition = (TH2F*)pfile->Get("hAllPosition");
693   if(hAllPosition == 0)
694     hAllPosition  = new TH2F("hAllPosition",  
695                              "Position of any RP with primary photon",100, 0., 5., 100, 0., 5.);
696   hPhotPosition= (TH2F*)pfile->Get("hPhotPosition");
697   if(hPhotPosition == 0)
698     hPhotPosition = new TH2F("hPhotPosition", 
699                              "Position of kGAMMA with primary photon",100, 0., 5., 100, 0., 5.);
700   hEMPosition= (TH2F*)pfile->Get("hEMPosition") ;
701   if(hEMPosition == 0)
702     hEMPosition   = new TH2F("hEMPosition",   
703                              "Position of EM with primary photon",    100, 0., 5., 100, 0., 5.);
704   hAllPositionX = (TH1F*)pfile->Get("hAllPositionX") ;                     
705   if(hAllPositionX == 0)
706     hAllPositionX = new TH1F("hAllPositionX", 
707                              "Delta X of any RP with primary photon",100, -2., 2.);
708   hAllPositionZ =(TH1F*) pfile->Get("hAllPositionZ") ;
709   if(hAllPositionZ == 0)
710     hAllPositionZ = new TH1F("hAllPositionZ", 
711                              "Delta X of any RP with primary photon",100, -2., 2.);
712
713
714   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(ffileName.Data(),branchTitle) ;
715   const AliPHOSGeometry * phosgeom = gime->PHOSGeometry() ; 
716
717   Int_t ievent;
718   Int_t maxevent = (Int_t)gAlice->TreeE()->GetEntries() ; 
719   for ( ievent=0; ievent < maxevent ; ievent++){
720     
721     //read the current event
722     gime->Event(ievent) ;
723     TClonesArray * rp = gime->RecParticles() ;
724     if(!rp) {
725       Error("PositionResolution", "Event %d,  Can't find RecParticles", ievent) ;
726       return ;
727     }
728     TClonesArray * ts = gime->TrackSegments() ;
729     if(!ts) {
730       Error("PositionResolution", "Event %d,  Can't find TrackSegments", ievent) ;
731       return ;
732     }
733     TObjArray * emcrp = gime->EmcRecPoints() ;
734     if(!emcrp){
735       Error("PositionResolution", "Event %d,  Can't find EmcRecPoints", ievent) ;
736       return ;
737     }
738     
739  
740     const AliPHOSRecParticle * recParticle ;
741     Int_t iRecParticle ;
742     for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast(); iRecParticle++ ){
743       recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
744       
745       //find the closest primary
746       Int_t moduleNumberRec ;
747       Double_t recX, recZ ;
748       phosgeom->ImpactOnEmc(recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
749       
750       Double_t minDistance  = 100. ;
751       Int_t closestPrimary = -1 ;
752       
753       //extract list of primaries: it is stored at EMC RecPoints
754       Int_t emcIndex = ((AliPHOSTrackSegment*) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
755       Int_t numberofprimaries ;
756       Int_t * listofprimaries  = ((AliPHOSEmcRecPoint *) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries)  ;
757
758       Int_t index ;
759       const TParticle * primary ;
760       Double_t distance = minDistance ;
761       Double_t dX = 1000; // incredible number
762       Double_t dZ = 1000; // for the case if no primary will be found
763       Double_t dXmin = 0.; 
764       Double_t dZmin = 0. ;
765       for ( index = 0 ; index < numberofprimaries ; index++){
766         primary = gime->Primary(listofprimaries[index]) ;
767         Int_t moduleNumber ;
768         Double_t primX, primZ ;
769         phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
770         if(moduleNumberRec == moduleNumber) {
771           dX = recX - primX;
772           dZ = recZ - primZ;
773           distance = TMath::Sqrt(dX*dX + dZ*dZ) ;
774           if(minDistance > distance) {
775             minDistance = distance ;
776             dXmin = dX;
777             dZmin = dZ;
778             closestPrimary = listofprimaries[index] ;
779           }
780         }
781       }
782       
783       //if found primary, fill histograms
784       if(closestPrimary >=0 ){
785         const TParticle * primary = gime->Primary(closestPrimary) ;
786         if(primary->GetPdgCode() == 22){
787           hAllPosition->Fill(primary->Energy(), minDistance) ;
788           hAllPositionX->Fill(primary->Energy(), dX) ;
789           hAllPositionZ->Fill(primary->Energy(), dZ) ;
790           if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST){
791             hPhotPosition->Fill(primary->Energy(), minDistance ) ; 
792             hEMPosition->Fill(primary->Energy(), minDistance ) ; 
793           }
794           else
795             if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW)
796               hEMPosition->Fill(primary->Energy(), minDistance ) ; 
797         }
798       }
799     }
800   }
801   
802   //Write output histgrams
803   pfile->cd() ;
804   hAllPosition->Write(0,kOverwrite) ;
805   hAllPositionX->Write(0,kOverwrite) ;
806   hAllPositionZ->Write(0,kOverwrite) ;
807   hPhotPosition->Write(0,kOverwrite) ;
808   hEMPosition->Write(0,kOverwrite) ;
809   pfile->Write() ;
810   pfile->Close() ;
811   delete pfile ;
812
813
814 }
815 //____________________________________________________________________________
816 void AliPHOSAnalyze::Contamination(const char* RecPointsTitle){
817 // fills spectra of primary photons and several kinds of 
818 // reconstructed particles, so that analyzing them one can 
819 // estimate conatmination, efficiency of registration etc.
820
821   //define several general histograms
822   TH1F * hPrimary = 0;   //spectrum (P_t distribution) of primary photons         
823   TH1F * hAllRP   = 0;   //spectrum of all RecParticles in PHOS
824   TH1F * hPhot    = 0;   //spectrum of kGAMMA RecParticles
825   TH1F * hShape   = 0;   //spectrum of all EM RecParticles
826   TH1F * hVeto    = 0;   //spectrum of all neutral RecParticles
827
828   //Now separate histograms in accoradance with primary
829   //primary - photon
830   TH1F * hPhotReg = 0;   //Registeres as photon
831   TH1F * hPhotEM  = 0;   //Registered as EM       
832
833   //primary - n
834   TH1F * hNReg = 0;   //Registeres as photon          
835   TH1F * hNEM  = 0;   //Registered as EM            
836
837   //primary - nBar
838   TH1F * hNBarReg = 0;   //Registeres as photon
839   TH1F * hNBarEM  = 0;   //Registered as EM          
840
841   //primary - charged hadron (pBar excluded)
842   TH1F * hChargedReg = 0;   //Registeres as photon  
843   TH1F * hChargedEM  = 0;   //Registered as EM           
844
845   //primary - pBar
846   TH1F * hPbarReg = 0;   //Registeres as photon  
847   TH1F * hPbarEM  = 0;   //Registered as EM 
848
849
850   //Reading histograms from the file
851   TFile * cfile = new TFile("contamination.root","update") ;
852
853   //read general histograms
854   hPrimary = (TH1F*) cfile->Get("hPrimary") ;
855   if(hPrimary == 0)
856     hPrimary= new TH1F("hPrimary", "Primary photon spectrum", 100, 0., 5.);
857   hAllRP = (TH1F*)cfile->Get("hAllRP") ;
858   if(hAllRP == 0)
859     hAllRP = new TH1F("hAllRP","All Reconstructed particles", 100, 0., 5.);
860   hPhot  = (TH1F*)cfile->Get("hPhot") ;
861   if(hPhot == 0)
862     hPhot = new TH1F("hPhot","All kGAMMA RecParticles",100, 0., 5.);
863   hShape = (TH1F*) cfile->Get("hShape") ;
864   if(hShape == 0)
865     hShape = new TH1F("hShape","All particles with EM shower",100, 0., 5.);
866   hVeto= (TH1F*)cfile->Get("hVeto") ;
867   if(hVeto == 0) 
868     hVeto  = new TH1F("hVeto", "All uncharged particles",      100, 0., 5.);
869
870
871   //primary - photon
872   hPhotReg = (TH1F*)cfile->Get("hPhotReg");
873   if(hPhotReg == 0)
874     hPhotReg   = new TH1F("hPhotReg","Photon registered as photon",100, 0., 5.);
875   hPhotEM  =(TH1F*)cfile->Get("hPhotEM");
876   if(hPhotEM== 0)
877     hPhotEM   = new TH1F("hPhotEM",  "Photon registered as EM", 100, 0., 5.);
878
879   //primary - n
880   hNReg = (TH1F*)cfile->Get("hNReg");
881   if(hNReg== 0)
882    hNReg      = new TH1F("hNReg",   "N registered as photon",              100, 0., 5.);
883   hNEM  = (TH1F*)cfile->Get("hNEM"); 
884   if(hNEM== 0)
885     hNEM      = new TH1F("hNEM",    "N registered as EM",      100, 0., 5.);
886
887   //primary - nBar
888   hNBarReg =(TH1F*)cfile->Get("hNBarReg");
889   if(hNBarReg== 0)
890    hNBarReg   = new TH1F("hNBarReg", "NBar registered as photon",           100, 0., 5.);
891   hNBarEM  =(TH1F*)cfile->Get("hNBarEM"); 
892   if(hNBarEM== 0)
893     hNBarEM   = new TH1F("hNBarEM",  "NBar registered as EM",   100, 0., 5.);
894
895   //primary - charged hadron (pBar excluded)
896   hChargedReg = (TH1F*)cfile->Get("hChargedReg");
897   if(hChargedReg== 0)
898     hChargedReg= new TH1F("hChargedReg", "Charged hadron registered as photon",100, 0., 5.);
899   hChargedEM  = (TH1F*)cfile->Get("hChargedEM"); 
900   if(hChargedEM== 0)
901     hChargedEM= new TH1F("hChargedEM","Charged registered as EM",100, 0., 5.);
902  
903   //primary - pBar
904   hPbarReg = (TH1F*)cfile->Get("hPbarReg");
905   if(hPbarReg== 0)
906     hPbarReg= new TH1F("hPbarReg", "pBar registered as photon",100, 0., 5.);
907   hPbarEM  = (TH1F*)cfile->Get("hPbarEM");
908   if(hPbarEM== 0)
909     hPbarEM= new TH1F("hPbarEM","Pbar registered as EM",100, 0., 5.);
910   
911
912   //Now make some initializations
913
914   Int_t counter[8][5] ;      //# of registered particles 
915   Int_t i1,i2 ;
916   for(i1 = 0; i1<8; i1++)
917     for(i2 = 0; i2<5; i2++)
918       counter[i1][i2] = 0 ;
919
920
921
922   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(ffileName.Data(),RecPointsTitle) ;
923   const AliPHOSGeometry * phosgeom = gime->PHOSGeometry() ; 
924   
925   Int_t ievent;
926   Int_t maxevent = (Int_t)gAlice->TreeE()->GetEntries() ; 
927   for ( ievent=0; ievent < maxevent ; ievent++){
928     
929     gime->Event(ievent) ;
930     
931     TClonesArray * rp = gime->RecParticles() ;
932     if(!rp) {
933       Error("Contamination", "Event %d,  Can't find RecParticles", ievent) ;
934       return ;
935     }
936     TClonesArray * ts = gime->TrackSegments() ;
937     if(!ts) {
938       Error("Contamination", "Event %d,  Can't find TrackSegments", ievent) ;
939       return ;
940     }
941     TObjArray * emcrp = gime->EmcRecPoints() ;
942     if(!emcrp){
943       Error("Contamination", "Event %d,  Can't find EmcRecPoints", ievent) ;
944       return ;
945     }
946     
947     
948     //=========== Make spectrum of the primary photons
949     const TParticle * primary ;
950     Int_t iPrimary ;
951     for( iPrimary = 0 ; iPrimary < gime->NPrimaries() ; iPrimary++){
952       primary = gime->Primary(iPrimary) ;
953       Int_t primaryType = primary->GetPdgCode() ;
954       if( primaryType == 22 ) {
955         //check, if photons folls onto PHOS
956         Int_t moduleNumber ;
957         Double_t primX, primZ ;
958         phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
959         if(moduleNumber)
960           hPrimary->Fill(primary->Energy()) ;
961         
962       }
963       
964     }
965     
966     //========== Now scan over RecParticles            
967     const AliPHOSRecParticle * recParticle ;
968     Int_t iRecParticle ;
969     for(iRecParticle = 0; iRecParticle < rp->GetEntriesFast(); iRecParticle++ ){
970       recParticle = (AliPHOSRecParticle *) rp->At(iRecParticle) ;
971       //fill histo spectrum of all RecParticles
972       hAllRP->Fill(CorrectedEnergy(recParticle->Energy())) ;
973       
974       //==========find the closest primary      
975       Int_t moduleNumberRec ;
976       Double_t recX, recZ ;
977       phosgeom->ImpactOnEmc(recParticle->Theta(), recParticle->Phi(), moduleNumberRec, recX, recZ) ;
978       
979       Double_t minDistance  = 100. ;
980       Int_t closestPrimary = -1 ;
981       
982       //extract list of primaries: it is stored at EMC RecPoints
983       Int_t emcIndex = ((AliPHOSTrackSegment *) ts->At(recParticle->GetPHOSTSIndex()))->GetEmcIndex() ;
984       Int_t numberofprimaries ;
985       Int_t * listofprimaries  = ((AliPHOSEmcRecPoint *) emcrp->At(emcIndex))->GetPrimaries(numberofprimaries)  ;
986       Int_t index ;
987       const TParticle * primary ;
988       Double_t distance = minDistance ;
989       Double_t dX, dZ; 
990       Double_t dXmin = 0.; 
991       Double_t dZmin = 0. ;
992       for ( index = 0 ; index < numberofprimaries ; index++){
993         primary = gime->Primary(listofprimaries[index]) ;
994         Int_t moduleNumber ;
995         Double_t primX, primZ ;
996         phosgeom->ImpactOnEmc(primary->Theta(), primary->Phi(), moduleNumber, primX, primZ) ;
997         if(moduleNumberRec == moduleNumber) {
998           dX = recX - primX;
999           dZ = recZ - primZ;
1000           distance = TMath::Sqrt(dX*dX + dZ*dZ) ;
1001           if(minDistance > distance) {
1002             minDistance = distance ;
1003             dXmin = dX;
1004             dZmin = dZ;
1005             closestPrimary = listofprimaries[index] ;
1006           }
1007         }
1008       }
1009       
1010       //===========define the "type" of closest primary
1011       if(closestPrimary >=0 ){
1012         Int_t primaryCode = -1;
1013         const TParticle * primary = gime->Primary(closestPrimary) ;
1014         Int_t primaryType = primary->GetPdgCode() ;
1015         if(primaryType == 22) // photon ?
1016           primaryCode = 0 ;
1017         else
1018           if(primaryType == 2112) // neutron
1019             primaryCode = 1 ; 
1020           else
1021             if(primaryType == -2112) // Anti neutron
1022               primaryCode = 2 ;
1023             else
1024               if(primaryType == -2122) //Anti proton
1025                 primaryCode = 4 ;
1026               else {
1027                 TParticle tempo(*primary) ; 
1028                 if(tempo.GetPDG()->Charge())
1029                   primaryCode = 3 ;
1030               }
1031
1032         //==========Now look at the type of RecParticle
1033         Float_t energy = CorrectedEnergy(recParticle->Energy()) ;
1034         if(recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST){
1035           hPhot->Fill(energy ) ;        
1036           switch(primaryCode){
1037           case 0:
1038             hPhotReg->Fill(energy ) ; 
1039             break ;
1040           case 1:
1041             hNReg->Fill(energy ) ; 
1042             break ;
1043           case 2:
1044             hNBarReg->Fill(energy ) ; 
1045             break ;
1046           case 3:
1047             hChargedReg->Fill(energy ) ;
1048             break ;
1049           case 4:
1050             hPbarReg->Fill(energy ) ;
1051             break ;
1052           default:
1053             break ;
1054           }
1055         }
1056         if((recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)||
1057            (recParticle->GetType() == AliPHOSFastRecParticle::kCHARGEDEMFAST)||
1058            (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW)||
1059            (recParticle->GetType() == AliPHOSFastRecParticle::kCHARGEDEMSLOW) ){ //with EM shower
1060           hShape->Fill(energy ) ;
1061           switch(primaryCode){
1062           case 0:
1063             hPhotEM->Fill(energy ) ; 
1064             break ;
1065           case 1:
1066             hNEM->Fill(energy ) ; 
1067             break ;
1068           case 2:
1069             hNBarEM->Fill(energy ) ; 
1070             break ;
1071           case 3:
1072             hChargedEM->Fill(energy ) ; 
1073             break ;
1074           case 4:
1075             hPbarEM->Fill(energy ) ; 
1076             break ;
1077           default:
1078             break ;
1079           }
1080         }
1081         
1082         if((recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMFAST)||
1083            (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALHAFAST) ||
1084            (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALEMSLOW) ||
1085            (recParticle->GetType() == AliPHOSFastRecParticle::kNEUTRALHASLOW) ) //nuetral
1086           hVeto->Fill(energy ) ;
1087         
1088         //fill number of primaries identified as ...
1089         if(primaryCode >= 0) // Primary code defined
1090           counter[recParticle->GetType()][primaryCode]++ ; 
1091         
1092       }
1093       
1094     } // no closest primary found
1095   }     
1096   
1097   
1098   //===================  SaveHistograms
1099   cfile->cd() ;
1100   hPrimary->Write(0,kOverwrite); 
1101   hAllRP->Write(0,kOverwrite);  
1102   hPhot->Write(0,kOverwrite);  
1103   hShape->Write(0,kOverwrite); 
1104   hVeto->Write(0,kOverwrite);  
1105   hPhotReg->Write(0,kOverwrite); 
1106   hPhotEM->Write(0,kOverwrite);   
1107   hNReg ->Write(0,kOverwrite);  
1108   hNEM  ->Write(0,kOverwrite); 
1109   hNBarReg ->Write(0,kOverwrite); 
1110   hNBarEM  ->Write(0,kOverwrite); 
1111   hChargedReg ->Write(0,kOverwrite); 
1112   hChargedEM  ->Write(0,kOverwrite); 
1113   hPbarReg ->Write(0,kOverwrite); 
1114   hPbarEM  ->Write(0,kOverwrite); 
1115   
1116   cfile->Write(0,kOverwrite); 
1117   cfile->Close();
1118   delete cfile ;
1119  
1120   
1121   //print Final Table
1122   maxevent = (Int_t)gAlice->TreeE()->GetEntries() ; 
1123
1124   TString message ; 
1125   message  = "Resolutions: Analyzed %d event(s)\n" ; 
1126  
1127   message += "        Primary:    Photon  Neutron  Antineutron  Charged hadron   AntiProton\n" ; 
1128   message += "--------------------------------------------------------------------------------\n" ;
1129   message += "         kGAMMA: " ; 
1130   message += "%d %d %d %d %d\n" ; 
1131   message += "       kGAMMAHA: " ;
1132   message += "%d %d %d %d %d\n" ; 
1133   message += "     kNEUTRALEM: " ; 
1134   message += "%d %d %d %d %d\n" ; 
1135   message += "     kNEUTRALHA: " ; 
1136   message += "%d %d %d %d %d\n" ;
1137   message += "      kABSURDEM: ";
1138   message += "%d %d %d %d %d\n" ;
1139   message += "      kABSURDHA: " ;
1140   message += "%d %d %d %d %d\n" ;
1141   message += "      kELECTRON: " ;
1142   message += "%d %d %d %d %d\n" ;
1143   message += "     kCHARGEDHA: " ;  
1144   message += "%d %d %d %d %d\n" ;
1145    
1146   message += "--------------------------------------------------------------------------------" ;
1147
1148  
1149   Int_t totalInd = 0 ;
1150   for(i1 = 0; i1<8; i1++)
1151     for(i2 = 0; i2<5; i2++)
1152       totalInd+=counter[i1][i2] ;
1153   message += "Indentified particles: %d" ; 
1154   
1155  Info("Contamination", message.Data(), maxevent, 
1156       counter[2][0], counter[2][1], counter[2][2], counter[2][3], counter[2][4], 
1157       counter[3][0], counter[3][1], counter[3][2], counter[3][3], counter[3][4], 
1158       counter[0][0], counter[0][1], counter[0][2], counter[0][3], counter[0][4], 
1159       counter[1][0], counter[1][1], counter[1][2], counter[1][3], counter[1][4], 
1160       counter[4][0], counter[4][1], counter[4][2], counter[4][3], counter[4][4], 
1161       counter[5][0], counter[5][1], counter[5][2], counter[5][3], counter[5][4], 
1162       counter[6][0], counter[6][1], counter[6][2], counter[6][3], counter[6][4], 
1163       counter[7][0], counter[7][1], counter[7][2], counter[7][3], counter[7][4], 
1164       totalInd ) ;
1165
1166 }
1167
1168
1169