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