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