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