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