]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDC.cxx
Time info in hits
[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        Entering 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 pmCZNC=0, pmCZPC=0, pmCZNA=0, pmCZPA=0, pmZEM1 = 0, pmZEM2 = 0;
261     Float_t pmQZNC[4], pmQZPC[4], pmQZNA[4], pmQZPA[4];
262     for(Int_t i = 0; i < 4; i++) pmQZNC[i] = pmQZPC[i] =  pmQZNA[i] = pmQZPA[i] = 0;
263
264     runLoader->GetEvent(iEvent);
265     TTree* treeH = fLoader->TreeH();
266     Int_t ntracks = (Int_t) treeH->GetEntries();
267     ResetHits();
268
269     // Tracks loop
270     Int_t sector[2];
271     for(Int_t itrack = 0; itrack < ntracks; itrack++) {
272       treeH->GetEntry(itrack);
273       for(AliZDCHit* zdcHit = (AliZDCHit*)FirstHit(-1); zdcHit;
274                       zdcHit = (AliZDCHit*)NextHit()) { 
275                       
276         sector[0] = zdcHit->GetVolume(0);
277         sector[1] = zdcHit->GetVolume(1);
278         if((sector[1] < 1) || (sector[1] > 5)) {
279           Error("Hits2SDigits", "sector[0] = %d, sector[1] = %d", 
280                 sector[0], sector[1]);
281           continue;
282         }
283         Float_t lightQ = zdcHit->GetLightPMQ();
284         Float_t lightC = zdcHit->GetLightPMC();
285      
286         if(sector[0] == 1) { //ZNC 
287           pmCZNC += lightC;
288           pmQZNC[sector[1]-1] += lightQ;
289         } 
290         else if(sector[0] == 2) { //ZPC 
291           pmCZPC += lightC;
292           pmQZPC[sector[1]-1] += lightQ;
293         } 
294         else if(sector[0] == 3) { //ZEM 
295           if(sector[1] == 1) pmZEM1 += lightC;
296           else pmZEM2 += lightQ;
297         }
298         if(sector[0] == 4) { //ZNA 
299           pmCZNA += lightC;
300           pmQZNA[sector[1]-1] += lightQ;
301         } 
302         else if(sector[0] == 5) { //ZPA 
303           pmCZPA += lightC;
304           pmQZPA[sector[1]-1] += lightQ;
305         } 
306       }//Hits loop
307     }
308
309     // create the output tree
310     fLoader->MakeTree("S");
311     TTree* treeS = fLoader->TreeS();
312     const Int_t kBufferSize = 4000;
313     treeS->Branch(GetName(), "AliZDCSDigit", &psdigit, kBufferSize);
314
315     // Create sdigits for ZNC
316     sector[0] = 1; // Detector = ZNC
317     sector[1] = 0; // Common PM ADC
318     new(psdigit) AliZDCSDigit(sector, pmCZNC);
319     if(pmCZNC > 0) treeS->Fill();
320     for(Int_t j = 0; j < 4; j++) {
321       sector[1] = j+1; // Towers PM ADCs
322       new(psdigit) AliZDCSDigit(sector, pmQZNC[j]);
323       if(pmQZNC[j] > 0) treeS->Fill();
324     }
325   
326     // Create sdigits for ZPC
327     sector[0] = 2; // Detector = ZPC
328     sector[1] = 0; // Common PM ADC
329     new(psdigit) AliZDCSDigit(sector, pmCZPC);
330     if(pmCZPC > 0) treeS->Fill();
331     for(Int_t j = 0; j < 4; j++) {
332       sector[1] = j+1; // Towers PM ADCs
333       new(psdigit) AliZDCSDigit(sector, pmQZPC[j]);
334       if(pmQZPC[j] > 0) treeS->Fill();
335     }
336
337     // Create sdigits for ZEM
338     sector[0] = 3; 
339     sector[1] = 1; // Detector = ZEM1
340     new(psdigit) AliZDCSDigit(sector, pmZEM1);
341     if(pmZEM1 > 0) treeS->Fill();
342     sector[1] = 2; // Detector = ZEM2
343     new(psdigit) AliZDCSDigit(sector, pmZEM2);
344     if(pmZEM2 > 0) treeS->Fill();
345
346     // Create sdigits for ZNA
347     sector[0] = 4; // Detector = ZNA
348     sector[1] = 0; // Common PM ADC
349     new(psdigit) AliZDCSDigit(sector, pmCZNA);
350     if(pmCZNA > 0) treeS->Fill();
351     for(Int_t j = 0; j < 4; j++) {
352       sector[1] = j+1; // Towers PM ADCs
353       new(psdigit) AliZDCSDigit(sector, pmQZNA[j]);
354       if(pmQZNA[j] > 0) treeS->Fill();
355     }
356   
357     // Create sdigits for ZPA
358     sector[0] = 5; // Detector = ZPA
359     sector[1] = 0; // Common PM ADC
360     new(psdigit) AliZDCSDigit(sector, pmCZPA);
361     if(pmCZPA > 0) treeS->Fill();
362     for(Int_t j = 0; j < 4; j++) {
363       sector[1] = j+1; // Towers PM ADCs
364       new(psdigit) AliZDCSDigit(sector, pmQZPA[j]);
365       if(pmQZPA[j] > 0) treeS->Fill();
366     }
367
368     // write the output tree
369     fLoader->WriteSDigits("OVERWRITE");
370   }
371
372   fLoader->UnloadHits();
373   fLoader->UnloadSDigits();
374 }
375
376 //_____________________________________________________________________________
377 AliDigitizer* AliZDC::CreateDigitizer(AliRunDigitizer* manager) const
378 {
379   // Create the digitizer for ZDC
380
381   return new AliZDCDigitizer(manager);
382 }
383
384 //_____________________________________________________________________________
385 void AliZDC::Digits2Raw()
386 {
387   // Convert ZDC digits to raw data
388
389   // 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
390   //         + 24 int values for the corresponding out of time channels
391   // For the CAEN module V965 we have an Header, the Data Words and an End Of Block
392   //    12 channels x 2 gain chains read from 1st ADC module
393   //    12 channels x 2 gain chains read from 2nd ADC module
394   //    12 channels x 2 gain chains read from 3rd ADC module (o.o.t.)
395   //    12 channels x 2 gain chains read from 4rth ADC module (o.o.t.)
396   //
397   const int knADCData1=24, knADCData2=24; // In principle the 2 numbers can be different!
398   UInt_t lADCHeader1; 
399   UInt_t lADCHeader2; 
400   UInt_t lADCData1[knADCData1];
401   UInt_t lADCData2[knADCData2];
402   UInt_t lADCData3[knADCData1];
403   UInt_t lADCData4[knADCData2];
404   //
405   UInt_t lADCEndBlock;
406
407   // load the digits
408   fLoader->LoadDigits("read");
409   AliZDCDigit digit;
410   AliZDCDigit* pdigit = &digit;
411   TTree* treeD = fLoader->TreeD();
412   if(!treeD) return;
413   treeD->SetBranchAddress("ZDC", &pdigit);
414   //printf("\t AliZDC::Digits2Raw -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
415
416   // Fill data array
417   // ADC header
418   UInt_t lADCHeaderGEO = 0;
419   UInt_t lADCHeaderCRATE = 0;
420   UInt_t lADCHeaderCNT1 = knADCData1;
421   UInt_t lADCHeaderCNT2 = knADCData2;
422     
423   lADCHeader1 = lADCHeaderGEO << 27 | 0x1 << 25 | lADCHeaderCRATE << 16 |
424                lADCHeaderCNT1 << 8 ;
425   lADCHeader2 = lADCHeaderGEO << 27 | 0x1 << 25 | lADCHeaderCRATE << 16 |
426                lADCHeaderCNT2 << 8 ;
427       
428   // ADC data word
429   UInt_t lADCDataGEO = lADCHeaderGEO;
430   //
431   UInt_t lADCDataValue1[knADCData1];
432   UInt_t lADCDataValue2[knADCData2];
433   UInt_t lADCDataValue3[knADCData1];
434   UInt_t lADCDataValue4[knADCData2];
435   //
436   UInt_t lADCDataOvFlw1[knADCData1];
437   UInt_t lADCDataOvFlw2[knADCData2];
438   UInt_t lADCDataOvFlw3[knADCData1];
439   UInt_t lADCDataOvFlw4[knADCData2];
440   //
441   for(Int_t i=0; i<knADCData1 ; i++){
442     lADCDataValue1[i] = 0;
443     lADCDataOvFlw1[i] = 0;
444     lADCDataValue3[i] = 0;
445     lADCDataOvFlw3[i] = 0;
446   }
447   for(Int_t i=0; i<knADCData2 ; i++){
448     lADCDataValue2[i] = 0;
449     lADCDataOvFlw2[i] = 0;
450     lADCDataValue4[i] = 0;
451     lADCDataOvFlw4[i] = 0;
452   }
453   //
454   UInt_t lADCDataChannel = 0;
455   
456   // loop over digits
457   for(Int_t iDigit=0; iDigit<treeD->GetEntries(); iDigit++){
458     treeD->GetEntry(iDigit);
459     if(!pdigit) continue;
460     //digit.Print("");
461     
462     // *** ADC data
463     Int_t index=0;
464     if(digit.GetSector(1)!=5){ // ZDC signal channels
465       // *** ADC1 (ZN1, ZP1, ZEM1,2) or ADC3 (ZN1, ZP1, ZEM1,2 o.o.t.)
466       if(digit.GetSector(0)==1 || digit.GetSector(0)==2 || digit.GetSector(0)==3){
467         if(digit.GetSector(0)==1 || digit.GetSector(0)==2){
468           index = (digit.GetSector(0)-1) + 4*digit.GetSector(1); // ZN1 or ZP1
469           lADCDataChannel = 8*(digit.GetSector(0)-1) + digit.GetSector(1);
470         }
471         else if(digit.GetSector(0)==3){ // ZEM 1,2
472           index = 20 + (digit.GetSector(1)-1);
473           lADCDataChannel = 5 + 8*(digit.GetSector(1)-1);
474         }
475         //
476         /*printf("\t AliZDC::Digits2Raw -> idig%d det %d quad %d index %d, ADCch %d ADCVal[%d, %d]\n",
477                 iDigit,digit.GetSector(0),digit.GetSector(1),index,lADCDataChannel,
478                 digit.GetADCValue(0),digit.GetADCValue(1));// Ch. debug
479         */
480         //
481         if(iDigit<knADCData1){ // *** In-time signals
482           lADCDataValue1[index] = digit.GetADCValue(0);   // High gain ADC ch.    
483           if(lADCDataValue1[index] > 2047) lADCDataOvFlw1[index] = 1; 
484           lADCDataValue1[index+2] = digit.GetADCValue(1); // Low gain ADC ch.
485           if(lADCDataValue1[index+2] > 2047) lADCDataOvFlw1[index+2] = 1; 
486         
487           lADCData1[index] = lADCDataGEO << 27 | 0x1 << 24 | lADCDataChannel << 17 | 
488                           lADCDataOvFlw1[index] << 12 | (lADCDataValue1[index] & 0xfff); 
489           lADCData1[index+2] = lADCDataGEO << 27 | 0x1 << 24  | lADCDataChannel << 17 | 0x1 << 16 |
490                           lADCDataOvFlw1[index+2] << 12 | (lADCDataValue1[index+2] & 0xfff);  
491         }
492         else{ // *** Out-of-time signals
493           lADCDataValue3[index] = digit.GetADCValue(0);   // High gain ADC ch.    
494           if(lADCDataValue3[index] > 2047) lADCDataOvFlw3[index] = 1; 
495           lADCDataValue3[index+2] = digit.GetADCValue(1); // Low gain ADC ch.
496           if(lADCDataValue3[index+2] > 2047) lADCDataOvFlw3[index+2] = 1; 
497       
498           lADCData3[index] = lADCDataGEO << 27 | lADCDataChannel << 17 | 
499                           lADCDataOvFlw3[index] << 12 | (lADCDataValue3[index] & 0xfff); 
500           lADCData3[index+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
501                           lADCDataOvFlw3[index+2] << 12 | (lADCDataValue3[index+2] & 0xfff);  
502         }                  
503       }
504       // *** ADC2 (ZN2, ZP2) or ADC4 (ZN2, ZP2 o.o.t.)
505       else if(digit.GetSector(0)==4 || digit.GetSector(0)==5){
506         index = (digit.GetSector(0)-4) + 4*digit.GetSector(1); // ZN2 or ZP2
507         lADCDataChannel = 8*(digit.GetSector(0)-4) + digit.GetSector(1);
508         //
509         /*printf("\t AliZDC::Digits2Raw -> idig%d det %d quad %d index %d, ADCch %d ADCVal[%d, %d]\n",
510                 iDigit,digit.GetSector(0),digit.GetSector(1),index,lADCDataChannel,
511                 digit.GetADCValue(0),digit.GetADCValue(1));// Ch. debug
512         */
513         //
514         if(iDigit<knADCData2){ // *** In-time signals
515           lADCDataValue2[index] = digit.GetADCValue(0);
516           if(lADCDataValue2[index] > 2047) lADCDataOvFlw2[index] = 1; 
517           lADCDataValue2[index+2] = digit.GetADCValue(1);
518           if(lADCDataValue2[index+2] > 2047) lADCDataOvFlw2[index+2] = 1; 
519           //
520           lADCData2[index] =   lADCDataGEO << 27 | lADCDataChannel << 17 | 
521                           lADCDataOvFlw2[index] << 12 | (lADCDataValue2[index] & 0xfff); 
522           lADCData2[index+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
523                           lADCDataOvFlw2[index+2] << 12 | (lADCDataValue2[index+2] & 0xfff);   
524         }                 
525         else{ // *** Out-of-time signals
526           lADCDataValue4[index] = digit.GetADCValue(0);
527           if(lADCDataValue4[index] > 2047) lADCDataOvFlw4[index] = 1; 
528           lADCDataValue4[index+2] = digit.GetADCValue(1);
529           if(lADCDataValue4[index+2] > 2047) lADCDataOvFlw4[index+2] = 1; 
530           //
531           lADCData4[index] =   lADCDataGEO << 27 | lADCDataChannel << 17 | 
532                         lADCDataOvFlw4[index] << 12 | (lADCDataValue4[index] & 0xfff); 
533           lADCData4[index+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
534                         lADCDataOvFlw4[index+2] << 12 | (lADCDataValue4[index+2] & 0xfff);   
535         }                 
536       }
537     }
538     // *** ADC2 (Reference PTMs) or ADC4 (Reference PTMs o.o.t.)
539     else if(digit.GetSector(1)==5){
540       index = 20 + (digit.GetSector(0)-1)/3; 
541       lADCDataChannel = 5 + 8*(digit.GetSector(0)-1)/3;
542       //
543       /*printf("\t AliZDC::Digits2Raw -> idig%d det %d quad %d index %d, ADCch %d ADCVal[%d, %d]\n",
544                 iDigit,digit.GetSector(0),digit.GetSector(1),index,lADCDataChannel,
545                 digit.GetADCValue(0),digit.GetADCValue(1));// Ch. debug
546       */
547       //
548       if(iDigit<knADCData2){ // *** In-time signals
549         lADCDataValue2[index] = digit.GetADCValue(0);
550         if(lADCDataValue2[index] > 2047) lADCDataOvFlw2[index] = 1; 
551         lADCDataValue2[index+2] = digit.GetADCValue(1);
552         if(lADCDataValue2[index+2] > 2047) lADCDataOvFlw2[index+2] = 1; 
553         //
554         lADCData2[index] =   lADCDataGEO << 27 | lADCDataChannel << 17 | 
555                         lADCDataOvFlw2[index] << 12 | (lADCDataValue2[index] & 0xfff); 
556         lADCData2[index+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
557                         lADCDataOvFlw2[index+2] << 12 | (lADCDataValue2[index+2] & 0xfff);   
558       }                 
559       else{ // *** Out-of-time signals
560         lADCDataValue4[index] = digit.GetADCValue(0);
561         if(lADCDataValue4[index] > 2047) lADCDataOvFlw4[index] = 1; 
562         lADCDataValue4[index+2] = digit.GetADCValue(1);
563         if(lADCDataValue4[index+2] > 2047) lADCDataOvFlw4[index+2] = 1; 
564         //
565         lADCData4[index] =   lADCDataGEO << 27 | lADCDataChannel << 17 | 
566                       lADCDataOvFlw4[index] << 12 | (lADCDataValue4[index] & 0xfff); 
567         lADCData4[index+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
568                       lADCDataOvFlw4[index+2] << 12 | (lADCDataValue4[index+2] & 0xfff);   
569       }                 
570            
571     }
572     if((index<0) || (index>23)) {
573       Error("Digits2Raw", "sector[0] = %d, sector[1] = %d", 
574             digit.GetSector(0), digit.GetSector(1));
575       continue;
576     }
577     
578     
579   }
580   //
581   /*
582   for(Int_t i=0;i<knADCData1;i++) printf("\t ADCData1[%d] = %x\n",i,lADCData1[i]);
583   for(Int_t i=0;i<knADCData2;i++) printf("\t ADCData2[%d] = %x\n",i,lADCData2[i]);
584   for(Int_t i=0;i<knADCData1;i++) printf("\t ADCData3[%d] = %x\n",i,lADCData3[i]);
585   for(Int_t i=0;i<knADCData2;i++) printf("\t ADCData4[%d] = %x\n",i,lADCData4[i]);
586   */
587  
588   // End of Block
589   UInt_t lADCEndBlockGEO = lADCHeaderGEO;
590   UInt_t lADCEndBlockEvCount = gAlice->GetEventNrInRun();
591   //  
592   lADCEndBlock = lADCEndBlockGEO << 27 | 0x1 << 26 | lADCEndBlockEvCount;
593   //printf("\t AliZDC::Digits2Raw -> ADCEndBlock = %d\n",lADCEndBlock);
594
595
596   // open the output file
597   char fileName[30];
598   strcpy(fileName,AliDAQ::DdlFileName("ZDC",0));
599
600   AliFstream* file = new AliFstream(fileName);
601
602   // write the DDL data header
603   AliRawDataHeaderSim header;
604   header.fSize = sizeof(header) + 
605                  sizeof(lADCHeader1) + sizeof(lADCData1) + sizeof(lADCEndBlock) +
606                  sizeof(lADCHeader2) + sizeof(lADCData2) + sizeof(lADCEndBlock) +
607                  sizeof(lADCHeader1) + sizeof(lADCData3) + sizeof(lADCEndBlock) +
608                  sizeof(lADCHeader2) + sizeof(lADCData4) + sizeof(lADCEndBlock);
609   //
610   /*printf("sizeof header = %d, ADCHeader1 = %d, ADCData1 = %d, ADCEndBlock = %d\n",
611           sizeof(header),sizeof(lADCHeader1),sizeof(lADCData1),sizeof(lADCEndBlock));
612   printf("sizeof header = %d, ADCHeader2 = %d, ADCData2 = %d, ADCEndBlock = %d\n",
613           sizeof(header),sizeof(lADCHeader2),sizeof(lADCData2),sizeof(lADCEndBlock));
614   */
615   //     
616   header.SetAttribute(0);  // valid data
617   file->WriteBuffer((char*)(&header), sizeof(header));
618
619   // write the raw data and close the file
620   file->WriteBuffer((char*) &lADCHeader1, sizeof (lADCHeader1));
621   file->WriteBuffer((char*)(lADCData1), sizeof(lADCData1));
622   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
623   file->WriteBuffer((char*) &lADCHeader2, sizeof (lADCHeader2));
624   file->WriteBuffer((char*)(lADCData2), sizeof(lADCData2));
625   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
626   file->WriteBuffer((char*) &lADCHeader1, sizeof (lADCHeader1));
627   file->WriteBuffer((char*)(lADCData3), sizeof(lADCData3));
628   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
629   file->WriteBuffer((char*) &lADCHeader2, sizeof (lADCHeader2));
630   file->WriteBuffer((char*)(lADCData4), sizeof(lADCData4));
631   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
632   delete file;
633
634   // unload the digits
635   fLoader->UnloadDigits();
636 }
637
638 //_____________________________________________________________________________
639 Bool_t AliZDC::Raw2SDigits(AliRawReader* rawReader)
640 {
641   // Convert ZDC raw data to Sdigits
642   
643   AliLoader* loader = (AliRunLoader::GetRunLoader())->GetLoader("ZDCLoader");
644   if(!loader) {
645     AliError("no ZDC loader found");
646     return kFALSE;
647   }
648
649 //  // Event loop
650   Int_t iEvent = 0;
651   while(rawReader->NextEvent()){
652     (AliRunLoader::GetRunLoader())->GetEvent(iEvent++);
653     // Create the output digit tree
654     TTree* treeS = loader->TreeS();
655     if(!treeS){
656       loader->MakeTree("S");
657       treeS = loader->TreeS();
658     }
659     //
660     AliZDCSDigit sdigit;
661     AliZDCSDigit* psdigit = &sdigit;
662     const Int_t kBufferSize = 4000;
663     treeS->Branch("ZDC", "AliZDCSDigit",  &psdigit, kBufferSize);
664     //
665     AliZDCRawStream rawStream(rawReader);
666     Int_t sector[2], resADC, rawADC, corrADC, nPheVal;
667     Int_t jcount = 0;
668     while(rawStream.Next()){
669       if(rawStream.IsADCDataWord()){
670         //For the moment only in-time SDigits are foreseen (1st 48 raw values)
671         if(jcount < 48){ 
672           for(Int_t j=0; j<2; j++) sector[j] = rawStream.GetSector(j);
673           rawADC = rawStream.GetADCValue();
674           resADC = rawStream.GetADCGain();
675           //printf("\t RAw2SDigits raw%d ->  RawADC[%d, %d, %d] read\n",
676           //    jcount, sector[0], sector[1], rawADC);
677           //
678           corrADC = rawADC - Pedestal(sector[0], sector[1], resADC);
679           if(corrADC<0) corrADC=0;
680           nPheVal = ADCch2Phe(sector[0], sector[1], corrADC, resADC);
681           //
682           //printf("\t \t ->  SDigit[%d, %d, %d] created\n",
683           //    sector[0], sector[1], nPheVal);
684           //
685           new(psdigit) AliZDCSDigit(sector, (Float_t) nPheVal);
686           treeS->Fill();
687           jcount++;
688         }
689       }//IsADCDataWord
690     }//rawStream.Next
691     // write the output tree
692     fLoader->WriteSDigits("OVERWRITE");
693     fLoader->UnloadSDigits();
694   }//Event loop 
695    
696   return kTRUE;
697 }
698
699 //_____________________________________________________________________________
700 Int_t AliZDC::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
701 {
702   // Returns a pedestal for detector det, PM quad, channel with res.
703   //
704   // Getting calibration object for ZDC set
705   AliCDBManager *man = AliCDBManager::Instance();
706   AliCDBEntry  *entry = man->Get("ZDC/Calib/Pedestals");
707   AliZDCPedestals *calibPed = (AliZDCPedestals*) entry->GetObject();
708   //
709   if(!calibPed){
710     printf("\t No calibration object found for ZDC!");
711     return -1;
712   }
713   //
714   Int_t index=0, kNch=24;
715   if(Quad!=5){
716     if(Det==1)        index = Quad+kNch*Res;     // ZN1
717     else if(Det==2)   index = Quad+5+kNch*Res;           // ZP1
718     else if(Det==3)   index = Quad+9+kNch*Res; // ZEM
719     else if(Det==4)   index = Quad+12+kNch*Res; // ZN2
720     else if(Det==5)   index = Quad+17+kNch*Res; // ZP2
721   }
722   else index = (Det-1)/3+22+kNch*Res; // Reference PMs
723   //
724   //
725   Float_t meanPed = calibPed->GetMeanPed(index);
726   Float_t pedWidth = calibPed->GetMeanPedWidth(index);
727   Float_t pedValue = gRandom->Gaus(meanPed,pedWidth);
728   //
729   //printf("\t AliZDC::Pedestal - det(%d, %d) - Ped[%d] = %d\n",Det, Quad, index,(Int_t) pedValue); // Chiara debugging!
730   
731   
732
733   return (Int_t) pedValue;
734 }
735
736
737 //_____________________________________________________________________________
738 Int_t AliZDC::ADCch2Phe(Int_t Det, Int_t Quad, Int_t ADCVal, Int_t Res) const
739 {
740   // Evaluation of the no. of phe produced
741   Float_t pmGain[6][5];
742   Float_t resADC[2];
743   for(Int_t j = 0; j < 5; j++){
744     pmGain[0][j] = 50000.;
745     pmGain[1][j] = 100000.;
746     pmGain[2][j] = 100000.;
747     pmGain[3][j] = 50000.;
748     pmGain[4][j] = 100000.;
749     pmGain[5][j] = 100000.;
750   }
751   // ADC Caen V965
752   resADC[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
753   resADC[1] = 0.0000064; // ADC Resolution low gain:  25  fC/adcCh
754   //
755   Int_t nPhe = (Int_t) (ADCVal * pmGain[Det-1][Quad] * resADC[Res]);
756   //
757   //printf("\t AliZDC::ADCch2Phe -> det(%d, %d) - ADC %d  phe %d\n",Det,Quad,ADCVal,nPhe);
758
759   return nPhe;
760 }
761
762 //______________________________________________________________________
763 void AliZDC::SetTreeAddress(){
764
765   // Set branch address for the Trees.
766   if(fLoader->TreeH() && (fHits == 0x0))
767     fHits   = new TClonesArray("AliZDCHit",1000);
768       
769   AliDetector::SetTreeAddress();
770 }