]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDC.cxx
Added AliL3Stopwatch.
[u/mrichter/AliRoot.git] / ZDC / AliZDC.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 /*
17 $Log$
18 Revision 1.29  2002/06/04 08:17:04  coppedis
19 Reconstruction method improved
20
21 Revision 1.28  2002/02/04 09:18:08  coppedis
22 Merging and reconstruction code review
23
24 Revision 1.26.2.2  2001/11/12 18:41:44  hristov
25 All the changes from the head are merged to the release
26
27 Revision 1.27  2001/10/21 18:27:45  hristov
28 Several pointers were set to zero in the default constructors to avoid memory management problems
29
30 Revision 1.26  2001/10/04 14:30:28  coppedis
31 Event merging for ZDC
32
33 Revision 1.25  2001/10/04 14:24:15  coppedis
34 Event merging for ZDC
35
36 Revision 1.24  2001/09/26 16:03:41  coppedis
37 Merging implemented
38
39 Revision 1.23  2001/05/15 13:44:57  coppedis
40 Changes in AddHit method
41
42 Revision 1.22  2001/05/14 09:53:32  coppedis
43 Adding functions ZMin and ZMax
44
45 Revision 1.21  2001/04/20 10:05:02  coppedis
46 Minor changes
47
48 Revision 1.20  2001/03/26 13:39:20  coppedis
49 Comment prints
50
51 Revision 1.19  2001/03/26 09:10:23  coppedis
52 Corrected bug in constructor (fIshunt has to be =1)
53
54 Revision 1.18  2001/03/20 08:21:55  coppedis
55 ZDC needs PIPE, ABSO, DIPO and SHIL
56
57 Revision 1.17  2001/03/16 16:18:03  coppedis
58 Correction for superposition of ZDC volumes with MUON arm one
59
60 Revision 1.16  2001/03/15 16:01:11  coppedis
61 Code review
62
63 Revision 1.15  2001/01/26 19:56:27  hristov
64 Major upgrade of AliRoot code
65
66 Revision 1.14  2000/12/12 13:17:01  coppedis
67 Minor corrections suggested by P. Hristov
68
69 Revision 1.12  2000/12/01 08:19:01  coppedis
70 Adding a message error if ZDC is constructed without DIPO
71
72 Revision 1.11  2000/11/30 17:21:03  coppedis
73 Introduce hit array fStHits reset only at the end of the event (for digitization)
74
75 Revision 1.10  2000/11/22 11:32:58  coppedis
76 Major code revision
77
78 Revision 1.9  2000/10/02 21:28:20  fca
79 Removal of useless dependecies via forward declarations
80
81 Revision 1.8  2000/07/10 13:58:01  fca
82 New version of ZDC from E.Scomparin & C.Oppedisano
83
84 Revision 1.7  2000/01/19 17:17:40  fca
85
86 Revision 1.6  1999/09/29 09:24:35  fca
87 Introduction of the Copyright and cvs Log
88
89 */
90
91 ///////////////////////////////////////////////////////////////////////////////
92 //                                                                           //
93 //                      Zero Degree Calorimeter                              //
94 //           This class contains the basic functions for the ZDCs;           //
95 //            functions specific to one particular geometry are              //
96 //                      contained in the derived classes                     //
97 //                                                                           //
98 ///////////////////////////////////////////////////////////////////////////////
99
100 #include <stdlib.h>
101 #include <iostream.h>
102
103 // --- ROOT system
104 #include <TBRIK.h>
105 #include <TNode.h>
106 #include <TGeometry.h>
107 #include <TFile.h>
108 #include <TTree.h>
109 #include <TDirectory.h>
110 #include <TF1.h>
111
112 // --- AliRoot header files
113 #include "AliZDC.h"
114 #include "AliZDCHit.h"
115 #include "AliZDCMergedHit.h"
116 #include "AliZDCMerger.h"
117 #include "AliZDCDigit.h"
118 #include "AliZDCReco.h"
119 #include "AliDetector.h"
120 #include "AliCallf77.h"
121 #include "AliConst.h"
122 #include "AliMC.h"
123 #include "AliRun.h"
124 #include "AliHeader.h"
125
126  
127 ClassImp(AliZDC)
128  
129 //_____________________________________________________________________________
130 AliZDC::AliZDC()
131 {
132   //
133   // Default constructor for the Zero Degree Calorimeter base class
134   //
135   
136   fIshunt     = 1;
137   fNoShower   = 0;
138
139   fMerger     = 0;
140   fHits       = 0;
141   fNhits      = 0;
142
143   fDigits     = 0;
144   fNdigits    = 0;
145
146   fMergedHits = 0;
147
148   fNRecPoints = 0;
149   fRecPoints  = 0;
150   
151 }
152  
153 //_____________________________________________________________________________
154 AliZDC::AliZDC(const char *name, const char *title)
155   : AliDetector(name,title)
156 {
157   //
158   // Standard constructor for the Zero Degree Calorimeter base class
159   //
160
161   fIshunt   = 1;
162   fNoShower = 0;
163   fMerger   = 0;
164
165   // Allocate the hits array  
166   fHits   = new TClonesArray("AliZDCHit",1000);
167   gAlice->AddHitList(fHits);
168   // Allocate the merged hits array  
169   fMergedHits = new TClonesArray("AliZDCMergedHit",1000);
170
171   // Allocate the digits array  
172   fDigits = new TClonesArray("AliZDCDigit",1000);
173   
174   fNRecPoints = 0;
175   fRecPoints = 0;  
176
177 }
178 //____________________________________________________________________________ 
179 AliZDC::~AliZDC()
180 {
181   //
182   // ZDC destructor
183   //
184
185   fIshunt   = 0;
186   
187   if(fMerger) delete fMerger;
188
189 }
190 //_____________________________________________________________________________
191 void AliZDC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
192 {
193   //
194   //            Add a ZDC hit to the hit list.
195   // -> We make use of 2 array of hits:
196   // [1]  fHits (the usual one) that contains hits for each PRIMARY
197   // [2]  fStHits that contains hits for each EVENT and is used to
198   //      obtain digits at the end of each event
199   //
200   
201   static Float_t primKinEn, xImpact, yImpact, sFlag;
202
203   AliZDCHit *newquad, *curprimquad;
204   newquad = new AliZDCHit(fIshunt, track, vol, hits);
205   TClonesArray &lhits = *fHits;
206   
207   if(fNhits==0){
208       // First hit -> setting flag for primary or secondary particle
209       Int_t primary = gAlice->GetPrimary(track);     
210       if(track != primary){
211         newquad->fSFlag = 1;  // SECONDARY particle entering the ZDC
212       }
213       else if(track == primary){
214         newquad->fSFlag = 0;  // PRIMARY particle entering the ZDC
215       }  
216       sFlag     = newquad->fSFlag;
217       primKinEn = newquad->fPrimKinEn;
218       xImpact   = newquad->fXImpact;
219       yImpact   = newquad->fYImpact;
220    }
221    else{       
222       newquad->fPrimKinEn = primKinEn;
223       newquad->fXImpact = xImpact;
224       newquad->fYImpact = yImpact;
225       newquad->fSFlag   = sFlag;
226    }
227  
228   Int_t j;
229   for(j=0; j<fNhits; j++){
230     // If hits are equal (same track, same volume), sum them.
231      curprimquad = (AliZDCHit*) lhits[j];
232      if(*curprimquad == *newquad){
233         *curprimquad = *curprimquad+*newquad;
234         delete newquad;
235         return;
236      } 
237   }
238
239     //Otherwise create a new hit
240     new(lhits[fNhits]) AliZDCHit(newquad);
241     fNhits++;
242     
243     delete newquad;
244 }
245
246 //_____________________________________________________________________________
247 void  AliZDC::AddDigit(Int_t *sect, Int_t digit)
248 {
249 //
250   AliZDCDigit *newdigit;
251   newdigit = new AliZDCDigit(sect, digit);
252   
253 //  printf("\n  AddDigit -> sector[0] = %d, sector[1] = %d, digit = %d",
254 //         sect[0], sect[1], digit);
255   new((*fDigits)[fNdigits]) AliZDCDigit(*newdigit);
256   fNdigits++;
257   delete newdigit;
258 }
259       
260 //_____________________________________________________________________________
261 void AliZDC::BuildGeometry()
262 {
263   //
264   // Build the ROOT TNode geometry for event display 
265   // in the Zero Degree Calorimeter
266   // This routine is dummy for the moment
267   //
268
269   TNode *node, *top;
270   TBRIK *brik;
271   const int kColorZDC  = kBlue;
272   
273   //
274   top=gAlice->GetGeometry()->GetNode("alice");
275   
276   // ZDC
277     brik = new TBRIK("S_ZDC","ZDC box","void",300,300,5);
278     top->cd();
279     node = new TNode("ZDC","ZDC","S_ZDC",0,0,600,"");
280     node->SetLineColor(kColorZDC);
281     fNodes->Add(node);
282 }
283
284 //_____________________________________________________________________________
285 Int_t AliZDC::DistancetoPrimitive(Int_t , Int_t )
286 {
287   //
288   // Distance from the mouse to the Zero Degree Calorimeter
289   // Dummy routine
290   //
291   return 9999;
292 }
293
294 //____________________________________________________________________________
295 Float_t AliZDC::ZMin(void) const
296 {
297   // Minimum dimension of the ZDC module in z
298   return 11600.;
299 }
300
301 //____________________________________________________________________________
302 Float_t AliZDC::ZMax(void) const
303 {
304   // Maximum dimension of the ZDC module in z
305   return  11750.;
306 }
307   
308
309 //_____________________________________________________________________________
310  void AliZDC::MakeBranch(Option_t *opt, const char *file)
311 {
312   //
313   // Create Tree branches for the ZDC
314   //
315
316   char branchname[10];
317   sprintf(branchname,"%s",GetName());
318   
319   AliDetector::MakeBranch(opt);
320
321   const char *cS = strstr(opt,"S");
322
323   if (gAlice->TreeS() && cS) {
324     if(fMergedHits!=0) fMergedHits->Clear();
325     else fMergedHits = new TClonesArray ("AliZDCMergedHit",1000);
326     MakeBranchInTree(gAlice->TreeS(), 
327                      branchname, &fMergedHits, fBufferSize, file) ;
328     printf("* AliZDC::MakeBranch    * Making Branch %s for SDigits\n\n",branchname);
329   }
330
331     
332   const char *cD = strstr(opt,"D");
333
334   if (gAlice->TreeD() && cD) {
335     if(fDigits!=0) fDigits->Clear();
336     else fDigits = new TClonesArray ("AliZDCDigit",1000);
337     MakeBranchInTree(gAlice->TreeD(), 
338                      branchname, &fDigits, fBufferSize, file) ;
339     printf("* AliZDC::MakeBranch    * Making Branch %s for Digits\n\n",branchname);
340   }
341
342   
343   const char *cR = strstr(opt,"R");
344
345   if (gAlice->TreeR() && cR) {
346     if(fRecPoints==0) fRecPoints = new TClonesArray("AliZDCReco",1000);
347     MakeBranchInTree(gAlice->TreeR(), 
348                      branchname, &fRecPoints, fBufferSize, file) ;
349     printf("* AliZDC::MakeBranch    * Making Branch %s for RecPoints\n\n",branchname);   }
350           
351 }
352
353 //_____________________________________________________________________________
354  void AliZDC::MakeBranchInTreeS(TTree *treeS, const char *file)
355 {
356   // MakeBranchInTree
357   const Int_t kBufferSize = 4000;
358   char  branchname[20];
359   sprintf(branchname,"%s",GetName());
360   MakeBranchInTree(treeS, branchname, &fMergedHits, kBufferSize, file) ;
361   printf("* AliZDC::MakeBranch    * Making Branch %s for SDigits\n\n",branchname);
362
363 }
364 //_____________________________________________________________________________
365  void AliZDC::MakeBranchInTreeD(TTree *treeD, const char *file)
366 {
367   // MakeBranchInTree
368   const Int_t kBufferSize = 4000;
369   char  branchname[20];
370   sprintf(branchname,"%s",GetName());
371   MakeBranchInTree(treeD, branchname, &fDigits, kBufferSize, file) ;
372   printf("* AliZDC::MakeBranch    * Making Branch %s for Digits\n\n",branchname);
373
374 }
375 //_____________________________________________________________________________
376  void AliZDC::MakeBranchInTreeR(TTree *treeR, const char *file)
377 {
378   // MakeBranchInTree
379   const Int_t kBufferSize = 4000;
380   char  branchname[20];
381   sprintf(branchname,"%s",GetName());
382   MakeBranchInTree(treeR, branchname, &fRecPoints, kBufferSize, file) ;
383   printf("* AliZDC::MakeBranch    * Making Branch %s for RecPoints\n\n",branchname);
384
385 }
386 //_____________________________________________________________________________
387 void AliZDC::Hits2SDigits()
388 {
389   printf("\n    Entering AliZDC::SDigits2Digits() ");
390   
391   //----------------------------------------------------------------
392   if(!fMerger){ 
393     printf("    ZDC digitization (without merging)\n");
394
395     AliZDCMergedHit *MHit;
396     Int_t j, sector[2];
397     Float_t MHits[7];
398     fNMergedhits = 0;
399
400     TTree *treeH = gAlice->TreeH();
401     Int_t ntracks = (Int_t) treeH->GetEntries();
402     gAlice->ResetHits();
403   
404     // Tracks loop
405     for(Int_t itrack=0; itrack<ntracks; itrack++){
406        treeH->GetEvent(itrack);
407        for(AliZDCHit* zdcHit=(AliZDCHit*)this->FirstHit(-1); zdcHit;
408                       zdcHit = (AliZDCHit*)this->NextHit()){ 
409                       
410            for(j=0; j<2; j++) sector[j] = zdcHit->GetVolume(j);
411            MHits[0] = zdcHit->GetPrimKinEn();
412            MHits[1] = zdcHit->GetXImpact();
413            MHits[2] = zdcHit->GetYImpact();
414            MHits[3] = zdcHit->GetSFlag();
415            MHits[4] = zdcHit->GetLightPMQ();
416            MHits[5] = zdcHit->GetLightPMC();
417            MHits[6] = zdcHit->GetEnergy();
418        }//Hits loop
419        
420           MHit = new AliZDCMergedHit(sector, MHits);
421           new((*fMergedHits)[fNMergedhits]) AliZDCMergedHit(*MHit);       
422           TClonesArray &sdigits = *fMergedHits;
423           new (sdigits[fNMergedhits]) AliZDCMergedHit(*MHit);
424           fNMergedhits++;
425           delete MHit;
426     }
427     gAlice->TreeS()->Fill();
428     gAlice->TreeS()->AutoSave(); 
429     gAlice->TreeS()->Reset();  
430   }
431   //----------------------------------------------------------------
432   else if(fMerger){
433     printf("    ZDC merging and digitization\n");
434     // ### Initialise merging
435     fMerger -> InitMerging();
436
437     // SDigits tree
438     TTree *treeS = 0;
439     char treeSName[20];
440     sprintf(treeSName,"TreeS%d",fMerger->EvNum());
441     if(gAlice->GetTreeSFile()){
442       gAlice->GetTreeSFile()->cd();
443       treeS = (TTree*)gAlice->GetTreeSFile()->Get(treeSName);
444     }
445     else {
446       treeS = gAlice->TreeS();
447     }
448     if(!treeS){
449       printf("\n ERROR -> Can't find TreeS%d in background file\n",fMerger->EvNum());
450     }    
451
452     // ### Get TCA of MergedHits from AliZDCMerger
453     fMergedHits  = fMerger->MergedHits();
454     fNMergedhits = fMerger->GetNMhits();
455
456     // Branch address
457     char branchSDname[20];
458     sprintf(branchSDname,"%s",GetName());
459     if(treeS && fMergedHits){
460       TBranch *branchSD = treeS->GetBranch(branchSDname);
461       if(branchSD) branchSD->SetAddress(&fMergedHits);
462       else if(!branchSD) MakeBranchInTreeS(treeS);
463     }
464     AliZDCMergedHit *MHit;
465     TClonesArray &sdigits = *fMergedHits;
466     Int_t imhit;
467     //Merged Hits loop
468     for(imhit=0; imhit<fNMergedhits; imhit++){
469        MHit = (AliZDCMergedHit*) fMergedHits->UncheckedAt(imhit);
470        new (sdigits[imhit]) AliZDCMergedHit(*MHit);
471     }
472     treeS->Fill();
473     treeS->AutoSave();
474   }
475   
476 }
477
478 //_____________________________________________________________________________
479 void AliZDC::SDigits2Digits()
480 {
481   if(!fMerger){ // Only digitization
482     printf("    ZDC digitization (no merging) \n");
483     fMerger = new AliZDCMerger();    
484     fMerger->Digitize(fNMergedhits, fMergedHits);
485
486     char hname[30];
487     sprintf(hname,"TreeD%d",gAlice->GetHeader()->GetEvent());
488     gAlice->TreeD()->Fill();
489     gAlice->TreeD()->AutoSave();
490     gAlice->TreeD()->Reset();  
491   }
492   else if(fMerger){     // Merging and digitization
493     printf("    ZDC merging and digitization\n");
494     fMerger->Digitize(fNMergedhits, fMergedHits);
495
496     // Digits tree
497     TTree *treeD = 0;
498     char treeDName[20];
499     sprintf(treeDName,"TreeD%d",fMerger->EvNum());  
500     if(gAlice->GetTreeDFile()){
501       treeD = (TTree*)gAlice->GetTreeDFile()->Get(treeDName);
502     }
503     else {
504       treeD = gAlice->TreeD();
505     }
506     if(!treeD){
507       printf("\n ERROR -> Can't find TreeD%d in background file\n",fMerger->EvNum());
508     }    
509     // Branch address
510     char branchDname[20];
511     sprintf(branchDname,"%s",GetName());
512     if(treeD && fDigits){
513       TBranch *branchD = treeD->GetBranch(branchDname);
514       if(branchD) branchD->SetAddress(&fDigits);
515       else if(!branchD) MakeBranchInTreeD(treeD);
516     }
517     treeD->Fill();
518     treeD->AutoSave();
519   }
520   
521   
522 }
523 //_____________________________________________________________________________
524 void AliZDC::Hits2Digits()
525 {
526     gAlice->Hits2SDigits();
527     gAlice->SDigits2Digits();
528 }
529
530 //_____________________________________________________________________________
531 void AliZDC::Digits2Reco()
532 {
533     printf("    Entering AliZDC::Digits2Reco\n");
534     AliDetector *ZDC  = gAlice->GetDetector("ZDC");
535     TClonesArray *ZDCdigits = ZDC->Digits();
536     
537     TTree *TD = 0;
538     char treeDame[20];
539     sprintf(treeDame,"TreeD%d",fMerger->EvNum());
540     if(gAlice->GetTreeDFile()){
541       TD = (TTree*)gAlice->GetTreeDFile()->Get(treeDame);
542     }
543     else {
544       TD = gAlice->TreeD();
545     }
546     if(TD){
547       char brname[20];
548       sprintf(brname,"%s",ZDC->GetName());
549       TBranch *br = TD->GetBranch(brname);
550       if(br) br->SetAddress(&ZDCdigits);
551     }
552     else if(!TD) printf("       ERROR -> TreeD NOT found in gAlice object\n");
553     
554     Int_t nt = (Int_t) (TD->GetEntries());
555     gAlice->ResetDigits();    
556     
557     AliZDCDigit *dig;
558     Int_t j, idig, ndigits, ZNraw=0, ZPraw=0, ZEMraw=0;
559     //  ---      Summing raw ADCs for each detector to obtain total light
560     for(j=0; j<nt; j++){
561       TD->GetEvent(j);
562       ndigits = ZDCdigits->GetEntries();
563       ZNraw=0;
564       ZPraw=0; 
565       ZEMraw=0;
566       //  ---  Loop over event digits
567       for(idig=0; idig<ndigits; idig++){
568          dig = (AliZDCDigit*) ZDCdigits->UncheckedAt(idig);
569          if(dig->GetSector(0) == 1)      ZNraw  += dig->GetADCValue();
570          else if(dig->GetSector(0) == 2) ZPraw  += dig->GetADCValue();
571          else if(dig->GetSector(0) == 3) ZEMraw += dig->GetADCValue();
572       } // Digits loop
573     } //  TreeD entries loop
574     printf("\n  ---     ZNraw = %d, ZPraw = %d, ZEMraw = %d\n",ZNraw, ZPraw, ZEMraw);
575     
576   //  ---      Pedestal subtraction
577   Int_t ZNcorr, ZPcorr, ZEMcorr, MeanPed=50;
578   ZNcorr  = ZNraw  - 5*MeanPed;
579   ZPcorr  = ZPraw  - 5*MeanPed;
580   ZEMcorr = ZEMraw - 2*MeanPed;
581   if(ZNcorr<0)  ZNcorr=0;
582   if(ZPcorr<0)  ZPcorr=0;
583   if(ZEMcorr<0) ZEMcorr=0;
584  printf("\n    ZNcorr = %d, ZPcorr = %d, ZEMcorr = %d\n",ZNcorr,ZPcorr,ZEMcorr);
585   
586   //  ---      ADCchannel -> photoelectrons
587   // NB-> PM gain = 10^(5), ADC resolution = 6.4*10^(-7)
588   Float_t ZNphe, ZPphe, ZEMphe, ConvFactor = 0.064;
589   ZNphe  = ZNcorr/ConvFactor;
590   ZPphe  = ZPcorr/ConvFactor;
591   ZEMphe = ZEMcorr/ConvFactor;
592   printf("\n    ZNphe = %f, ZPphe = %f, ZEMphe = %f\n",ZNphe, ZPphe, ZEMphe);
593   
594   //  ---      Energy calibration
595   // Conversion factors for hadronic ZDCs goes from phe yield to TRUE incident 
596   //  energy (conversion from GeV to TeV is included); while for EM calos 
597   // conversion is from light yield to detected energy calculated by GEANT
598   // NB -> ZN and ZP conversion factors are constant since incident spectators
599   // have all the same energy, ZEM energy is obtained through a fit over the whole
600   // range of incident particle energies (obtained with full HIJING simulations) 
601   Float_t ZNenergy, ZPenergy, ZEMenergy, ZDCenergy;
602   Float_t ZNphexTeV=329., ZPphexTeV=369.;
603   ZNenergy  = ZNphe/ZNphexTeV;
604   ZPenergy  = ZPphe/ZPphexTeV;
605   ZDCenergy = ZNenergy+ZPenergy;
606   ZEMenergy = -4.81+0.3238*ZEMphe;
607   if(ZEMenergy<0) ZEMenergy=0;
608   printf("    ZNenergy = %f TeV, ZPenergy = %f TeV, ZDCenergy = %f GeV, "
609          "\n            ZEMenergy = %f TeV\n", ZNenergy, ZPenergy, 
610          ZDCenergy, ZEMenergy);
611   
612   if(ZDCenergy==0)
613     printf("\n\n        ###     ATTENZIONE!!! -> ev# %d: ZNenergy = %f TeV, ZPenergy = %f TeV, ZDCenergy = %f GeV, "
614          " ZEMenergy = %f TeV\n\n", fMerger->EvNum(), ZNenergy, ZPenergy, ZDCenergy, ZEMenergy); 
615   
616   //  ---      Number of incident spectator nucleons
617   Int_t NDetSpecN, NDetSpecP;
618   NDetSpecN = (Int_t) (ZNenergy/2.760);
619   NDetSpecP = (Int_t) (ZPenergy/2.760);
620   printf("\n    NDetSpecN = %d, NDetSpecP = %d\n",NDetSpecN, NDetSpecP);
621   
622   //  ---      Number of generated spectator nucleons and impact parameter
623   // --------------------------------------------------------------------------------------------------
624   // [1] ### Results in Chiara's PhD thesis -> 0<b<15 fm (Dec 2001)
625   /*// Fit results for neutrons (Nspectator n true vs. EZN)
626   TF1 *fZNCen = new TF1("fZNCen",
627       "(-2.116909+sqrt(2.116909*2.116909-4*(-0.00651)*(14.556798-x)))/(2*(-0.00651))",0.,158.5);
628   TF1 *fZNPer = new TF1("fZNPer",
629       "(-34.695134-sqrt(34.695134*34.695134-4*(-0.174780)*(-1562.283443-x)))/(2*(-0.174780))",0.,158.5);
630   // Fit results for protons (Nspectator p true vs. EZP)
631   TF1 *fZPCen = new TF1("fZPCen",
632       "(-1.3217+sqrt(1.3217*1.3217-4*(-0.007934)*(4.742873-x)))/(2*(-0.007934))",0.,58.91);
633   TF1 *fZPPer = new TF1("fZPPer",
634       "(-15.788267-sqrt(15.788267*15.788267-4*(-0.133359)*(-383.800673-x)))/(2*(-0.133359))",0.,58.91);
635   // Fit results for total number of spectators (Nspectators true vs. EZDC)
636   TF1 *fZDCCen = new TF1("fZDCCen",
637       "(-1.867335+sqrt(1.867335*1.867335-4*(-0.004119)*(19.100289-x)))/(2*(-0.004119))",0.,220.4);
638   TF1 *fZDCPer = new TF1("fZDCPer",
639       "(-22.429097-sqrt(22.429097*22.429097-4*(-0.072435)*(-1482.034526-x)))/(2*(-0.072435))",0.,220.4);*/
640   // --------------------------------------------------------------------------------------------------
641   // [1] ### Results from a new production  -> 0<b<18 fm (Apr 2002)
642   // Fit results for neutrons (Nspectator n true vs. EZN)
643   TF1 *fZNCen = new TF1("fZNCen",
644       "(-2.287920+sqrt(2.287920*2.287920-4*(-0.007629)*(11.921710-x)))/(2*(-0.007629))",0.,164.);
645   TF1 *fZNPer = new TF1("fZNPer",
646       "(-37.812280-sqrt(37.812280*37.812280-4*(-0.190932)*(-1709.249672-x)))/(2*(-0.190932))",0.,164.);
647   // Fit results for protons (Nspectator p true vs. EZP)
648   TF1 *fZPCen = new TF1("fZPCen",
649        "(-1.321353+sqrt(1.321353*1.321353-4*(-0.007283)*(3.550697-x)))/(2*(-0.007283))",0.,60.);
650   TF1 *fZPPer = new TF1("fZPPer",
651       "(-42.643308-sqrt(42.643308*42.643308-4*(-0.310786)*(-1402.945615-x)))/(2*(-0.310786))",0.,60.);
652   // Fit results for total number of spectators (Nspectators true vs. EZDC)
653   TF1 *fZDCCen = new TF1("fZDCCen",
654       "(-1.934991+sqrt(1.934991*1.934991-4*(-0.004080)*(15.111124-x)))/(2*(-0.004080))",0.,225.);
655   TF1 *fZDCPer = new TF1("fZDCPer",
656       "(-34.380639-sqrt(34.380639*34.380639-4*(-0.104251)*(-2612.189017-x)))/(2*(-0.104251))",0.,225.);
657   // --------------------------------------------------------------------------------------------------
658   // [1] ### Results in Chiara's PhD thesis -> 0<b<15 fm (Dec 2001)
659   /*// Fit results for b (b vs. EZDC)
660   //TF1 *fbCen = new TF1("fbCen","0.611543+0.052231*x-0.000112*x*x+0.000000374*x*x*x",0.,222.);
661   //TF1 *fbPer = new TF1("fbPer","16.552010-0.023866*x-0.00001*x*x",0.,222.);
662   TF1 *fbCen = new TF1("fbCen","0.612769+0.051929*x-0.0001074*x*x+0.0000003724*x*x*x",0.,225.);
663   TF1 *fbPer = new TF1("fbPer","16.6131016-0.026053*x+0.000006893*x*x",0.,225.);*/
664   // --------------------------------------------------------------------------------------------------
665   // [2] ### Results from a new production  -> 0<b<18 fm (Apr 2002)
666   TF1 *fbCen = new TF1("fbCen","-0.056923+0.079703*x-0.0004301*x*x+0.000001366*x*x*x",0.,220.);
667   TF1 *fbPer = new TF1("fbPer","17.943998-0.046846*x+0.000074*x*x",0.,220.);
668   // --------------------------------------------------------------------------------------------------
669   // Evaluating Nspectators and b from ZEM energy
670   // [1] ### Results in Chiara's PhD thesis -> 0<b<15 fm (Dec 2001)
671   /*TF1 *fZEMn  = new TF1("fZEMn","124.2-0.0566*x+0.000006014*x*x",0.,3500.);
672   TF1 *fZEMp  = new TF1("fZEMp","81.3-0.03834*x+0.000004359*x*x",0.,3500.);
673   TF1 *fZEMsp = new TF1("fZEMsp","205.6-0.09567*x+0.00001056*x*x",0.,3500.);
674   TF1 *fZEMb  = new TF1("fZEMb","15.8-0.02084*x+2.802e-5*x*x-2.007e-8*x*x*x+6.586e-12*x*x*x*x-8.042e-16*x*x*x*x*x",0.,3500.);*/
675   // --------------------------------------------------------------------------------------------------
676   // [2] ### Results from a new production  -> 0<b<18 fm (Apr 2002)
677   TF1 *fZEMn  = new TF1("fZEMn","126.2-0.05399*x+0.000005679*x*x",0.,4000.);
678   TF1 *fZEMp  = new TF1("fZEMp","82.49-0.03611*x+0.00000385*x*x",0.,4000.);
679   TF1 *fZEMsp = new TF1("fZEMsp","208.7-0.09006*x+0.000009526*x*x",0.,4000.);
680   TF1 *fZEMb  = new TF1("fZEMb","16.06-0.01633*x+1.44e-5*x*x-6.778e-9*x*x*x+1.438e-12*x*x*x*x-1.112e-16*x*x*x*x*x",0.,4000.);
681   
682   Int_t NGenSpecN=0, NGenSpecP=0, NGenSpec=0;
683   Double_t ImpPar=0;
684   // Cut value for Ezem (GeV)
685   // [1] ### Results in Chiara's PhD thesis -> 0<b<15 fm (Dec 2001)
686   //Float_t EZEMCut = 360.; 
687   // [2] ### Results from a new production  -> 0<b<18 fm (Apr 2002)
688   Float_t EZEMCut = 420.;
689   Float_t DeltaEZEMSup = 690.; 
690   Float_t DeltaEZEMInf = 270.; 
691   if(ZEMenergy > (EZEMCut+DeltaEZEMSup)){
692     NGenSpecN = (Int_t) (fZNCen->Eval(ZNenergy));
693     NGenSpecP = (Int_t) (fZPCen->Eval(ZPenergy));
694     NGenSpec  = (Int_t) (fZDCCen->Eval(ZDCenergy));
695     ImpPar    = fbCen->Eval(ZDCenergy);
696     //printf("    fZNCen = %f, fZPCen = %f, fZDCCen = %f\n",fZNCen->Eval(ZNenergy),
697     //            fZPCen->Eval(ZPenergy),fZDCCen->Eval(ZDCenergy));
698   }
699   else if(ZEMenergy < (EZEMCut-DeltaEZEMInf)){
700     NGenSpecN = (Int_t) (fZNPer->Eval(ZNenergy)); 
701     NGenSpecP = (Int_t) (fZPPer->Eval(ZPenergy));
702     NGenSpec  = (Int_t) (fZDCPer->Eval(ZDCenergy));
703     ImpPar    = fbPer->Eval(ZDCenergy);
704     //printf("    fZNPer = %f, fZPPer = %f, fZDCPer = %f\n",fZNPer->Eval(ZNenergy),
705     //            fZPPer->Eval(ZPenergy),fZDCPer->Eval(ZDCenergy));
706   }
707   else if(ZEMenergy >= (EZEMCut-DeltaEZEMInf) && ZEMenergy <= (EZEMCut+DeltaEZEMSup)){
708     NGenSpecN = (Int_t) (fZEMn->Eval(ZEMenergy));
709     NGenSpecP = (Int_t) (fZEMp->Eval(ZEMenergy));
710     NGenSpec  = (Int_t)(fZEMsp->Eval(ZEMenergy));
711     ImpPar    =  fZEMb->Eval(ZEMenergy);
712     //printf("    Nspec ZEM = %f, Nspec ZDC = %f\n",fZEMsp->Eval(ZNenergy),fZDCPer->Eval(ZDCenergy));
713   }
714   // [1] ### Results in Chiara's PhD thesis -> 0<b<15 fm (Dec 2001)
715   /*if(ZNenergy>158.5)  NGenSpecN = (Int_t) (fZEMn->Eval(ZEMenergy));
716   if(ZPenergy>58.91)  NGenSpecP = (Int_t) (fZEMp->Eval(ZEMenergy));
717   if(ZDCenergy>220.4) NGenSpec  = (Int_t)(fZEMsp->Eval(ZEMenergy));
718   if(ZDCenergy>225.)  ImpPar    =          fZEMb->Eval(ZEMenergy);*/
719   // [2] ### Results from a new production  -> 0<b<18 fm (Apr 2002)
720   if(ZNenergy>162.)  NGenSpecN = (Int_t) (fZEMn->Eval(ZEMenergy));
721   if(ZPenergy>59.75)  NGenSpecP = (Int_t) (fZEMp->Eval(ZEMenergy));
722   if(ZDCenergy>221.5) NGenSpec  = (Int_t)(fZEMsp->Eval(ZEMenergy));
723   if(ZDCenergy>220.)  ImpPar    =  fZEMb->Eval(ZEMenergy);
724   
725   if(NGenSpecN>125)    NGenSpecN=125;
726   else if(NGenSpecN<0) NGenSpecN=0;
727   if(NGenSpecP>82)     NGenSpecP=82;
728   else if(NGenSpecP<0) NGenSpecP=0;
729   if(NGenSpec>207)     NGenSpec=207;
730   else if(NGenSpec<0)  NGenSpec=0;
731   //printf("    NRecSpecN = %d, NRecSpecP = %d, NRecSpec = %d\n",NGenSpecN,NGenSpecP,NGenSpec);
732   
733   //  ---      Number of participants
734   Int_t NPart, NPartTot;
735   NPart = 207-NGenSpecN-NGenSpecP;
736   NPartTot = 207-NGenSpec;
737   //printf("    ###     NPart(ZP+ZN) = %d, NPart(ZDC) = %d, b = %f fm\n",NPart,NPartTot,ImpPar);
738   printf("      ###     NPart = %d, b = %f fm\n",NPartTot,ImpPar);
739   
740   //  ---     Writing RecPoints TCA
741   // Allocate the RecPoints TCA 
742   fRecPoints = new TClonesArray("AliZDCReco",1000);
743   AliZDCReco *reco = new AliZDCReco(ZNenergy,ZPenergy,ZDCenergy,ZEMenergy,
744               NDetSpecN,NDetSpecP,NGenSpecN,NGenSpecP,NGenSpec,NPartTot,ImpPar);
745   new((*fRecPoints)[fNRecPoints]) AliZDCReco(*reco);
746   //fNRecPoints++;
747   //fRecPoints->Dump();
748   delete reco;
749   
750   // TreeR
751   TTree *treeR = gAlice->TreeR();
752   if(!treeR) printf("\n ERROR -> Can't find TreeR%d in background file\n",fMerger->EvNum());
753   // Branch address
754   char branchRname[20];
755   sprintf(branchRname,"%s",GetName());
756   if(fRecPoints){
757     TBranch *branchR = treeR->GetBranch(branchRname);
758     if(branchR) branchR->SetAddress(&fRecPoints);
759     else if(!branchR) MakeBranchInTreeR(treeR);
760   }
761   treeR->Fill();
762   treeR->AutoSave();
763   treeR->Reset();
764 }
765
766