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