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