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