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