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