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