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