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