]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDC.cxx
Changes to allow NO spectator signal when specattors are tracked in sim
[u/mrichter/AliRoot.git] / ZDC / AliZDC.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 //                      Zero Degree Calorimeter                              //
21 //           This class contains the basic functions for the ZDCs;           //
22 //            functions specific to one particular geometry are              //
23 //                      contained in the derived classes                     //
24 //                                                                           //
25 ///////////////////////////////////////////////////////////////////////////////
26
27 // --- ROOT system
28 #include <TClonesArray.h>
29 #include <TTree.h>
30 #include <TFile.h>
31 #include <TSystem.h>
32 #include <TRandom.h>
33
34 // --- AliRoot header files
35 #include "AliDetector.h"
36 #include "AliRawDataHeaderSim.h"
37 #include "AliRawReader.h"
38 #include "AliLoader.h"
39 #include "AliRun.h"
40 #include "AliMC.h"
41 #include "AliLog.h"
42 #include "AliDAQ.h"
43 #include "AliZDC.h"
44 #include "AliZDCHit.h"
45 #include "AliZDCSDigit.h"
46 #include "AliZDCDigit.h"
47 #include "AliZDCDigitizer.h"
48 #include "AliZDCRawStream.h"
49 #include "AliZDCPedestals.h"
50 #include "AliZDCCalib.h"
51 #include "AliFstream.h"
52
53  
54 ClassImp(AliZDC)
55
56 //_____________________________________________________________________________
57 AliZDC::AliZDC() :
58   AliDetector(),
59   fNoShower(0),
60   fPedCalib(0),
61   fCalibData(0),
62   fZDCCalibFName(""),
63   fSpectatorTracked(1)
64 {
65   //
66   // Default constructor for the Zero Degree Calorimeter base class
67   //
68   
69   fIshunt = 1;
70   fNhits  = 0;
71   fHits = 0;
72   fDigits = 0;
73   fNdigits = 0;
74   
75 }
76  
77 //_____________________________________________________________________________
78 AliZDC::AliZDC(const char *name, const char *title) : 
79   AliDetector(name,title),
80   fNoShower  (0),
81   fPedCalib(0),
82   fCalibData(0),
83   fZDCCalibFName(""),
84   fSpectatorTracked(1)
85 {
86   //
87   // Standard constructor for the Zero Degree Calorimeter base class
88   //
89     
90   fIshunt = 1;
91   fNhits  = 0;
92   fDigits = 0;
93   fNdigits = 0;
94  
95   fHits = new TClonesArray("AliZDCHit",1000);
96   gAlice->GetMCApp()->AddHitList(fHits);
97   
98   char sensname[5],senstitle[25];
99   sprintf(sensname,"ZDC");
100   sprintf(senstitle,"ZDC dummy");
101   SetName(sensname); SetTitle(senstitle);
102
103 }
104
105 //____________________________________________________________________________ 
106 AliZDC::~AliZDC()
107 {
108   //
109   // ZDC destructor
110   //
111
112   fIshunt = 0;
113   delete fPedCalib;
114   delete fCalibData;
115
116 }
117
118 //_____________________________________________________________________________
119 AliZDC::AliZDC(const AliZDC& ZDC) :
120 AliDetector("ZDC","ZDC"),
121 fNoShower(ZDC.fNoShower),
122 fPedCalib(ZDC.fPedCalib),
123 fCalibData(ZDC.fCalibData),
124 fZDCCalibFName(ZDC.fZDCCalibFName),
125 fSpectatorTracked(ZDC.fSpectatorTracked)
126 {
127   // copy constructor
128 }
129
130 //_____________________________________________________________________________
131 AliZDC& AliZDC::operator=(const AliZDC& ZDC)
132 {
133   // assignement operator
134   if(this!=&ZDC){
135     fNoShower = ZDC.fNoShower;
136     fPedCalib = ZDC.fPedCalib;
137     fCalibData = ZDC.fCalibData;
138     fZDCCalibFName = ZDC.fZDCCalibFName;
139   } return *this;
140 }
141
142 //_____________________________________________________________________________
143 void AliZDC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
144 {
145   //
146   //            Add a ZDC hit to the hit list.
147   
148   static Float_t trackTime=0., primKinEn=0., xImpact=0., yImpact=0., sFlag=0.;
149   static Int_t   pcPDGcode;
150
151   AliZDCHit *newquad, *curprimquad;
152   newquad = new AliZDCHit(fIshunt, track, vol, hits);
153   TClonesArray &lhits = *fHits;
154   
155   if(fNhits==0){
156       // First hit -> setting flag for primary or secondary particle
157       Int_t primary = gAlice->GetMCApp()->GetPrimary(track);     
158       //
159       if(track != primary){
160         newquad->SetSFlag(1);  // SECONDARY particle entering the ZDC
161       }
162       else if(track == primary){
163         newquad->SetSFlag(0);  // PRIMARY particle entering the ZDC
164       }  
165       sFlag     = newquad->GetSFlag();
166       primKinEn = newquad->GetPrimKinEn();
167       xImpact   = newquad->GetXImpact();
168       yImpact   = newquad->GetYImpact();
169       pcPDGcode = newquad->GetPDGCode();
170       trackTime = newquad->GetTrackTOF();
171    }
172    else{       
173       newquad->SetPrimKinEn(primKinEn);
174       newquad->SetXImpact(xImpact);
175       newquad->SetYImpact(yImpact);
176       newquad->SetSFlag(sFlag);
177       newquad->SetPDGCode(pcPDGcode);
178       newquad->SetTrackTOF(trackTime);
179    }
180  
181   Int_t j;
182   for(j=0; j<fNhits; j++){
183     // If hits are equal (same track, same volume), sum them.
184      curprimquad = (AliZDCHit*) lhits[j];
185      if(*curprimquad == *newquad){
186         *curprimquad = *curprimquad+*newquad;
187         // Ch. debug
188         //printf("\n\t Summing hits **************** \n", fNhits);
189         //curprimquad->Print("");
190         //
191         delete newquad;
192         return;
193      } 
194   }
195
196     //Otherwise create a new hit
197     new(lhits[fNhits]) AliZDCHit(*newquad);
198     fNhits++;
199     // Ch. debug
200     //printf("\n\t New ZDC hit added! fNhits = %d\n", fNhits);
201     //newquad->Print("");
202     
203     delete newquad;
204 }
205
206 //____________________________________________________________________________
207 Float_t AliZDC::ZMin(void) const
208 {
209   // Minimum dimension of the ZDC module in z
210   return -11600.;
211 }
212
213 //____________________________________________________________________________
214 Float_t AliZDC::ZMax(void) const
215 {
216   // Maximum dimension of the ZDC module in z
217   return 11750.;
218 }
219   
220
221 //_____________________________________________________________________________
222 void AliZDC::MakeBranch(Option_t *opt)
223 {
224   //
225   // Create Tree branches for the ZDC
226   //
227
228   char branchname[10];
229   sprintf(branchname,"%s",GetName());
230
231   const char *cH = strstr(opt,"H");
232   
233   if(cH && fLoader->TreeH()) {
234     if (fHits) {
235       fHits->Clear();
236       fNhits = 0;
237     }
238     else {
239       fHits   = new TClonesArray("AliZDCHit",1000); 
240       if (gAlice && gAlice->GetMCApp())
241         gAlice->GetMCApp()->AddHitList(fHits);
242     }
243   }
244   
245   AliDetector::MakeBranch(opt);
246 }
247
248 //_____________________________________________________________________________
249 void AliZDC::Hits2SDigits()
250 {
251   // Create summable digits from hits
252   
253   AliDebug(1,"\n        AliZDC::Hits2SDigits() ");
254   
255   fLoader->LoadHits("read");
256   fLoader->LoadSDigits("recreate");
257   AliRunLoader* runLoader = fLoader->GetRunLoader();
258   AliZDCSDigit sdigit;
259   AliZDCSDigit* psdigit = &sdigit;
260
261   // Event loop
262   for(Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
263     Float_t pmZNC[5], pmZPC[5], pmZNA[5], pmZPA[5], pmZEM1=0., pmZEM2=0.;
264     for(Int_t i=0; i<4; i++) pmZNC[i] = pmZPC[i] =  pmZNA[i] = pmZPA[i] = 0;
265
266     runLoader->GetEvent(iEvent);
267     TTree* treeH = fLoader->TreeH();
268     Int_t ntracks = (Int_t) treeH->GetEntries();
269     ResetHits();
270
271     // Tracks loop
272     Int_t sector[2]; Float_t trackTime = 0.;
273     for(Int_t itrack = 0; itrack < ntracks; itrack++) {
274       treeH->GetEntry(itrack);
275       for(AliZDCHit* zdcHit = (AliZDCHit*)FirstHit(-1); zdcHit;
276           zdcHit = (AliZDCHit*)NextHit()) { 
277                       
278         sector[0] = zdcHit->GetVolume(0);
279         sector[1] = zdcHit->GetVolume(1);
280         if((sector[1] < 1) || (sector[1]>5)) {
281           Error("Hits2SDigits", "sector[0] = %d, sector[1] = %d", sector[0], sector[1]);
282           continue;
283         }
284         Float_t lightQ = zdcHit->GetLightPMQ();
285         Float_t lightC = zdcHit->GetLightPMC();
286         trackTime = zdcHit->GetTrackTOF();
287         // Signals from ZEM are delayed to arrive in time with ZDC signals
288         if(sector[0] == 3)  trackTime += 320;
289         // Ch. debug
290         //printf("\t det %d vol %d trackTOF %f lightQ %1.0f lightC %1.0f\n",
291         //      sector[0], sector[1], trackTime, lightQ, lightC);
292      
293         if(sector[0] == 1) { //ZNC 
294           pmZNC[0] += lightC;
295           pmZNC[sector[1]] += lightQ;
296         } 
297         else if(sector[0] == 2) { //ZPC 
298           pmZPC[0] += lightC;
299           pmZPC[sector[1]] += lightQ;
300         } 
301         else if(sector[0] == 3) { //ZEM 
302           if(sector[1] == 1) pmZEM1 += lightC;
303           else pmZEM2 += lightQ;
304         }
305         if(sector[0] == 4) { //ZNA 
306           pmZNA[0] += lightC;
307           pmZNA[sector[1]] += lightQ;
308         } 
309         else if(sector[0] == 5) { //ZPA 
310           pmZPA[0] += lightC;
311           pmZPA[sector[1]] += lightQ;
312         } 
313       }//Hits loop
314     }//Tracks loop
315
316     // create the output tree
317     fLoader->MakeTree("S");
318     TTree* treeS = fLoader->TreeS();
319     const Int_t kBufferSize = 4000;
320     treeS->Branch(GetName(), "AliZDCSDigit", &psdigit, kBufferSize);
321
322     // Create sdigits for ZNC
323     sector[0] = 1; // Detector = ZNC
324     for(Int_t j = 0; j < 5; j++) {
325       sector[1] = j; 
326       if(pmZNC[j]>0){
327         new(psdigit) AliZDCSDigit(sector, pmZNC[j], trackTime);
328         treeS->Fill();
329         // Ch. debug
330         //printf("\t SDigit created: det %d quad %d pmZNC[%d] %1.0f trackTOF %f\n",
331         //      sector[0], sector[1], j, pmZNC[j], trackTime);
332       }
333     }
334   
335     // Create sdigits for ZPC
336     sector[0] = 2; // Detector = ZPC
337     for(Int_t j = 0; j < 5; j++) {
338       sector[1] = j; // Towers PM ADCs
339       if(pmZPC[j]>0){
340         new(psdigit) AliZDCSDigit(sector, pmZPC[j], trackTime);
341         treeS->Fill();
342         // Ch. debug
343         //printf("\t SDigit created: det %d quad %d pmZPC[%d] %1.0f trackTOF %f\n",
344         //      sector[0], sector[1], j, pmZPC[j], trackTime);
345       }
346     }
347
348     // Create sdigits for ZEM
349     sector[0] = 3; 
350     sector[1] = 1; // Detector = ZEM1
351     if(pmZEM1>0){ 
352       new(psdigit) AliZDCSDigit(sector, pmZEM1, trackTime);
353       treeS->Fill();
354       // Ch. debug
355       //printf("\t SDigit created: det %d quad %d pmZEM1 %1.0f trackTOF %f\n",
356       //        sector[0], sector[1], pmZEM1, trackTime);
357     }
358     sector[1] = 2; // Detector = ZEM2
359     if(pmZEM2>0){
360       new(psdigit) AliZDCSDigit(sector, pmZEM2, trackTime);
361       treeS->Fill();
362       // Ch. debug
363       //printf("\t SDigit created: det %d quad %d pmZEM2 %1.0f trackTOF %f\n",
364       //        sector[0], sector[1], pmZEM2, trackTime);
365     }
366
367     // Create sdigits for ZNA
368     sector[0] = 4; // Detector = ZNA
369     for(Int_t j = 0; j < 5; j++) {
370       sector[1] = j; // Towers PM ADCs
371       if(pmZNA[j]>0){
372         new(psdigit) AliZDCSDigit(sector, pmZNA[j], trackTime);
373         treeS->Fill();
374         // Ch. debug
375         //printf("\t SDigit created: det %d quad %d pmZNA[%d] %1.0f trackTOF %f\n",
376         //      sector[0], sector[1], j, pmZNA[j], trackTime);
377       }
378     }
379   
380     // Create sdigits for ZPA
381     sector[0] = 5; // Detector = ZPA
382     sector[1] = 0; // Common PM ADC
383     for(Int_t j = 0; j < 5; j++) {
384       sector[1] = j; // Towers PM ADCs
385       if(pmZPA[j]>0){
386         new(psdigit) AliZDCSDigit(sector, pmZPA[j], trackTime);
387         treeS->Fill();
388         // Ch. debug
389         //printf("\t SDigit created: det %d quad %d pmZPA[%d] %1.0f trackTOF %f\n",
390         //      sector[0], sector[1], j, pmZPA[j], trackTime);
391       }
392     }
393
394     // write the output tree
395     fLoader->WriteSDigits("OVERWRITE");
396   }
397
398   fLoader->UnloadHits();
399   fLoader->UnloadSDigits();
400 }
401
402 //_____________________________________________________________________________
403 AliDigitizer* AliZDC::CreateDigitizer(AliRunDigitizer* manager) const
404 {
405   // Create the digitizer for ZDC
406
407   return new AliZDCDigitizer(manager);
408 }
409
410 //_____________________________________________________________________________
411 void AliZDC::Digits2Raw()
412 {
413   // Convert ZDC digits to raw data
414
415   // Format: 24 int values -> ZN1(C+Q1-4), ZP1(C+Q1-4), ZEM1, ZEM2, ZN(C+Q1-4), ZP2(C+Q1-4), 2 Ref PMs
416   //         + 24 int values for the corresponding out of time channels
417   // For the CAEN module V965 we have an Header, the Data Words and an End Of Block
418   //    12 channels x 2 gain chains read from 1st ADC module
419   //    12 channels x 2 gain chains read from 2nd ADC module
420   //    12 channels x 2 gain chains read from 3rd ADC module (o.o.t.)
421   //    12 channels x 2 gain chains read from 4rth ADC module (o.o.t.)
422   //
423   const int knADCData1=24, knADCData2=24; // In principle the 2 numbers can be different!
424   UInt_t lADCHeader1; 
425   UInt_t lADCHeader2; 
426   UInt_t lADCData1[knADCData1];
427   UInt_t lADCData2[knADCData2];
428   UInt_t lADCData3[knADCData1];
429   UInt_t lADCData4[knADCData2];
430   //
431   UInt_t lADCEndBlock;
432
433   // load the digits
434   fLoader->LoadDigits("read");
435   AliZDCDigit digit;
436   AliZDCDigit* pdigit = &digit;
437   TTree* treeD = fLoader->TreeD();
438   if(!treeD) return;
439   treeD->SetBranchAddress("ZDC", &pdigit);
440   //printf("\t AliZDC::Digits2Raw -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
441
442   // Fill data array
443   // ADC header
444   UInt_t lADCHeaderGEO = 0;
445   UInt_t lADCHeaderCRATE = 0;
446   UInt_t lADCHeaderCNT1 = knADCData1;
447   UInt_t lADCHeaderCNT2 = knADCData2;
448     
449   lADCHeader1 = lADCHeaderGEO << 27 | 0x1 << 25 | lADCHeaderCRATE << 16 |
450                lADCHeaderCNT1 << 8 ;
451   lADCHeader2 = lADCHeaderGEO << 27 | 0x1 << 25 | lADCHeaderCRATE << 16 |
452                lADCHeaderCNT2 << 8 ;
453       
454   // ADC data word
455   UInt_t lADCDataGEO = lADCHeaderGEO;
456   //
457   UInt_t lADCDataValue1[knADCData1];
458   UInt_t lADCDataValue2[knADCData2];
459   UInt_t lADCDataValue3[knADCData1];
460   UInt_t lADCDataValue4[knADCData2];
461   //
462   UInt_t lADCDataOvFlw1[knADCData1];
463   UInt_t lADCDataOvFlw2[knADCData2];
464   UInt_t lADCDataOvFlw3[knADCData1];
465   UInt_t lADCDataOvFlw4[knADCData2];
466   //
467   for(Int_t i=0; i<knADCData1 ; i++){
468     lADCDataValue1[i] = 0;
469     lADCDataOvFlw1[i] = 0;
470     lADCDataValue3[i] = 0;
471     lADCDataOvFlw3[i] = 0;
472   }
473   for(Int_t i=0; i<knADCData2 ; i++){
474     lADCDataValue2[i] = 0;
475     lADCDataOvFlw2[i] = 0;
476     lADCDataValue4[i] = 0;
477     lADCDataOvFlw4[i] = 0;
478   }
479   //
480   UInt_t lADCDataChannel = 0;
481   
482   // loop over digits
483   for(Int_t iDigit=0; iDigit<treeD->GetEntries(); iDigit++){
484     treeD->GetEntry(iDigit);
485     if(!pdigit) continue;
486     //digit.Print("");
487     
488     // *** ADC data
489     Int_t index=0;
490     if(digit.GetSector(1)!=5){ // ZDC signal channels
491       // *** ADC1 (ZN1, ZP1, ZEM1,2) or ADC3 (ZN1, ZP1, ZEM1,2 o.o.t.)
492       if(digit.GetSector(0)==1 || digit.GetSector(0)==2 || digit.GetSector(0)==3){
493         if(digit.GetSector(0)==1 || digit.GetSector(0)==2){
494           index = (digit.GetSector(0)-1) + 4*digit.GetSector(1); // ZN1 or ZP1
495           lADCDataChannel = 8*(digit.GetSector(0)-1) + digit.GetSector(1);
496         }
497         else if(digit.GetSector(0)==3){ // ZEM 1,2
498           index = 20 + (digit.GetSector(1)-1);
499           lADCDataChannel = 5 + 8*(digit.GetSector(1)-1);
500         }
501         //
502         /*printf("\t AliZDC::Digits2Raw -> idig%d det %d quad %d index %d, ADCch %d ADCVal[%d, %d]\n",
503                 iDigit,digit.GetSector(0),digit.GetSector(1),index,lADCDataChannel,
504                 digit.GetADCValue(0),digit.GetADCValue(1));// Ch. debug
505         */
506         //
507         if(iDigit<knADCData1){ // *** In-time signals
508           lADCDataValue1[index] = digit.GetADCValue(0);   // High gain ADC ch.    
509           if(lADCDataValue1[index] > 2047) lADCDataOvFlw1[index] = 1; 
510           lADCDataValue1[index+2] = digit.GetADCValue(1); // Low gain ADC ch.
511           if(lADCDataValue1[index+2] > 2047) lADCDataOvFlw1[index+2] = 1; 
512         
513           lADCData1[index] = lADCDataGEO << 27 | 0x1 << 24 | lADCDataChannel << 17 | 
514                           lADCDataOvFlw1[index] << 12 | (lADCDataValue1[index] & 0xfff); 
515           lADCData1[index+2] = lADCDataGEO << 27 | 0x1 << 24  | lADCDataChannel << 17 | 0x1 << 16 |
516                           lADCDataOvFlw1[index+2] << 12 | (lADCDataValue1[index+2] & 0xfff);  
517         }
518         else{ // *** Out-of-time signals
519           lADCDataValue3[index] = digit.GetADCValue(0);   // High gain ADC ch.    
520           if(lADCDataValue3[index] > 2047) lADCDataOvFlw3[index] = 1; 
521           lADCDataValue3[index+2] = digit.GetADCValue(1); // Low gain ADC ch.
522           if(lADCDataValue3[index+2] > 2047) lADCDataOvFlw3[index+2] = 1; 
523       
524           lADCData3[index] = lADCDataGEO << 27 | lADCDataChannel << 17 | 
525                           lADCDataOvFlw3[index] << 12 | (lADCDataValue3[index] & 0xfff); 
526           lADCData3[index+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
527                           lADCDataOvFlw3[index+2] << 12 | (lADCDataValue3[index+2] & 0xfff);  
528         }                  
529       }
530       // *** ADC2 (ZN2, ZP2) or ADC4 (ZN2, ZP2 o.o.t.)
531       else if(digit.GetSector(0)==4 || digit.GetSector(0)==5){
532         index = (digit.GetSector(0)-4) + 4*digit.GetSector(1); // ZN2 or ZP2
533         lADCDataChannel = 8*(digit.GetSector(0)-4) + digit.GetSector(1);
534         //
535         /*printf("\t AliZDC::Digits2Raw -> idig%d det %d quad %d index %d, ADCch %d ADCVal[%d, %d]\n",
536                 iDigit,digit.GetSector(0),digit.GetSector(1),index,lADCDataChannel,
537                 digit.GetADCValue(0),digit.GetADCValue(1));// Ch. debug
538         */
539         //
540         if(iDigit<knADCData2){ // *** In-time signals
541           lADCDataValue2[index] = digit.GetADCValue(0);
542           if(lADCDataValue2[index] > 2047) lADCDataOvFlw2[index] = 1; 
543           lADCDataValue2[index+2] = digit.GetADCValue(1);
544           if(lADCDataValue2[index+2] > 2047) lADCDataOvFlw2[index+2] = 1; 
545           //
546           lADCData2[index] =   lADCDataGEO << 27 | lADCDataChannel << 17 | 
547                           lADCDataOvFlw2[index] << 12 | (lADCDataValue2[index] & 0xfff); 
548           lADCData2[index+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
549                           lADCDataOvFlw2[index+2] << 12 | (lADCDataValue2[index+2] & 0xfff);   
550         }                 
551         else{ // *** Out-of-time signals
552           lADCDataValue4[index] = digit.GetADCValue(0);
553           if(lADCDataValue4[index] > 2047) lADCDataOvFlw4[index] = 1; 
554           lADCDataValue4[index+2] = digit.GetADCValue(1);
555           if(lADCDataValue4[index+2] > 2047) lADCDataOvFlw4[index+2] = 1; 
556           //
557           lADCData4[index] =   lADCDataGEO << 27 | lADCDataChannel << 17 | 
558                         lADCDataOvFlw4[index] << 12 | (lADCDataValue4[index] & 0xfff); 
559           lADCData4[index+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
560                         lADCDataOvFlw4[index+2] << 12 | (lADCDataValue4[index+2] & 0xfff);   
561         }                 
562       }
563     }
564     // *** ADC2 (Reference PTMs) or ADC4 (Reference PTMs o.o.t.)
565     else if(digit.GetSector(1)==5){
566       index = 20 + (digit.GetSector(0)-1)/3; 
567       lADCDataChannel = 5 + 8*(digit.GetSector(0)-1)/3;
568       //
569       /*printf("\t AliZDC::Digits2Raw -> idig%d det %d quad %d index %d, ADCch %d ADCVal[%d, %d]\n",
570                 iDigit,digit.GetSector(0),digit.GetSector(1),index,lADCDataChannel,
571                 digit.GetADCValue(0),digit.GetADCValue(1));// Ch. debug
572       */
573       //
574       if(iDigit<knADCData2){ // *** In-time signals
575         lADCDataValue2[index] = digit.GetADCValue(0);
576         if(lADCDataValue2[index] > 2047) lADCDataOvFlw2[index] = 1; 
577         lADCDataValue2[index+2] = digit.GetADCValue(1);
578         if(lADCDataValue2[index+2] > 2047) lADCDataOvFlw2[index+2] = 1; 
579         //
580         lADCData2[index] =   lADCDataGEO << 27 | lADCDataChannel << 17 | 
581                         lADCDataOvFlw2[index] << 12 | (lADCDataValue2[index] & 0xfff); 
582         lADCData2[index+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
583                         lADCDataOvFlw2[index+2] << 12 | (lADCDataValue2[index+2] & 0xfff);   
584       }                 
585       else{ // *** Out-of-time signals
586         lADCDataValue4[index] = digit.GetADCValue(0);
587         if(lADCDataValue4[index] > 2047) lADCDataOvFlw4[index] = 1; 
588         lADCDataValue4[index+2] = digit.GetADCValue(1);
589         if(lADCDataValue4[index+2] > 2047) lADCDataOvFlw4[index+2] = 1; 
590         //
591         lADCData4[index] =   lADCDataGEO << 27 | lADCDataChannel << 17 | 
592                       lADCDataOvFlw4[index] << 12 | (lADCDataValue4[index] & 0xfff); 
593         lADCData4[index+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
594                       lADCDataOvFlw4[index+2] << 12 | (lADCDataValue4[index+2] & 0xfff);   
595       }                 
596            
597     }
598     if((index<0) || (index>23)) {
599       Error("Digits2Raw", "sector[0] = %d, sector[1] = %d", 
600             digit.GetSector(0), digit.GetSector(1));
601       continue;
602     }
603     
604     
605   }
606   //
607   /*
608   for(Int_t i=0;i<knADCData1;i++) printf("\t ADCData1[%d] = %x\n",i,lADCData1[i]);
609   for(Int_t i=0;i<knADCData2;i++) printf("\t ADCData2[%d] = %x\n",i,lADCData2[i]);
610   for(Int_t i=0;i<knADCData1;i++) printf("\t ADCData3[%d] = %x\n",i,lADCData3[i]);
611   for(Int_t i=0;i<knADCData2;i++) printf("\t ADCData4[%d] = %x\n",i,lADCData4[i]);
612   */
613  
614   // End of Block
615   UInt_t lADCEndBlockGEO = lADCHeaderGEO;
616   UInt_t lADCEndBlockEvCount = gAlice->GetEventNrInRun();
617   //  
618   lADCEndBlock = lADCEndBlockGEO << 27 | 0x1 << 26 | lADCEndBlockEvCount;
619   //printf("\t AliZDC::Digits2Raw -> ADCEndBlock = %d\n",lADCEndBlock);
620
621
622   // open the output file
623   char fileName[30];
624   strcpy(fileName,AliDAQ::DdlFileName("ZDC",0));
625
626   AliFstream* file = new AliFstream(fileName);
627
628   // write the DDL data header
629   AliRawDataHeaderSim header;
630   header.fSize = sizeof(header) + 
631                  sizeof(lADCHeader1) + sizeof(lADCData1) + sizeof(lADCEndBlock) +
632                  sizeof(lADCHeader2) + sizeof(lADCData2) + sizeof(lADCEndBlock) +
633                  sizeof(lADCHeader1) + sizeof(lADCData3) + sizeof(lADCEndBlock) +
634                  sizeof(lADCHeader2) + sizeof(lADCData4) + sizeof(lADCEndBlock);
635   //
636   /*printf("sizeof header = %d, ADCHeader1 = %d, ADCData1 = %d, ADCEndBlock = %d\n",
637           sizeof(header),sizeof(lADCHeader1),sizeof(lADCData1),sizeof(lADCEndBlock));
638   printf("sizeof header = %d, ADCHeader2 = %d, ADCData2 = %d, ADCEndBlock = %d\n",
639           sizeof(header),sizeof(lADCHeader2),sizeof(lADCData2),sizeof(lADCEndBlock));
640   */
641   //     
642   header.SetAttribute(0);  // valid data
643   file->WriteBuffer((char*)(&header), sizeof(header));
644
645   // write the raw data and close the file
646   file->WriteBuffer((char*) &lADCHeader1, sizeof (lADCHeader1));
647   file->WriteBuffer((char*)(lADCData1), sizeof(lADCData1));
648   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
649   file->WriteBuffer((char*) &lADCHeader2, sizeof (lADCHeader2));
650   file->WriteBuffer((char*)(lADCData2), sizeof(lADCData2));
651   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
652   file->WriteBuffer((char*) &lADCHeader1, sizeof (lADCHeader1));
653   file->WriteBuffer((char*)(lADCData3), sizeof(lADCData3));
654   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
655   file->WriteBuffer((char*) &lADCHeader2, sizeof (lADCHeader2));
656   file->WriteBuffer((char*)(lADCData4), sizeof(lADCData4));
657   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
658   delete file;
659
660   // unload the digits
661   fLoader->UnloadDigits();
662 }
663
664 //_____________________________________________________________________________
665 Bool_t AliZDC::Raw2SDigits(AliRawReader* rawReader)
666 {
667   // Convert ZDC raw data to Sdigits
668   
669   AliLoader* loader = (AliRunLoader::Instance())->GetLoader("ZDCLoader");
670   if(!loader) {
671     AliError("no ZDC loader found");
672     return kFALSE;
673   }
674
675 //  // Event loop
676   Int_t iEvent = 0;
677   while(rawReader->NextEvent()){
678     (AliRunLoader::Instance())->GetEvent(iEvent++);
679     // Create the output digit tree
680     TTree* treeS = loader->TreeS();
681     if(!treeS){
682       loader->MakeTree("S");
683       treeS = loader->TreeS();
684     }
685     //
686     AliZDCSDigit sdigit;
687     AliZDCSDigit* psdigit = &sdigit;
688     const Int_t kBufferSize = 4000;
689     treeS->Branch("ZDC", "AliZDCSDigit",  &psdigit, kBufferSize);
690     //
691     AliZDCRawStream rawStream(rawReader);
692     Int_t sector[2], resADC, rawADC, corrADC, nPheVal;
693     Int_t jcount = 0;
694     while(rawStream.Next()){
695       if(rawStream.IsADCDataWord()){
696         //For the moment only in-time SDigits are foreseen (1st 48 raw values)
697         if(jcount < 48){ 
698           for(Int_t j=0; j<2; j++) sector[j] = rawStream.GetSector(j);
699           rawADC = rawStream.GetADCValue();
700           resADC = rawStream.GetADCGain();
701           //printf("\t RAw2SDigits raw%d ->  RawADC[%d, %d, %d] read\n",
702           //    jcount, sector[0], sector[1], rawADC);
703           //
704           corrADC = rawADC - Pedestal(sector[0], sector[1], resADC);
705           if(corrADC<0) corrADC=0;
706           nPheVal = ADCch2Phe(sector[0], sector[1], corrADC, resADC);
707           //
708           //printf("\t \t ->  SDigit[%d, %d, %d] created\n",
709           //    sector[0], sector[1], nPheVal);
710           //
711           new(psdigit) AliZDCSDigit(sector, (Float_t) nPheVal, 0.);
712           treeS->Fill();
713           jcount++;
714         }
715       }//IsADCDataWord
716     }//rawStream.Next
717     // write the output tree
718     fLoader->WriteSDigits("OVERWRITE");
719     fLoader->UnloadSDigits();
720   }//Event loop 
721    
722   return kTRUE;
723 }
724
725 //_____________________________________________________________________________
726 Int_t AliZDC::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
727 {
728   // Returns a pedestal for detector det, PM quad, channel with res.
729   //
730   // Getting calibration object for ZDC set
731   AliCDBManager *man = AliCDBManager::Instance();
732   AliCDBEntry  *entry = man->Get("ZDC/Calib/Pedestals");
733   AliZDCPedestals *calibPed = (AliZDCPedestals*) entry->GetObject();
734   //
735   if(!calibPed){
736     printf("\t No calibration object found for ZDC!");
737     return -1;
738   }
739   //
740   Int_t index=0, kNch=24;
741   if(Quad!=5){
742     if(Det==1)        index = Quad+kNch*Res;     // ZN1
743     else if(Det==2)   index = Quad+5+kNch*Res;           // ZP1
744     else if(Det==3)   index = Quad+9+kNch*Res; // ZEM
745     else if(Det==4)   index = Quad+12+kNch*Res; // ZN2
746     else if(Det==5)   index = Quad+17+kNch*Res; // ZP2
747   }
748   else index = (Det-1)/3+22+kNch*Res; // Reference PMs
749   //
750   //
751   Float_t meanPed = calibPed->GetMeanPed(index);
752   Float_t pedWidth = calibPed->GetMeanPedWidth(index);
753   Float_t pedValue = gRandom->Gaus(meanPed,pedWidth);
754   //
755   //printf("\t AliZDC::Pedestal - det(%d, %d) - Ped[%d] = %d\n",Det, Quad, index,(Int_t) pedValue); // Chiara debugging!
756   
757   
758
759   return (Int_t) pedValue;
760 }
761
762
763 //_____________________________________________________________________________
764 Int_t AliZDC::ADCch2Phe(Int_t Det, Int_t Quad, Int_t ADCVal, Int_t Res) const
765 {
766   // Evaluation of the no. of phe produced
767   Float_t pmGain[6][5];
768   Float_t resADC[2];
769   for(Int_t j = 0; j < 5; j++){
770     pmGain[0][j] = 50000.;
771     pmGain[1][j] = 100000.;
772     pmGain[2][j] = 100000.;
773     pmGain[3][j] = 50000.;
774     pmGain[4][j] = 100000.;
775     pmGain[5][j] = 100000.;
776   }
777   // ADC Caen V965
778   resADC[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
779   resADC[1] = 0.0000064; // ADC Resolution low gain:  25  fC/adcCh
780   //
781   Int_t nPhe = (Int_t) (ADCVal * pmGain[Det-1][Quad] * resADC[Res]);
782   //
783   //printf("\t AliZDC::ADCch2Phe -> det(%d, %d) - ADC %d  phe %d\n",Det,Quad,ADCVal,nPhe);
784
785   return nPhe;
786 }
787
788 //______________________________________________________________________
789 void AliZDC::SetTreeAddress(){
790
791   // Set branch address for the Trees.
792   if(fLoader->TreeH() && (fHits == 0x0))
793     fHits   = new TClonesArray("AliZDCHit",1000);
794       
795   AliDetector::SetTreeAddress();
796 }