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