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