]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDC.cxx
M
[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 <TBRIK.h>
29 #include <TClonesArray.h>
30 #include <TGeometry.h>
31 #include <TNode.h>
32 #include <TTree.h>
33 #include <TFile.h>
34 #include <TSystem.h>
35 #include <TRandom.h>
36
37 // --- AliRoot header files
38 #include "AliDetector.h"
39 #include "AliRawDataHeaderSim.h"
40 #include "AliRawReader.h"
41 #include "AliLoader.h"
42 #include "AliRun.h"
43 #include "AliMC.h"
44 #include "AliLog.h"
45 #include "AliDAQ.h"
46 #include "AliZDC.h"
47 #include "AliZDCHit.h"
48 #include "AliZDCSDigit.h"
49 #include "AliZDCDigit.h"
50 #include "AliZDCDigitizer.h"
51 #include "AliZDCRawStream.h"
52 #include "AliZDCPedestals.h"
53 #include "AliZDCCalib.h"
54 #include "AliFstream.h"
55
56  
57 ClassImp(AliZDC)
58
59 //_____________________________________________________________________________
60 AliZDC::AliZDC() :
61   AliDetector(),
62   fNoShower(0),
63   fPedCalib(0),
64   fCalibData(0),
65   fZDCCalibFName("")
66 {
67   //
68   // Default constructor for the Zero Degree Calorimeter base class
69   //
70   
71   fIshunt = 1;
72   fNhits  = 0;
73   fHits = 0;
74   fDigits = 0;
75   fNdigits = 0;
76   
77 }
78  
79 //_____________________________________________________________________________
80 AliZDC::AliZDC(const char *name, const char *title) : 
81   AliDetector(name,title),
82   fNoShower  (0),
83   fPedCalib(0),
84   fCalibData(0),
85   fZDCCalibFName("")
86 {
87   //
88   // Standard constructor for the Zero Degree Calorimeter base class
89   //
90     
91   fIshunt = 1;
92   fNhits  = 0;
93   fDigits = 0;
94   fNdigits = 0;
95  
96   fHits = new TClonesArray("AliZDCHit",1000);
97   gAlice->GetMCApp()->AddHitList(fHits);
98   
99   char sensname[5],senstitle[25];
100   sprintf(sensname,"ZDC");
101   sprintf(senstitle,"ZDC dummy");
102   SetName(sensname); SetTitle(senstitle);
103
104 }
105
106 //____________________________________________________________________________ 
107 AliZDC::~AliZDC()
108 {
109   //
110   // ZDC destructor
111   //
112
113   fIshunt = 0;
114   delete fPedCalib;
115   delete fCalibData;
116
117 }
118
119 //_____________________________________________________________________________
120 AliZDC::AliZDC(const AliZDC& ZDC) :
121 AliDetector("ZDC","ZDC"),
122 fNoShower(ZDC.fNoShower),
123 fPedCalib(ZDC.fPedCalib),
124 fCalibData(ZDC.fCalibData),
125 fZDCCalibFName(ZDC.fZDCCalibFName)
126 {
127   // copy constructor
128 }
129
130 //_____________________________________________________________________________
131 AliZDC& AliZDC::operator=(const AliZDC& ZDC)
132 {
133   // assignement operator
134   if(this!=&ZDC){
135     fNoShower = ZDC.fNoShower;
136     fPedCalib = ZDC.fPedCalib;
137     fCalibData = ZDC.fCalibData;
138     fZDCCalibFName = ZDC.fZDCCalibFName;
139   } return *this;
140 }
141
142 //_____________________________________________________________________________
143 void AliZDC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
144 {
145   //
146   //            Add a ZDC hit to the hit list.
147   // -> We make use of 2 array of hits:
148   // [1]  fHits (the usual one) that contains hits for each PRIMARY
149   // [2]  fStHits that contains hits for each EVENT and is used to
150   //      obtain digits at the end of each event
151   //
152   
153   static Float_t primKinEn=0., xImpact=0., yImpact=0., sFlag=0.;
154   static Int_t   pcPDGcode;
155
156   AliZDCHit *newquad, *curprimquad;
157   newquad = new AliZDCHit(fIshunt, track, vol, hits);
158   TClonesArray &lhits = *fHits;
159   
160   if(fNhits==0){
161       // First hit -> setting flag for primary or secondary particle
162       Int_t primary = gAlice->GetMCApp()->GetPrimary(track);     
163       //
164       if(track != primary){
165         newquad->SetSFlag(1);  // SECONDARY particle entering the ZDC
166       }
167       else if(track == primary){
168         newquad->SetSFlag(0);  // PRIMARY particle entering the ZDC
169       }  
170       sFlag     = newquad->GetSFlag();
171       primKinEn = newquad->GetPrimKinEn();
172       xImpact   = newquad->GetXImpact();
173       yImpact   = newquad->GetYImpact();
174       pcPDGcode = newquad->GetPDGCode();
175    }
176    else{       
177       newquad->SetPrimKinEn(primKinEn);
178       newquad->SetXImpact(xImpact);
179       newquad->SetYImpact(yImpact);
180       newquad->SetSFlag(sFlag);
181       newquad->SetPDGCode(pcPDGcode);
182    }
183  
184   Int_t j;
185   for(j=0; j<fNhits; j++){
186     // If hits are equal (same track, same volume), sum them.
187      curprimquad = (AliZDCHit*) lhits[j];
188      if(*curprimquad == *newquad){
189         *curprimquad = *curprimquad+*newquad;
190         // CH. debug
191         /*if(newquad->GetEnergy() != 0. || newquad->GetLightPMC() != 0. || 
192            newquad->GetLightPMQ() != 0.){
193           printf("\n\t --- Equal hits found\n");
194           curprimquad->Print("");
195           newquad->Print("");
196           printf("\t --- Det. %d, Quad. %d: X = %f, E = %f, LightPMC = %f, LightPMQ = %f\n",
197           curprimquad->GetVolume(0),curprimquad->GetVolume(1),curprimquad->GetXImpact(),
198           curprimquad->GetEnergy(), curprimquad->GetLightPMC(), curprimquad->GetLightPMQ());
199         }*/
200         //
201         delete newquad;
202         return;
203      } 
204   }
205
206     //Otherwise create a new hit
207     new(lhits[fNhits]) AliZDCHit(*newquad);
208     fNhits++;
209     // CH. debug
210     /*printf("\n\t New ZDC hit added! fNhits = %d\n", fNhits);
211     printf("\t Det. %d, Quad.t %d: X = %f, E = %f, LightPMC = %f, LightPMQ = %f\n",
212     newquad->GetVolume(0),newquad->GetVolume(1),newquad->GetXImpact(),
213     newquad->GetEnergy(), newquad->GetLightPMC(), newquad->GetLightPMQ());
214     */
215     delete newquad;
216 }
217
218 //_____________________________________________________________________________
219 void AliZDC::BuildGeometry()
220 {
221   //
222   // Build the ROOT TNode geometry for event display 
223   // in the Zero Degree Calorimeter
224   // This routine is dummy for the moment
225   //
226
227   TNode *node, *top;
228   TBRIK *brik;
229   const int kColorZDC  = kBlue;
230   
231   //
232   top=gAlice->GetGeometry()->GetNode("alice");
233   
234   // ZDC
235   brik = new TBRIK("S_ZDC","ZDC box","void",300,300,5);
236   top->cd();
237   node = new TNode("ZDC","ZDC","S_ZDC",0,0,600,"");
238   node->SetLineColor(kColorZDC);
239   fNodes->Add(node);
240 }
241
242 //____________________________________________________________________________
243 Float_t AliZDC::ZMin(void) const
244 {
245   // Minimum dimension of the ZDC module in z
246   return -11600.;
247 }
248
249 //____________________________________________________________________________
250 Float_t AliZDC::ZMax(void) const
251 {
252   // Maximum dimension of the ZDC module in z
253   return 11750.;
254 }
255   
256
257 //_____________________________________________________________________________
258 void AliZDC::MakeBranch(Option_t *opt)
259 {
260   //
261   // Create Tree branches for the ZDC
262   //
263
264   char branchname[10];
265   sprintf(branchname,"%s",GetName());
266
267   const char *cH = strstr(opt,"H");
268   
269   if(cH && fLoader->TreeH()) {
270     if (fHits) {
271       fHits->Clear();
272       fNhits = 0;
273     }
274     else {
275       fHits   = new TClonesArray("AliZDCHit",1000); 
276       if (gAlice && gAlice->GetMCApp())
277         gAlice->GetMCApp()->AddHitList(fHits);
278     }
279   }
280   
281   AliDetector::MakeBranch(opt);
282 }
283
284 //_____________________________________________________________________________
285 void AliZDC::Hits2SDigits()
286 {
287   // Create summable digits from hits
288   
289   AliDebug(1,"\n        Entering AliZDC::Hits2SDigits() ");
290   
291   fLoader->LoadHits("read");
292   fLoader->LoadSDigits("recreate");
293   AliRunLoader* runLoader = fLoader->GetRunLoader();
294   AliZDCSDigit sdigit;
295   AliZDCSDigit* psdigit = &sdigit;
296
297   // Event loop
298   for(Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
299     Float_t pmCZNC=0, pmCZPC=0, pmCZNA=0, pmCZPA=0, pmZEM1 = 0, pmZEM2 = 0;
300     Float_t pmQZNC[4], pmQZPC[4], pmQZNA[4], pmQZPA[4];
301     for(Int_t i = 0; i < 4; i++) pmQZNC[i] = pmQZPC[i] =  pmQZNA[i] = pmQZPA[i] = 0;
302
303     runLoader->GetEvent(iEvent);
304     TTree* treeH = fLoader->TreeH();
305     Int_t ntracks = (Int_t) treeH->GetEntries();
306     ResetHits();
307
308     // Tracks loop
309     Int_t sector[2];
310     for(Int_t itrack = 0; itrack < ntracks; itrack++) {
311       treeH->GetEntry(itrack);
312       for(AliZDCHit* zdcHit = (AliZDCHit*)FirstHit(-1); zdcHit;
313                       zdcHit = (AliZDCHit*)NextHit()) { 
314                       
315         sector[0] = zdcHit->GetVolume(0);
316         sector[1] = zdcHit->GetVolume(1);
317         if((sector[1] < 1) || (sector[1] > 5)) {
318           Error("Hits2SDigits", "sector[0] = %d, sector[1] = %d", 
319                 sector[0], sector[1]);
320           continue;
321         }
322         Float_t lightQ = zdcHit->GetLightPMQ();
323         Float_t lightC = zdcHit->GetLightPMC();
324      
325         if(sector[0] == 1) { //ZNC 
326           pmCZNC += lightC;
327           pmQZNC[sector[1]-1] += lightQ;
328         } 
329         else if(sector[0] == 2) { //ZPC 
330           pmCZPC += lightC;
331           pmQZPC[sector[1]-1] += lightQ;
332         } 
333         else if(sector[0] == 3) { //ZEM 
334           if(sector[1] == 1) pmZEM1 += lightC;
335           else pmZEM2 += lightQ;
336         }
337         if(sector[0] == 4) { //ZNA 
338           pmCZNA += lightC;
339           pmQZNA[sector[1]-1] += lightQ;
340         } 
341         else if(sector[0] == 5) { //ZPA 
342           pmCZPA += lightC;
343           pmQZPA[sector[1]-1] += lightQ;
344         } 
345       }//Hits loop
346     }
347
348     // create the output tree
349     fLoader->MakeTree("S");
350     TTree* treeS = fLoader->TreeS();
351     const Int_t kBufferSize = 4000;
352     treeS->Branch(GetName(), "AliZDCSDigit", &psdigit, kBufferSize);
353
354     // Create sdigits for ZNC
355     sector[0] = 1; // Detector = ZNC
356     sector[1] = 0; // Common PM ADC
357     new(psdigit) AliZDCSDigit(sector, pmCZNC);
358     if(pmCZNC > 0) treeS->Fill();
359     for(Int_t j = 0; j < 4; j++) {
360       sector[1] = j+1; // Towers PM ADCs
361       new(psdigit) AliZDCSDigit(sector, pmQZNC[j]);
362       if(pmQZNC[j] > 0) treeS->Fill();
363     }
364   
365     // Create sdigits for ZPC
366     sector[0] = 2; // Detector = ZPC
367     sector[1] = 0; // Common PM ADC
368     new(psdigit) AliZDCSDigit(sector, pmCZPC);
369     if(pmCZPC > 0) treeS->Fill();
370     for(Int_t j = 0; j < 4; j++) {
371       sector[1] = j+1; // Towers PM ADCs
372       new(psdigit) AliZDCSDigit(sector, pmQZPC[j]);
373       if(pmQZPC[j] > 0) treeS->Fill();
374     }
375
376     // Create sdigits for ZEM
377     sector[0] = 3; 
378     sector[1] = 1; // Detector = ZEM1
379     new(psdigit) AliZDCSDigit(sector, pmZEM1);
380     if(pmZEM1 > 0) treeS->Fill();
381     sector[1] = 2; // Detector = ZEM2
382     new(psdigit) AliZDCSDigit(sector, pmZEM2);
383     if(pmZEM2 > 0) treeS->Fill();
384
385     // Create sdigits for ZNA
386     sector[0] = 4; // Detector = ZNA
387     sector[1] = 0; // Common PM ADC
388     new(psdigit) AliZDCSDigit(sector, pmCZNA);
389     if(pmCZNA > 0) treeS->Fill();
390     for(Int_t j = 0; j < 4; j++) {
391       sector[1] = j+1; // Towers PM ADCs
392       new(psdigit) AliZDCSDigit(sector, pmQZNA[j]);
393       if(pmQZNA[j] > 0) treeS->Fill();
394     }
395   
396     // Create sdigits for ZPA
397     sector[0] = 5; // Detector = ZPA
398     sector[1] = 0; // Common PM ADC
399     new(psdigit) AliZDCSDigit(sector, pmCZPA);
400     if(pmCZPA > 0) treeS->Fill();
401     for(Int_t j = 0; j < 4; j++) {
402       sector[1] = j+1; // Towers PM ADCs
403       new(psdigit) AliZDCSDigit(sector, pmQZPA[j]);
404       if(pmQZPA[j] > 0) treeS->Fill();
405     }
406
407     // write the output tree
408     fLoader->WriteSDigits("OVERWRITE");
409   }
410
411   fLoader->UnloadHits();
412   fLoader->UnloadSDigits();
413 }
414
415 //_____________________________________________________________________________
416 AliDigitizer* AliZDC::CreateDigitizer(AliRunDigitizer* manager) const
417 {
418   // Create the digitizer for ZDC
419
420   return new AliZDCDigitizer(manager);
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   //
620   /*
621   for(Int_t i=0;i<knADCData1;i++) printf("\t ADCData1[%d] = %x\n",i,lADCData1[i]);
622   for(Int_t i=0;i<knADCData2;i++) printf("\t ADCData2[%d] = %x\n",i,lADCData2[i]);
623   for(Int_t i=0;i<knADCData1;i++) printf("\t ADCData3[%d] = %x\n",i,lADCData3[i]);
624   for(Int_t i=0;i<knADCData2;i++) printf("\t ADCData4[%d] = %x\n",i,lADCData4[i]);
625   */
626  
627   // End of Block
628   UInt_t lADCEndBlockGEO = lADCHeaderGEO;
629   UInt_t lADCEndBlockEvCount = gAlice->GetEventNrInRun();
630   //  
631   lADCEndBlock = lADCEndBlockGEO << 27 | 0x1 << 26 | lADCEndBlockEvCount;
632   //printf("\t AliZDC::Digits2Raw -> ADCEndBlock = %d\n",lADCEndBlock);
633
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 = (gAlice->GetRunLoader())->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     (gAlice->GetRunLoader())->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);
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 }