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