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