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