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