]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDC.cxx
Hits, SDigits, Digits updated to AliZDCv3 (2nd ZDC set included)
[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 "AliZDCRecParam.h"
55 #include "AliFstream.h"
56
57  
58 ClassImp(AliZDC)
59
60 //_____________________________________________________________________________
61 AliZDC::AliZDC() :
62   AliDetector(),
63   fNoShower  (0),
64   fPedCalib(0),
65   fCalibData(0),
66   fRecParam(0)
67 {
68   //
69   // Default constructor for the Zero Degree Calorimeter base class
70   //
71   
72   fIshunt = 1;
73   fNhits  = 0;
74   fHits = 0;
75   fDigits = 0;
76   fNdigits = 0;
77   
78 }
79  
80 //_____________________________________________________________________________
81 AliZDC::AliZDC(const char *name, const char *title) : 
82   AliDetector(name,title),
83   fNoShower  (0),
84   fPedCalib(0),
85   fCalibData(0),
86   fRecParam(0)
87 {
88   //
89   // Standard constructor for the Zero Degree Calorimeter base class
90   //
91     
92   fIshunt = 1;
93   fNhits  = 0;
94   fDigits = 0;
95   fNdigits = 0;
96  
97   fHits = new TClonesArray("AliZDCHit",1000);
98   gAlice->GetMCApp()->AddHitList(fHits);
99   
100   char sensname[5],senstitle[25];
101   sprintf(sensname,"ZDC");
102   sprintf(senstitle,"ZDC dummy");
103   SetName(sensname); SetTitle(senstitle);
104
105 }
106
107 //____________________________________________________________________________ 
108 AliZDC::~AliZDC()
109 {
110   //
111   // ZDC destructor
112   //
113
114   fIshunt = 0;
115   delete fPedCalib;
116   delete fCalibData;
117   delete fRecParam;
118
119 }
120
121 //_____________________________________________________________________________
122 AliZDC::AliZDC(const AliZDC& ZDC) :
123   AliDetector("ZDC","ZDC")
124 {
125   // copy constructor
126     fNoShower = ZDC.fNoShower;
127     fPedCalib = ZDC.fPedCalib;
128     fCalibData = ZDC.fCalibData;
129     fRecParam = ZDC.fRecParam;
130     fZDCCalibFName = ZDC.fZDCCalibFName;
131 }
132
133 //_____________________________________________________________________________
134 AliZDC& AliZDC::operator=(const AliZDC& ZDC)
135 {
136   // assignement operator
137   if(this!=&ZDC){
138     fNoShower = ZDC.fNoShower;
139     fPedCalib = ZDC.fPedCalib;
140     fCalibData = ZDC.fCalibData;
141     fRecParam = ZDC.fRecParam;
142     fZDCCalibFName = ZDC.fZDCCalibFName;
143   } return *this;
144 }
145
146 //_____________________________________________________________________________
147 void AliZDC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
148 {
149   //
150   //            Add a ZDC hit to the hit list.
151   // -> We make use of 2 array of hits:
152   // [1]  fHits (the usual one) that contains hits for each PRIMARY
153   // [2]  fStHits that contains hits for each EVENT and is used to
154   //      obtain digits at the end of each event
155   //
156   
157   static Float_t primKinEn, xImpact, yImpact, sFlag;
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       Int_t primary = gAlice->GetMCApp()->GetPrimary(track);     
166       if(track != primary){
167         newquad->SetSFlag(1);  // SECONDARY particle entering the ZDC
168       }
169       else if(track == primary){
170         newquad->SetSFlag(0);  // PRIMARY particle entering the ZDC
171       }  
172       sFlag     = newquad->GetSFlag();
173       primKinEn = newquad->GetPrimKinEn();
174       xImpact   = newquad->GetXImpact();
175       yImpact   = newquad->GetYImpact();
176    }
177    else{       
178       newquad->SetPrimKinEn(primKinEn);
179       newquad->SetXImpact(xImpact);
180       newquad->SetYImpact(yImpact);
181       newquad->SetSFlag(sFlag);
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    fHits   = new TClonesArray("AliZDCHit",1000); 
271   
272   AliDetector::MakeBranch(opt);
273 }
274
275 //_____________________________________________________________________________
276 void AliZDC::Hits2SDigits()
277 {
278   // Create summable digits from hits
279   
280   AliDebug(1,"\n        Entering AliZDC::Hits2SDigits() ");
281   
282   fLoader->LoadHits("read");
283   fLoader->LoadSDigits("recreate");
284   AliRunLoader* runLoader = fLoader->GetRunLoader();
285   AliZDCSDigit sdigit;
286   AliZDCSDigit* psdigit = &sdigit;
287
288   // Event loop
289   for(Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
290     Float_t pmCZNC=0, pmCZPC=0, pmCZNA=0, pmCZPA=0, pmZEM1 = 0, pmZEM2 = 0;
291     Float_t pmQZNC[4], pmQZPC[4], pmQZNA[4], pmQZPA[4];
292     for(Int_t i = 0; i < 4; i++) pmQZNC[i] = pmQZPC[i] =  pmQZNA[i] = pmQZPA[i] = 0;
293
294     runLoader->GetEvent(iEvent);
295     TTree* treeH = fLoader->TreeH();
296     Int_t ntracks = (Int_t) treeH->GetEntries();
297     ResetHits();
298
299     // Tracks loop
300     Int_t sector[2];
301     for(Int_t itrack = 0; itrack < ntracks; itrack++) {
302       treeH->GetEntry(itrack);
303       for(AliZDCHit* zdcHit = (AliZDCHit*)FirstHit(-1); zdcHit;
304                       zdcHit = (AliZDCHit*)NextHit()) { 
305                       
306         sector[0] = zdcHit->GetVolume(0);
307         sector[1] = zdcHit->GetVolume(1);
308         if((sector[1] < 1) || (sector[1] > 5)) {
309           Error("Hits2SDigits", "sector[0] = %d, sector[1] = %d", 
310                 sector[0], sector[1]);
311           continue;
312         }
313         Float_t lightQ = zdcHit->GetLightPMQ();
314         Float_t lightC = zdcHit->GetLightPMC();
315      
316         if(sector[0] == 1) { //ZNC 
317           pmCZNC += lightC;
318           pmQZNC[sector[1]-1] += lightQ;
319         } 
320         else if(sector[0] == 2) { //ZPC 
321           pmCZPC += lightC;
322           pmQZPC[sector[1]-1] += lightQ;
323         } 
324         else if(sector[0] == 3) { //ZEM 
325           if(sector[1] == 1) pmZEM1 += lightC;
326           else pmZEM2 += lightQ;
327         }
328         if(sector[0] == 4) { //ZNA 
329           pmCZNA += lightC;
330           pmQZNA[sector[1]-1] += lightQ;
331         } 
332         else if(sector[0] == 5) { //ZPA 
333           pmCZPA += lightC;
334           pmQZPA[sector[1]-1] += lightQ;
335         } 
336       }//Hits loop
337     }
338
339     // create the output tree
340     fLoader->MakeTree("S");
341     TTree* treeS = fLoader->TreeS();
342     const Int_t kBufferSize = 4000;
343     treeS->Branch(GetName(), "AliZDCSDigit", &psdigit, kBufferSize);
344
345     // Create sdigits for ZNC
346     sector[0] = 1; // Detector = ZNC
347     sector[1] = 0; // Common PM ADC
348     new(psdigit) AliZDCSDigit(sector, pmCZNC);
349     if(pmCZNC > 0) treeS->Fill();
350     for(Int_t j = 0; j < 4; j++) {
351       sector[1] = j+1; // Towers PM ADCs
352       new(psdigit) AliZDCSDigit(sector, pmQZNC[j]);
353       if(pmQZNC[j] > 0) treeS->Fill();
354     }
355   
356     // Create sdigits for ZPC
357     sector[0] = 2; // Detector = ZPC
358     sector[1] = 0; // Common PM ADC
359     new(psdigit) AliZDCSDigit(sector, pmCZPC);
360     if(pmCZPC > 0) treeS->Fill();
361     for(Int_t j = 0; j < 4; j++) {
362       sector[1] = j+1; // Towers PM ADCs
363       new(psdigit) AliZDCSDigit(sector, pmQZPC[j]);
364       if(pmQZPC[j] > 0) treeS->Fill();
365     }
366
367     // Create sdigits for ZEM
368     sector[0] = 3; 
369     sector[1] = 1; // Detector = ZEM1
370     new(psdigit) AliZDCSDigit(sector, pmZEM1);
371     if(pmZEM1 > 0) treeS->Fill();
372     sector[1] = 2; // Detector = ZEM2
373     new(psdigit) AliZDCSDigit(sector, pmZEM2);
374     if(pmZEM2 > 0) treeS->Fill();
375
376     // Create sdigits for ZNA
377     sector[0] = 4; // Detector = ZNA
378     sector[1] = 0; // Common PM ADC
379     new(psdigit) AliZDCSDigit(sector, pmCZNA);
380     if(pmCZNA > 0) treeS->Fill();
381     for(Int_t j = 0; j < 4; j++) {
382       sector[1] = j+1; // Towers PM ADCs
383       new(psdigit) AliZDCSDigit(sector, pmQZNA[j]);
384       if(pmQZNA[j] > 0) treeS->Fill();
385     }
386   
387     // Create sdigits for ZPA
388     sector[0] = 5; // Detector = ZPA
389     sector[1] = 0; // Common PM ADC
390     new(psdigit) AliZDCSDigit(sector, pmCZPA);
391     if(pmCZPA > 0) treeS->Fill();
392     for(Int_t j = 0; j < 4; j++) {
393       sector[1] = j+1; // Towers PM ADCs
394       new(psdigit) AliZDCSDigit(sector, pmQZPA[j]);
395       if(pmQZPA[j] > 0) treeS->Fill();
396     }
397
398     // write the output tree
399     fLoader->WriteSDigits("OVERWRITE");
400   }
401
402   fLoader->UnloadHits();
403   fLoader->UnloadSDigits();
404 }
405
406 //_____________________________________________________________________________
407 AliDigitizer* AliZDC::CreateDigitizer(AliRunDigitizer* manager) const
408 {
409   // Create the digitizer for ZDC
410
411   return new AliZDCDigitizer(manager);
412 }
413
414 //_____________________________________________________________________________
415 void AliZDC::Digits2Raw()
416 {
417   // Convert ZDC digits to raw data
418
419   // 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
420   //         + 24 int values for the corresponding out of time channels
421   // For the CAEN module V965 we have an Header, the Data Words and an End Of Block
422   //    12 channels x 2 gain chains read from 1st ADC module
423   //    12 channels x 2 gain chains read from 2nd ADC module
424   //    12 channels x 2 gain chains read from 3rd ADC module (o.o.t.)
425   //    12 channels x 2 gain chains read from 4rth ADC module (o.o.t.)
426   //
427   const int knADCData1=24, knADCData2=24; // In principle the 2 numbers can be different!
428   UInt_t lADCHeader1; 
429   UInt_t lADCHeader2; 
430   UInt_t lADCData1[knADCData1];
431   UInt_t lADCData2[knADCData2];
432   UInt_t lADCData3[knADCData1];
433   UInt_t lADCData4[knADCData2];
434   //
435   UInt_t lADCEndBlock;
436
437   // load the digits
438   fLoader->LoadDigits("read");
439   AliZDCDigit digit;
440   AliZDCDigit* pdigit = &digit;
441   TTree* treeD = fLoader->TreeD();
442   if(!treeD) return;
443   treeD->SetBranchAddress("ZDC", &pdigit);
444   //printf("\t AliZDC::Digits2Raw -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
445
446   // Fill data array
447   // ADC header
448   UInt_t lADCHeaderGEO = 0;
449   UInt_t lADCHeaderCRATE = 0;
450   UInt_t lADCHeaderCNT1 = knADCData1;
451   UInt_t lADCHeaderCNT2 = knADCData2;
452     
453   lADCHeader1 = lADCHeaderGEO << 27 | 0x1 << 25 | lADCHeaderCRATE << 16 |
454                lADCHeaderCNT1 << 8 ;
455   lADCHeader2 = lADCHeaderGEO << 27 | 0x1 << 25 | lADCHeaderCRATE << 16 |
456                lADCHeaderCNT2 << 8 ;
457       
458   // ADC data word
459   UInt_t lADCDataGEO = lADCHeaderGEO;
460   //
461   UInt_t lADCDataValue1[knADCData1];
462   UInt_t lADCDataValue2[knADCData2];
463   UInt_t lADCDataValue3[knADCData1];
464   UInt_t lADCDataValue4[knADCData2];
465   //
466   UInt_t lADCDataOvFlw1[knADCData1];
467   UInt_t lADCDataOvFlw2[knADCData2];
468   UInt_t lADCDataOvFlw3[knADCData1];
469   UInt_t lADCDataOvFlw4[knADCData2];
470   //
471   for(Int_t i=0; i<knADCData1 ; i++){
472     lADCDataValue1[i] = 0;
473     lADCDataOvFlw1[i] = 0;
474     lADCDataValue3[i] = 0;
475     lADCDataOvFlw3[i] = 0;
476   }
477   for(Int_t i=0; i<knADCData2 ; i++){
478     lADCDataValue2[i] = 0;
479     lADCDataOvFlw2[i] = 0;
480     lADCDataValue4[i] = 0;
481     lADCDataOvFlw4[i] = 0;
482   }
483   //
484   UInt_t lADCDataChannel = 0;
485   
486   // loop over digits
487   for(Int_t iDigit=0; iDigit<treeD->GetEntries(); iDigit++){
488     treeD->GetEntry(iDigit);
489     if(!pdigit) continue;
490     //digit.Print("");
491     
492     // *** ADC data
493     Int_t index=0;
494     if(digit.GetSector(1)!=5){ // ZDC signal channels
495       // *** ADC1 (ZN1, ZP1, ZEM1,2) or ADC3 (ZN1, ZP1, ZEM1,2 o.o.t.)
496       if(digit.GetSector(0)==1 || digit.GetSector(0)==2 || digit.GetSector(0)==3){
497         if(digit.GetSector(0)==1 || digit.GetSector(0)==2){
498           index = (digit.GetSector(0)-1) + 4*digit.GetSector(1); // ZN1 or ZP1
499           lADCDataChannel = 8*(digit.GetSector(0)-1) + digit.GetSector(1);
500         }
501         else if(digit.GetSector(0)==3){ // ZEM 1,2
502           index = 20 + (digit.GetSector(1)-1);
503           lADCDataChannel = 5 + 8*(digit.GetSector(1)-1);
504         }
505         //
506         /*printf("\t AliZDC::Digits2Raw -> idig%d det %d quad %d index %d, ADCch %d ADCVal[%d, %d]\n",
507                 iDigit,digit.GetSector(0),digit.GetSector(1),index,lADCDataChannel,
508                 digit.GetADCValue(0),digit.GetADCValue(1));// Ch. debug
509         */
510         //
511         if(iDigit<knADCData1){ // *** In-time signals
512           lADCDataValue1[index] = digit.GetADCValue(0);   // High gain ADC ch.    
513           if(lADCDataValue1[index] > 2047) lADCDataOvFlw1[index] = 1; 
514           lADCDataValue1[index+2] = digit.GetADCValue(1); // Low gain ADC ch.
515           if(lADCDataValue1[index+2] > 2047) lADCDataOvFlw1[index+2] = 1; 
516         
517           lADCData1[index] = lADCDataGEO << 27 | 0x1 << 24 | lADCDataChannel << 17 | 
518                           lADCDataOvFlw1[index] << 12 | (lADCDataValue1[index] & 0xfff); 
519           lADCData1[index+2] = lADCDataGEO << 27 | 0x1 << 24  | lADCDataChannel << 17 | 0x1 << 16 |
520                           lADCDataOvFlw1[index+2] << 12 | (lADCDataValue1[index+2] & 0xfff);  
521         }
522         else{ // *** Out-of-time signals
523           lADCDataValue3[index] = digit.GetADCValue(0);   // High gain ADC ch.    
524           if(lADCDataValue3[index] > 2047) lADCDataOvFlw3[index] = 1; 
525           lADCDataValue3[index+2] = digit.GetADCValue(1); // Low gain ADC ch.
526           if(lADCDataValue3[index+2] > 2047) lADCDataOvFlw3[index+2] = 1; 
527       
528           lADCData3[index] = lADCDataGEO << 27 | lADCDataChannel << 17 | 
529                           lADCDataOvFlw3[index] << 12 | (lADCDataValue3[index] & 0xfff); 
530           lADCData3[index+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
531                           lADCDataOvFlw3[index+2] << 12 | (lADCDataValue3[index+2] & 0xfff);  
532         }                  
533       }
534       // *** ADC2 (ZN2, ZP2) or ADC4 (ZN2, ZP2 o.o.t.)
535       else if(digit.GetSector(0)==4 || digit.GetSector(0)==5){
536         index = (digit.GetSector(0)-4) + 4*digit.GetSector(1); // ZN2 or ZP2
537         lADCDataChannel = 8*(digit.GetSector(0)-4) + digit.GetSector(1);
538         //
539         /*printf("\t AliZDC::Digits2Raw -> idig%d det %d quad %d index %d, ADCch %d ADCVal[%d, %d]\n",
540                 iDigit,digit.GetSector(0),digit.GetSector(1),index,lADCDataChannel,
541                 digit.GetADCValue(0),digit.GetADCValue(1));// Ch. debug
542         */
543         //
544         if(iDigit<knADCData2){ // *** In-time signals
545           lADCDataValue2[index] = digit.GetADCValue(0);
546           if(lADCDataValue2[index] > 2047) lADCDataOvFlw2[index] = 1; 
547           lADCDataValue2[index+2] = digit.GetADCValue(1);
548           if(lADCDataValue2[index+2] > 2047) lADCDataOvFlw2[index+2] = 1; 
549           //
550           lADCData2[index] =   lADCDataGEO << 27 | lADCDataChannel << 17 | 
551                           lADCDataOvFlw2[index] << 12 | (lADCDataValue2[index] & 0xfff); 
552           lADCData2[index+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
553                           lADCDataOvFlw2[index+2] << 12 | (lADCDataValue2[index+2] & 0xfff);   
554         }                 
555         else{ // *** Out-of-time signals
556           lADCDataValue4[index] = digit.GetADCValue(0);
557           if(lADCDataValue4[index] > 2047) lADCDataOvFlw4[index] = 1; 
558           lADCDataValue4[index+2] = digit.GetADCValue(1);
559           if(lADCDataValue4[index+2] > 2047) lADCDataOvFlw4[index+2] = 1; 
560           //
561           lADCData4[index] =   lADCDataGEO << 27 | lADCDataChannel << 17 | 
562                         lADCDataOvFlw4[index] << 12 | (lADCDataValue4[index] & 0xfff); 
563           lADCData4[index+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
564                         lADCDataOvFlw4[index+2] << 12 | (lADCDataValue4[index+2] & 0xfff);   
565         }                 
566       }
567     }
568     // *** ADC2 (Reference PTMs) or ADC4 (Reference PTMs o.o.t.)
569     else if(digit.GetSector(1)==5){
570       index = 20 + (digit.GetSector(0)-1)/3; 
571       lADCDataChannel = 5 + 8*(digit.GetSector(0)-1)/3;
572       //
573       /*printf("\t AliZDC::Digits2Raw -> idig%d det %d quad %d index %d, ADCch %d ADCVal[%d, %d]\n",
574                 iDigit,digit.GetSector(0),digit.GetSector(1),index,lADCDataChannel,
575                 digit.GetADCValue(0),digit.GetADCValue(1));// Ch. debug
576       */
577       //
578       if(iDigit<knADCData2){ // *** In-time signals
579         lADCDataValue2[index] = digit.GetADCValue(0);
580         if(lADCDataValue2[index] > 2047) lADCDataOvFlw2[index] = 1; 
581         lADCDataValue2[index+2] = digit.GetADCValue(1);
582         if(lADCDataValue2[index+2] > 2047) lADCDataOvFlw2[index+2] = 1; 
583         //
584         lADCData2[index] =   lADCDataGEO << 27 | lADCDataChannel << 17 | 
585                         lADCDataOvFlw2[index] << 12 | (lADCDataValue2[index] & 0xfff); 
586         lADCData2[index+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
587                         lADCDataOvFlw2[index+2] << 12 | (lADCDataValue2[index+2] & 0xfff);   
588       }                 
589       else{ // *** Out-of-time signals
590         lADCDataValue4[index] = digit.GetADCValue(0);
591         if(lADCDataValue4[index] > 2047) lADCDataOvFlw4[index] = 1; 
592         lADCDataValue4[index+2] = digit.GetADCValue(1);
593         if(lADCDataValue4[index+2] > 2047) lADCDataOvFlw4[index+2] = 1; 
594         //
595         lADCData4[index] =   lADCDataGEO << 27 | lADCDataChannel << 17 | 
596                       lADCDataOvFlw4[index] << 12 | (lADCDataValue4[index] & 0xfff); 
597         lADCData4[index+2] = lADCDataGEO << 27 | lADCDataChannel << 17 | 0x1 << 16 |
598                       lADCDataOvFlw4[index+2] << 12 | (lADCDataValue4[index+2] & 0xfff);   
599       }                 
600            
601     }
602     if((index<0) || (index>23)) {
603       Error("Digits2Raw", "sector[0] = %d, sector[1] = %d", 
604             digit.GetSector(0), digit.GetSector(1));
605       continue;
606     }
607     
608     
609   }
610   //
611   /*
612   for(Int_t i=0;i<knADCData1;i++) printf("\t ADCData1[%d] = %x\n",i,lADCData1[i]);
613   for(Int_t i=0;i<knADCData2;i++) printf("\t ADCData2[%d] = %x\n",i,lADCData2[i]);
614   for(Int_t i=0;i<knADCData1;i++) printf("\t ADCData3[%d] = %x\n",i,lADCData3[i]);
615   for(Int_t i=0;i<knADCData2;i++) printf("\t ADCData4[%d] = %x\n",i,lADCData4[i]);
616   */
617  
618   // End of Block
619   UInt_t lADCEndBlockGEO = lADCHeaderGEO;
620   UInt_t lADCEndBlockEvCount = gAlice->GetEventNrInRun();
621   //  
622   lADCEndBlock = lADCEndBlockGEO << 27 | 0x1 << 26 | lADCEndBlockEvCount;
623   //printf("\t AliZDC::Digits2Raw -> ADCEndBlock = %d\n",lADCEndBlock);
624
625
626   // open the output file
627   char fileName[30];
628   strcpy(fileName,AliDAQ::DdlFileName("ZDC",0));
629
630   AliFstream* file = new AliFstream(fileName);
631
632   // write the DDL data header
633   AliRawDataHeaderSim header;
634   header.fSize = sizeof(header) + 
635                  sizeof(lADCHeader1) + sizeof(lADCData1) + sizeof(lADCEndBlock) +
636                  sizeof(lADCHeader2) + sizeof(lADCData2) + sizeof(lADCEndBlock) +
637                  sizeof(lADCHeader1) + sizeof(lADCData3) + sizeof(lADCEndBlock) +
638                  sizeof(lADCHeader2) + sizeof(lADCData4) + sizeof(lADCEndBlock);
639   //
640   /*printf("sizeof header = %d, ADCHeader1 = %d, ADCData1 = %d, ADCEndBlock = %d\n",
641           sizeof(header),sizeof(lADCHeader1),sizeof(lADCData1),sizeof(lADCEndBlock));
642   printf("sizeof header = %d, ADCHeader2 = %d, ADCData2 = %d, ADCEndBlock = %d\n",
643           sizeof(header),sizeof(lADCHeader2),sizeof(lADCData2),sizeof(lADCEndBlock));
644   */
645   //     
646   header.SetAttribute(0);  // valid data
647   file->WriteBuffer((char*)(&header), sizeof(header));
648
649   // write the raw data and close the file
650   file->WriteBuffer((char*) &lADCHeader1, sizeof (lADCHeader1));
651   file->WriteBuffer((char*)(lADCData1), sizeof(lADCData1));
652   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
653   file->WriteBuffer((char*) &lADCHeader2, sizeof (lADCHeader2));
654   file->WriteBuffer((char*)(lADCData2), sizeof(lADCData2));
655   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
656   file->WriteBuffer((char*) &lADCHeader1, sizeof (lADCHeader1));
657   file->WriteBuffer((char*)(lADCData3), sizeof(lADCData3));
658   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
659   file->WriteBuffer((char*) &lADCHeader2, sizeof (lADCHeader2));
660   file->WriteBuffer((char*)(lADCData4), sizeof(lADCData4));
661   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
662   delete file;
663
664   // unload the digits
665   fLoader->UnloadDigits();
666 }
667
668 //_____________________________________________________________________________
669 Bool_t AliZDC::Raw2SDigits(AliRawReader* rawReader)
670 {
671   // Convert ZDC raw data to Sdigits
672   
673   AliLoader* loader = (gAlice->GetRunLoader())->GetLoader("ZDCLoader");
674   if(!loader) {
675     AliError("no ZDC loader found");
676     return kFALSE;
677   }
678
679 //  // Event loop
680   Int_t iEvent = 0;
681   while(rawReader->NextEvent()){
682     (gAlice->GetRunLoader())->GetEvent(iEvent++);
683     // Create the output digit tree
684     TTree* treeS = loader->TreeS();
685     if(!treeS){
686       loader->MakeTree("S");
687       treeS = loader->TreeS();
688     }
689     //
690     AliZDCSDigit sdigit;
691     AliZDCSDigit* psdigit = &sdigit;
692     const Int_t kBufferSize = 4000;
693     treeS->Branch("ZDC", "AliZDCSDigit",  &psdigit, kBufferSize);
694     //
695     AliZDCRawStream rawStream(rawReader);
696     Int_t sector[2], resADC, rawADC, corrADC, nPheVal;
697     Int_t jcount = 0;
698     while(rawStream.Next()){
699       if(rawStream.IsADCDataWord()){
700         //For the moment only in-time SDigits are foreseen (1st 48 raw values)
701         if(jcount < 48){ 
702           for(Int_t j=0; j<2; j++) sector[j] = rawStream.GetSector(j);
703           rawADC = rawStream.GetADCValue();
704           resADC = rawStream.GetADCGain();
705           //printf("\t RAw2SDigits raw%d ->  RawADC[%d, %d, %d] read\n",
706           //    jcount, sector[0], sector[1], rawADC);
707           //
708           corrADC = rawADC - Pedestal(sector[0], sector[1], resADC);
709           if(corrADC<0) corrADC=0;
710           nPheVal = ADCch2Phe(sector[0], sector[1], corrADC, resADC);
711           //
712           //printf("\t \t ->  SDigit[%d, %d, %d] created\n",
713           //    sector[0], sector[1], nPheVal);
714           //
715           new(psdigit) AliZDCSDigit(sector, (Float_t) nPheVal);
716           treeS->Fill();
717           jcount++;
718         }
719       }//IsADCDataWord
720     }//rawStream.Next
721     // write the output tree
722     fLoader->WriteSDigits("OVERWRITE");
723     fLoader->UnloadSDigits();
724   }//Event loop 
725    
726   return kTRUE;
727 }
728
729 //_____________________________________________________________________________
730 Int_t AliZDC::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
731 {
732   // Returns a pedestal for detector det, PM quad, channel with res.
733   //
734   // Getting calibration object for ZDC set
735   AliCDBManager *man = AliCDBManager::Instance();
736   AliCDBEntry  *entry = man->Get("ZDC/Calib/Pedestals");
737   AliZDCPedestals *calibPed = (AliZDCPedestals*) entry->GetObject();
738   //
739   if(!calibPed){
740     printf("\t No calibration object found for ZDC!");
741     return -1;
742   }
743   //
744   Int_t index=0, kNch=24;
745   if(Quad!=5){
746     if(Det==1)        index = Quad+kNch*Res;     // ZN1
747     else if(Det==2)   index = Quad+5+kNch*Res;           // ZP1
748     else if(Det==3)   index = Quad+9+kNch*Res; // ZEM
749     else if(Det==4)   index = Quad+12+kNch*Res; // ZN2
750     else if(Det==5)   index = Quad+17+kNch*Res; // ZP2
751   }
752   else index = (Det-1)/3+22+kNch*Res; // Reference PMs
753   //
754   //
755   Float_t meanPed = calibPed->GetMeanPed(index);
756   Float_t pedWidth = calibPed->GetMeanPedWidth(index);
757   Float_t pedValue = gRandom->Gaus(meanPed,pedWidth);
758   //
759   //printf("\t AliZDC::Pedestal - det(%d, %d) - Ped[%d] = %d\n",Det, Quad, index,(Int_t) pedValue); // Chiara debugging!
760   
761   
762
763   return (Int_t) pedValue;
764 }
765
766
767 //_____________________________________________________________________________
768 Int_t AliZDC::ADCch2Phe(Int_t Det, Int_t Quad, Int_t ADCVal, Int_t Res) const
769 {
770   // Evaluation of the no. of phe produced
771   Float_t pmGain[6][5];
772   Float_t resADC[2];
773   for(Int_t j = 0; j < 5; j++){
774     pmGain[0][j] = 50000.;
775     pmGain[1][j] = 100000.;
776     pmGain[2][j] = 100000.;
777     pmGain[3][j] = 50000.;
778     pmGain[4][j] = 100000.;
779     pmGain[5][j] = 100000.;
780   }
781   // ADC Caen V965
782   resADC[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
783   resADC[1] = 0.0000064; // ADC Resolution low gain:  25  fC/adcCh
784   //
785   Int_t nPhe = (Int_t) (ADCVal * pmGain[Det-1][Quad] * resADC[Res]);
786   //
787   //printf("\t AliZDC::ADCch2Phe -> det(%d, %d) - ADC %d  phe %d\n",Det,Quad,ADCVal,nPhe);
788
789   return nPhe;
790 }
791
792 //______________________________________________________________________
793 void AliZDC::SetTreeAddress(){
794
795   // Set branch address for the Trees.
796   if(fLoader->TreeH() && (fHits == 0x0))
797     fHits   = new TClonesArray("AliZDCHit",1000);
798       
799   AliDetector::SetTreeAddress();
800 }