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