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