]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDC.cxx
Opening more the TDC check window
[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 <TClonesArray.h>
29 #include <TTree.h>
30 #include <TFile.h>
31 #include <TSystem.h>
32 #include <TRandom.h>
33 #include <TParticle.h>
34
35 // --- AliRoot header files
36 #include "AliDetector.h"
37 #include "AliRawDataHeaderSim.h"
38 #include "AliRawReader.h"
39 #include "AliLoader.h"
40 #include "AliRun.h"
41 #include "AliMC.h"
42 #include "AliLog.h"
43 #include "AliDAQ.h"
44 #include "AliZDC.h"
45 #include "AliZDCHit.h"
46 #include "AliZDCSDigit.h"
47 #include "AliZDCDigit.h"
48 #include "AliZDCDigitizer.h"
49 #include "AliZDCRawStream.h"
50 #include "AliZDCPedestals.h"
51 #include "AliZDCEnCalib.h"
52 #include "AliZDCTowerCalib.h"
53 #include "AliFstream.h"
54
55  
56 ClassImp(AliZDC)
57
58 //_____________________________________________________________________________
59 AliZDC::AliZDC() :
60   AliDetector(),
61   fNoShower(0),
62   fPedCalib(0),
63   fEnCalibData(0),
64   fTowCalibData(0),
65   fZDCCalibFName(""),
66   fSpectatorTracked(1)
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   fEnCalibData(0),
86   fTowCalibData(0),
87   fZDCCalibFName(""),
88   fSpectatorTracked(1)
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   SetName("ZDC"); SetTitle("ZDC");
103
104 }
105
106 //____________________________________________________________________________ 
107 AliZDC::~AliZDC()
108 {
109   //
110   // ZDC destructor
111   //
112
113   fIshunt = 0;
114   if(fPedCalib) delete fPedCalib;
115   if(fEnCalibData) delete fEnCalibData;
116   if(fEnCalibData) delete fEnCalibData;
117
118 }
119
120 //_____________________________________________________________________________
121 AliZDC::AliZDC(const AliZDC& ZDC) :
122 AliDetector("ZDC","ZDC"),
123 fNoShower(ZDC.fNoShower),
124 fPedCalib(ZDC.fPedCalib),
125 fEnCalibData(ZDC.fEnCalibData),
126 fTowCalibData(ZDC.fTowCalibData),
127 fZDCCalibFName(ZDC.fZDCCalibFName),
128 fSpectatorTracked(ZDC.fSpectatorTracked)
129 {
130   // copy constructor
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     fEnCalibData = ZDC.fEnCalibData;
141     fTowCalibData = ZDC.fTowCalibData;
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   
152   static Float_t trackTime=0., primKinEn=0., xImpact=0., yImpact=0., sFlag=0.;
153   static Int_t   pcPDGcode, motPDGcode;
154
155   AliZDCHit *newquad, *curprimquad;
156   newquad = new AliZDCHit(fIshunt, track, vol, hits);
157   TClonesArray &lhits = *fHits;
158   
159   if(fNhits==0){
160       // First hit -> setting flag for primary or secondary particle
161       TParticle * p = gAlice->GetMCApp()->Particle(track);
162       Int_t imo = p->GetFirstMother();
163       //
164       if(track != imo){
165         newquad->SetSFlag(1);  // SECONDARY particle entering the ZDC
166       }
167       else if(track == imo){
168         newquad->SetSFlag(0);  // PRIMARY particle entering the ZDC
169       }
170       //  
171       sFlag      = newquad->GetSFlag();
172       primKinEn  = newquad->GetPrimKinEn();
173       xImpact    = newquad->GetXImpact();
174       yImpact    = newquad->GetYImpact();
175       pcPDGcode  = newquad->GetPDGCode();
176       motPDGcode = newquad->GetMotherPDGCode();
177       trackTime  = newquad->GetTrackTOF();
178    }
179    else{       
180       newquad->SetPrimKinEn(primKinEn);
181       newquad->SetXImpact(xImpact);
182       newquad->SetYImpact(yImpact);
183       newquad->SetSFlag(sFlag);
184       newquad->SetPDGCode(pcPDGcode);
185       newquad->SetMotherPDGCode(motPDGcode);
186       newquad->SetTrackTOF(trackTime);
187    }
188  
189   Int_t j;
190   for(j=0; j<fNhits; j++){
191     // If hits are equal (same track, same volume), sum them.
192      curprimquad = (AliZDCHit*) lhits[j];
193      if(*curprimquad == *newquad){
194         *curprimquad = *curprimquad+*newquad;
195         // Ch. debug
196         //printf("\n\t Summing hits **************** \n", fNhits);
197         //curprimquad->Print("");
198         //
199         delete newquad;
200         return;
201      } 
202   }
203
204     //Otherwise create a new hit
205     new(lhits[fNhits]) AliZDCHit(*newquad);
206     fNhits++;
207     // Ch. debug
208     //printf("\n\t New ZDC hit added! fNhits = %d\n", fNhits);
209     //newquad->Print("");
210     
211     delete newquad;
212 }
213
214 //____________________________________________________________________________
215 Float_t AliZDC::ZMin(void) const
216 {
217   // Minimum dimension of the ZDC module in z
218   return -11600.;
219 }
220
221 //____________________________________________________________________________
222 Float_t AliZDC::ZMax(void) const
223 {
224   // Maximum dimension of the ZDC module in z
225   return 11750.;
226 }
227   
228
229 //_____________________________________________________________________________
230 void AliZDC::MakeBranch(Option_t *opt)
231 {
232   //
233   // Create Tree branches for the ZDC
234   //
235
236   char branchname[10];
237   snprintf(branchname, 10, "%s", GetName());
238
239   const char *cH = strstr(opt,"H");
240   
241   if(cH && fLoader->TreeH()) {
242     if (fHits) {
243       fHits->Clear();
244       fNhits = 0;
245     }
246     else {
247       fHits   = new TClonesArray("AliZDCHit",1000); 
248       if (gAlice && gAlice->GetMCApp())
249         gAlice->GetMCApp()->AddHitList(fHits);
250     }
251   }
252   
253   AliDetector::MakeBranch(opt);
254 }
255
256 //_____________________________________________________________________________
257 void AliZDC::Hits2SDigits()
258 {
259   // Create summable digits from hits
260   
261   AliDebug(1,"\n        AliZDC::Hits2SDigits() ");
262   
263   fLoader->LoadHits("read");
264   fLoader->LoadSDigits("recreate");
265   AliRunLoader* runLoader = fLoader->GetRunLoader();
266   AliZDCSDigit sdigit;
267   AliZDCSDigit* psdigit = &sdigit;
268
269   // Event loop
270   for(Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
271     Float_t pmZNC[5], pmZPC[5], pmZNA[5], pmZPA[5], pmZEM1=0., pmZEM2=0.;
272     for(Int_t i=0; i<5; i++) pmZNC[i] = pmZPC[i] =  pmZNA[i] = pmZPA[i] = 0;
273
274     runLoader->GetEvent(iEvent);
275     TTree* treeH = fLoader->TreeH();
276     Int_t ntracks = (Int_t) treeH->GetEntries();
277     ResetHits();
278
279     // Tracks loop
280     Int_t sector[2]; Float_t trackTime = 0.;
281     for(Int_t itrack = 0; itrack < ntracks; itrack++) {
282       treeH->GetEntry(itrack);
283       for(AliZDCHit* zdcHit = (AliZDCHit*)FirstHit(-1); zdcHit;
284           zdcHit = (AliZDCHit*)NextHit()) { 
285                       
286         sector[0] = zdcHit->GetVolume(0);
287         sector[1] = zdcHit->GetVolume(1);
288         if((sector[1] < 1) || (sector[1]>5)) {
289           Error("Hits2SDigits", "sector[0] = %d, sector[1] = %d", sector[0], sector[1]);
290           continue;
291         }
292         Float_t lightQ = zdcHit->GetLightPMQ();
293         Float_t lightC = zdcHit->GetLightPMC();
294         trackTime = zdcHit->GetTrackTOF();
295         // Signals from ZEM are delayed to arrive in time with ZDC signals
296         if(sector[0] == 3)  trackTime += 320;
297         // Ch. debug
298         //printf("\t det %d vol %d trackTOF %f lightQ %1.0f lightC %1.0f\n",
299         //      sector[0], sector[1], trackTime, lightQ, lightC);
300      
301         if(sector[0] == 1) { //ZNC 
302           pmZNC[0] += lightC;
303           pmZNC[sector[1]] += lightQ;
304         } 
305         else if(sector[0] == 2) { //ZPC 
306           pmZPC[0] += lightC;
307           pmZPC[sector[1]] += lightQ;
308         } 
309         else if(sector[0] == 3) { //ZEM 
310           if(sector[1] == 1) pmZEM1 += lightC;
311           else pmZEM2 += lightQ;
312         }
313         if(sector[0] == 4) { //ZNA 
314           pmZNA[0] += lightC;
315           pmZNA[sector[1]] += lightQ;
316         } 
317         else if(sector[0] == 5) { //ZPA 
318           pmZPA[0] += lightC;
319           pmZPA[sector[1]] += lightQ;
320         } 
321       }//Hits loop
322     }//Tracks loop
323
324     // create the output tree
325     fLoader->MakeTree("S");
326     TTree* treeS = fLoader->TreeS();
327     const Int_t kBufferSize = 4000;
328     treeS->Branch(GetName(), "AliZDCSDigit", &psdigit, kBufferSize);
329
330     // Create sdigits for ZNC
331     sector[0] = 1; // Detector = ZNC
332     for(Int_t j = 0; j < 5; j++) {
333       sector[1] = j; 
334       if(pmZNC[j]>0){
335         new(psdigit) AliZDCSDigit(sector, pmZNC[j], trackTime);
336         treeS->Fill();
337         // Ch. debug
338         //printf("\t SDigit created: det %d quad %d pmZNC[%d] %1.0f trackTOF %f\n",
339         //      sector[0], sector[1], j, pmZNC[j], trackTime);
340       }
341     }
342   
343     // Create sdigits for ZPC
344     sector[0] = 2; // Detector = ZPC
345     for(Int_t j = 0; j < 5; j++) {
346       sector[1] = j; // Towers PM ADCs
347       if(pmZPC[j]>0){
348         new(psdigit) AliZDCSDigit(sector, pmZPC[j], trackTime);
349         treeS->Fill();
350         // Ch. debug
351         //printf("\t SDigit created: det %d quad %d pmZPC[%d] %1.0f trackTOF %f\n",
352         //      sector[0], sector[1], j, pmZPC[j], trackTime);
353       }
354     }
355
356     // Create sdigits for ZEM
357     sector[0] = 3; 
358     sector[1] = 1; // Detector = ZEM1
359     if(pmZEM1>0){ 
360       new(psdigit) AliZDCSDigit(sector, pmZEM1, trackTime);
361       treeS->Fill();
362       // Ch. debug
363       //printf("\t SDigit created: det %d quad %d pmZEM1 %1.0f trackTOF %f\n",
364       //        sector[0], sector[1], pmZEM1, trackTime);
365     }
366     sector[1] = 2; // Detector = ZEM2
367     if(pmZEM2>0){
368       new(psdigit) AliZDCSDigit(sector, pmZEM2, trackTime);
369       treeS->Fill();
370       // Ch. debug
371       //printf("\t SDigit created: det %d quad %d pmZEM2 %1.0f trackTOF %f\n",
372       //        sector[0], sector[1], pmZEM2, trackTime);
373     }
374
375     // Create sdigits for ZNA
376     sector[0] = 4; // Detector = ZNA
377     for(Int_t j = 0; j < 5; j++) {
378       sector[1] = j; // Towers PM ADCs
379       if(pmZNA[j]>0){
380         new(psdigit) AliZDCSDigit(sector, pmZNA[j], trackTime);
381         treeS->Fill();
382         // Ch. debug
383         //printf("\t SDigit created: det %d quad %d pmZNA[%d] %1.0f trackTOF %f\n",
384         //      sector[0], sector[1], j, pmZNA[j], trackTime);
385       }
386     }
387   
388     // Create sdigits for ZPA
389     sector[0] = 5; // Detector = ZPA
390     sector[1] = 0; // Common PM ADC
391     for(Int_t j = 0; j < 5; j++) {
392       sector[1] = j; // Towers PM ADCs
393       if(pmZPA[j]>0){
394         new(psdigit) AliZDCSDigit(sector, pmZPA[j], trackTime);
395         treeS->Fill();
396         // Ch. debug
397         //printf("\t SDigit created: det %d quad %d pmZPA[%d] %1.0f trackTOF %f\n",
398         //      sector[0], sector[1], j, pmZPA[j], trackTime);
399       }
400     }
401
402     // write the output tree
403     fLoader->WriteSDigits("OVERWRITE");
404   }
405
406   fLoader->UnloadHits();
407   fLoader->UnloadSDigits();
408 }
409
410 //_____________________________________________________________________________
411 AliDigitizer* AliZDC::CreateDigitizer(AliRunDigitizer* manager) const
412 {
413   // Create the digitizer for ZDC
414   AliZDCDigitizer *zdcDigitizer = new AliZDCDigitizer(manager);
415   if(fSpectatorTracked==0) zdcDigitizer->SetSpectators2Track();
416   //printf("\n**************************ZDC digitizer created with Spectators2Track = %d\n\n", fSpectatorTracked);
417   return zdcDigitizer;
418 }
419
420 //_____________________________________________________________________________
421 void AliZDC::Digits2Raw()
422 {
423   // Convert ZDC digits to raw data
424
425   // 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
426   //         + 24 int values for the corresponding out of time channels
427   // For the CAEN module V965 we have an Header, the Data Words and an End Of Block
428   //    12 channels x 2 gain chains read from 1st ADC module
429   //    12 channels x 2 gain chains read from 2nd ADC module
430   //    12 channels x 2 gain chains read from 3rd ADC module (o.o.t.)
431   //    12 channels x 2 gain chains read from 4rth ADC module (o.o.t.)
432   //
433   const int knADCData1=12, knADCData2=12; 
434   const int knADCData3=12, knADCData4=12; 
435   //
436   UInt_t lADCHeader1; 
437   UInt_t lADCHeader2; 
438   UInt_t lADCHeader3; 
439   UInt_t lADCHeader4; 
440   //
441   UInt_t lADCData1[2*knADCData1];
442   UInt_t lADCData2[2*knADCData2];
443   UInt_t lADCData3[2*knADCData3];
444   UInt_t lADCData4[2*knADCData4];
445   //
446   UInt_t lADCEndBlock;
447
448   // load the digits
449   fLoader->LoadDigits("read");
450   AliZDCDigit digit;
451   AliZDCDigit* pdigit = &digit;
452   TTree* treeD = fLoader->TreeD();
453   if(!treeD) return;
454   treeD->SetBranchAddress("ZDC", &pdigit);
455   //printf("\t AliZDC::Digits2Raw -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
456
457   // Reading channel map
458   //printf("\n\t Reading ADC mapping from OCDB\n");
459   AliZDCChMap * chMap = GetChMap();
460   const int nCh = knADCData1+knADCData2+knADCData3+knADCData4;
461   Int_t  mapADC[nCh][4]; 
462   for(Int_t i=0; i<nCh; i++){
463     mapADC[i][0] = chMap->GetADCModule(i);
464     mapADC[i][1] = chMap->GetADCChannel(i);
465     mapADC[i][2] = chMap->GetDetector(i);
466     mapADC[i][3] = chMap->GetSector(i);
467     // Ch. debug
468     //printf("  mapADC[%d] = (%d %d %d %d)\n", i,
469     //  mapADC[i][0],mapADC[i][1],mapADC[i][2],mapADC[i][3]);
470   }
471
472   // *** Fill data array
473   // ** ADC header
474   UInt_t lADCHeaderGEO1 = 0;
475   UInt_t lADCHeaderGEO2 = 1;
476   UInt_t lADCHeaderGEO3 = 2;
477   UInt_t lADCHeaderGEO4 = 3;
478   UInt_t lADCHeaderCRATE = 0;
479   UInt_t lADCHeaderCNT1 = knADCData1;
480   UInt_t lADCHeaderCNT2 = knADCData2;
481   UInt_t lADCHeaderCNT3 = knADCData3;
482   UInt_t lADCHeaderCNT4 = knADCData4;
483     
484   lADCHeader1 = lADCHeaderGEO1 << 27 | 0x1 << 25 | lADCHeaderCRATE << 16 |
485                lADCHeaderCNT1 << 8 ;
486   lADCHeader2 = lADCHeaderGEO2 << 27 | 0x1 << 25 | lADCHeaderCRATE << 16 |
487                lADCHeaderCNT2 << 8 ;
488   lADCHeader3 = lADCHeaderGEO3 << 27 | 0x1 << 25 | lADCHeaderCRATE << 16 |
489                lADCHeaderCNT3 << 8 ;
490   lADCHeader4 = lADCHeaderGEO4 << 27 | 0x1 << 25 | lADCHeaderCRATE << 16 |
491                lADCHeaderCNT4 << 8 ;
492       
493   // ** ADC data word
494   UInt_t lADCDataGEO = 0;
495   //
496   UInt_t lADCDataValue1[2*knADCData1];
497   UInt_t lADCDataValue2[2*knADCData2];
498   UInt_t lADCDataValue3[2*knADCData3];
499   UInt_t lADCDataValue4[2*knADCData4];
500   //
501   UInt_t lADCDataOvFlwHG = 0;
502   UInt_t lADCDataOvFlwLG = 0;
503   //
504   for(Int_t i=0; i<2*knADCData1 ; i++) lADCDataValue1[i] = 0;
505   for(Int_t i=0; i<2*knADCData2 ; i++) lADCDataValue2[i] = 0;
506   for(Int_t i=0; i<2*knADCData3 ; i++) lADCDataValue3[i] = 0;
507   for(Int_t i=0; i<2*knADCData4 ; i++) lADCDataValue4[i] = 0;
508   //
509   UInt_t lADCDataChannel = 0;
510   
511   Int_t indADC0=0, indADC1=0, indADC2=0, indADC3=0;
512   
513   // loop over digits
514   for(Int_t iDigit=0; iDigit<(Int_t) (treeD->GetEntries()); iDigit++){
515     treeD->GetEntry(iDigit);
516     if(!pdigit) continue;
517     //digit.Print("");
518    
519     // *** ADC data
520     // Scan of the map to assign the correct ADC module-channel
521     for(Int_t k=0; k<nCh; k++){
522       if(iDigit<knADCData1+knADCData2){ 
523        if(digit.GetSector(0)==mapADC[k][2] && digit.GetSector(1)==mapADC[k][3]){
524          lADCDataGEO = (UInt_t) mapADC[k][0];
525          lADCDataChannel = (UInt_t) mapADC[k][1];
526          break;
527        } 
528       }
529       else if(iDigit>knADCData1+knADCData2){
530        if(digit.GetSector(0)==mapADC[k][2] && digit.GetSector(1)==mapADC[k][3]){
531          lADCDataGEO = (UInt_t) mapADC[k][0];
532          lADCDataChannel = (UInt_t) mapADC[k][1];
533          if(k>knADCData1+knADCData2) break;
534        } 
535       }
536     }
537     // Ch. debug
538     //printf(" det %d sec %d -> lADCDataGEO %d  lADCDataChannel %d\n",
539     //  digit.GetSector(0),digit.GetSector(1),lADCDataGEO,lADCDataChannel);
540      
541     if(lADCDataGEO==0){ 
542       if(indADC0>=knADCData1){
543         AliError(" Problem with digit index 4 ADC0\n");
544         return;
545       }
546       Int_t indLG = indADC0+knADCData1;
547       // High gain ADC ch.       
548       if(digit.GetADCValue(0) > 2047) lADCDataOvFlwHG = 1; 
549       lADCDataValue1[indADC0] = digit.GetADCValue(0);    
550       lADCData1[indADC0] = lADCDataGEO << 27 |  lADCDataChannel << 17 | 
551                     lADCDataOvFlwHG << 12 | (lADCDataValue1[indADC0] & 0xfff); 
552       // Low gain ADC ch.
553       if(digit.GetADCValue(1) > 2047) lADCDataOvFlwLG = 1; 
554       lADCDataValue1[indLG-1] = digit.GetADCValue(1); 
555       lADCData1[indLG] = lADCDataGEO << 27 |  lADCDataChannel << 17 | 0x1 << 16 |
556                     lADCDataOvFlwLG << 12 | (lADCDataValue1[indLG] & 0xfff);  
557                     
558       indADC0++;
559       //
560       // Ch. debug
561       //printf(" lADCDataGEO %d lADCDataValue1[%d] = %d  lADCDataValue1[%d] = %d\n", 
562       //  lADCDataGEO,iDigit,lADCDataValue1[indLG],indLG,lADCDataValue1[indLG]);
563     }
564     else if(lADCDataGEO==1){ 
565       if(indADC1>=knADCData2){
566          AliError(" Problem with digit index 4 ADC1\n");
567          return;
568       }
569       Int_t indLG = indADC1+knADCData2;
570       // High gain ADC ch.       
571       if(digit.GetADCValue(0) > 2047) lADCDataOvFlwHG = 1; 
572       lADCDataValue2[indADC1] = digit.GetADCValue(0);    
573       lADCData2[indADC1] = lADCDataGEO << 27 | lADCDataChannel << 17 | 
574                     lADCDataOvFlwHG << 12 | (lADCDataValue2[indADC1] & 0xfff); 
575       // Low gain ADC ch.
576       if(digit.GetADCValue(1) > 2047) lADCDataOvFlwLG = 1; 
577       lADCDataValue2[indLG] = digit.GetADCValue(1); 
578       lADCData2[indLG] = lADCDataGEO << 27 |  lADCDataChannel << 17 | 0x1 << 16 |
579                     lADCDataOvFlwLG << 12 | (lADCDataValue2[indLG] & 0xfff);  
580                   
581       indADC1++;
582       //
583       //Ch. debug
584       //printf(" lADCDataGEO %d  lADCDataValue2[%d] = %d  lADCDataValue2[%d] = %d\n", 
585       //  lADCDataGEO,indHG,lADCDataValue2[indHG],indLG,lADCDataValue2[indLG]);
586     }
587     else if(lADCDataGEO==2){ 
588       if(indADC2>=knADCData3){
589         AliError(" Problem with digit index 4 ADC2\n");
590         return;
591       }
592       Int_t indLG = indADC2+knADCData3;
593       // High gain ADC ch.       
594       if(digit.GetADCValue(0) > 2047) lADCDataOvFlwHG = 1; 
595       lADCDataValue3[indADC1] = digit.GetADCValue(0);    
596       lADCData3[indADC1] = lADCDataGEO << 27 | lADCDataChannel << 17 | 
597                     lADCDataOvFlwHG << 12 | (lADCDataValue3[indADC2] & 0xfff); 
598       // Low gain ADC ch.
599       if(digit.GetADCValue(1) > 2047) lADCDataOvFlwLG = 1; 
600       lADCDataValue3[indLG] = digit.GetADCValue(1); 
601       lADCData3[indLG] = lADCDataGEO << 27 |  lADCDataChannel << 17 | 0x1 << 16 |
602                     lADCDataOvFlwLG << 12 | (lADCDataValue3[indLG] & 0xfff);  
603                   
604       indADC2++;
605       //
606       //Ch. debug
607       //printf(" lADCDataGEO %d   lADCDataValue3[%d] = %d  lADCDataValue3[%d] = %d\n", 
608       //  lADCDataGEO,indHG,lADCDataValue3[indHG],indLG,lADCDataValue3[indLG]);
609     }
610     else if(lADCDataGEO==3){ 
611       if(indADC3>=knADCData4){
612          AliError(" Problem with digit index 4 ADC2\n");
613          return;
614       }
615       Int_t indLG = indADC3+knADCData4;
616       // High gain ADC ch.       
617       if(digit.GetADCValue(0) > 2047) lADCDataOvFlwHG = 1; 
618       lADCDataValue4[indADC3] = digit.GetADCValue(0);    
619       lADCData4[indADC3] = lADCDataGEO << 27 | lADCDataChannel << 17 | 
620                     lADCDataOvFlwHG << 12 | (lADCDataValue4[indADC3] & 0xfff); 
621       // Low gain ADC ch.
622       if(digit.GetADCValue(1) > 2047) lADCDataOvFlwLG = 1; 
623       lADCDataValue4[indLG] = digit.GetADCValue(1); 
624       lADCData4[indLG] = lADCDataGEO << 27 |  lADCDataChannel << 17 | 0x1 << 16 |
625                     lADCDataOvFlwLG << 12 | (lADCDataValue4[indLG] & 0xfff);  
626                   
627       indADC3++;
628       //
629       // Ch. debug
630       //printf(" lADCDataGEO %d lADCDataValue4[%d] = %d  lADCDataValue4[%d] = %d\n", 
631       //  lADCDataGEO,indHG,lADCDataValue4[indHG],indLG,lADCDataValue4[indLG]);
632     }             
633
634   }
635   //
636   /*for(Int_t i=0;i<2*knADCData1;i++) printf("\t ADCData1[%d] = %x\n",i,lADCData1[i]);
637   for(Int_t i=0;i<2*knADCData2;i++) printf("\t ADCData2[%d] = %x\n",i,lADCData2[i]);
638   for(Int_t i=0;i<2*knADCData3;i++) printf("\t ADCData3[%d] = %x\n",i,lADCData3[i]);
639   for(Int_t i=0;i<2*knADCData4;i++) printf("\t ADCData4[%d] = %x\n",i,lADCData4[i]);*/
640    
641   // End of Block
642   UInt_t lADCEndBlockGEO = 0;
643   // Event counter in ADC EOB -> getting no. of events in run from AliRunLoader
644   // get run loader
645   AliRunLoader* runLoader = fLoader->GetRunLoader(); 
646   UInt_t lADCEndBlockEvCount = runLoader->GetEventNumber();
647   //  
648   lADCEndBlock = lADCEndBlockGEO << 27 | 0x1 << 26 | lADCEndBlockEvCount;
649   //printf("\t AliZDC::Digits2Raw -> ADCEndBlock = %d\n",lADCEndBlock);
650
651   // open the output file
652   TString fileName;
653   fileName.Form("%s",AliDAQ::DdlFileName("ZDC",0)); 
654
655   AliFstream* file = new AliFstream(fileName.Data());
656
657   // write the DDL data header
658   AliRawDataHeaderSim header;
659   header.fSize = sizeof(header) + 
660                  sizeof(lADCHeader1) + sizeof(lADCData1) + sizeof(lADCEndBlock) +
661                  sizeof(lADCHeader2) + sizeof(lADCData2) + sizeof(lADCEndBlock) +
662                  sizeof(lADCHeader3) + sizeof(lADCData3) + sizeof(lADCEndBlock) +
663                  sizeof(lADCHeader4) + sizeof(lADCData4) + sizeof(lADCEndBlock);
664   //
665   /*printf("sizeof header = %d, ADCHeader1 = %d, ADCData1 = %d, ADCEndBlock = %d\n",
666           sizeof(header),sizeof(lADCHeader1),sizeof(lADCData1),sizeof(lADCEndBlock));
667   printf("sizeof header = %d, ADCHeader2 = %d, ADCData2 = %d, ADCEndBlock = %d\n",
668           sizeof(header),sizeof(lADCHeader2),sizeof(lADCData2),sizeof(lADCEndBlock));
669   printf("sizeof header = %d, ADCHeader3 = %d, ADCData3 = %d, ADCEndBlock = %d\n",
670           sizeof(header),sizeof(lADCHeader1),sizeof(lADCData1),sizeof(lADCEndBlock));
671   printf("sizeof header = %d, ADCHeader4 = %d, ADCData4 = %d, ADCEndBlock = %d\n",
672           sizeof(header),sizeof(lADCHeader2),sizeof(lADCData2),sizeof(lADCEndBlock));*/
673   
674   header.SetAttribute(0);  // valid data
675   file->WriteBuffer((char*)(&header), sizeof(header));
676   // write the raw data and close the file
677   file->WriteBuffer((char*) &lADCHeader1,  sizeof (lADCHeader1));
678   file->WriteBuffer((char*) &lADCData1,   sizeof(lADCData1));
679   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
680   file->WriteBuffer((char*) &lADCHeader2,  sizeof (lADCHeader2));
681   file->WriteBuffer((char*) (lADCData2),   sizeof(lADCData2));
682   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
683   file->WriteBuffer((char*) &lADCHeader3,  sizeof (lADCHeader3));
684   file->WriteBuffer((char*) (lADCData3),   sizeof(lADCData3));
685   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
686   file->WriteBuffer((char*) &lADCHeader4,  sizeof (lADCHeader4));
687   file->WriteBuffer((char*) (lADCData4),   sizeof(lADCData4));
688   file->WriteBuffer((char*) &lADCEndBlock, sizeof(lADCEndBlock));
689   delete file;
690
691   // unload the digits
692   fLoader->UnloadDigits();
693 }
694
695 //_____________________________________________________________________________
696 Bool_t AliZDC::Raw2SDigits(AliRawReader* rawReader)
697 {
698   // Convert ZDC raw data to Sdigits
699   const int kNch = 48;
700   AliLoader* loader = (AliRunLoader::Instance())->GetLoader("ZDCLoader");
701   if(!loader) {
702     AliError("no ZDC loader found");
703     return kFALSE;
704   }
705
706   // Create the output digit tree
707   TTree* treeS = loader->TreeS();
708   if(!treeS){
709     loader->MakeTree("S");
710     treeS = loader->TreeS();
711   }
712   //
713   AliZDCSDigit sdigit;
714   AliZDCSDigit* psdigit = &sdigit;
715   const Int_t kBufferSize = 4000;
716   treeS->Branch("ZDC", "AliZDCSDigit",  &psdigit, kBufferSize);
717   //
718   AliZDCRawStream rawStream(rawReader);
719   Int_t sector[2], resADC, rawADC, corrADC, nPheVal;
720   Int_t jcount = 0;
721   while(rawStream.Next()){
722     if(rawStream.IsADCDataWord()){
723       //For the moment only in-time SDigits are foreseen (1st 48 raw values)
724       if(jcount < kNch){ 
725         for(Int_t j=0; j<2; j++) sector[j] = rawStream.GetSector(j);
726         rawADC = rawStream.GetADCValue();
727         resADC = rawStream.GetADCGain();
728         //printf("\t RAw2SDigits raw%d ->  RawADC[%d, %d, %d] read\n",
729         //    jcount, sector[0], sector[1], rawADC);
730         //
731         corrADC = rawADC - Pedestal(sector[0], sector[1], resADC);
732         if(corrADC<0) corrADC=0;
733         nPheVal = ADCch2Phe(sector[0], sector[1], corrADC, resADC);
734         //
735         //printf("\t \t ->  SDigit[%d, %d, %d] created\n",
736         //    sector[0], sector[1], nPheVal);
737         //
738         new(psdigit) AliZDCSDigit(sector, (Float_t) nPheVal, 0.);
739         treeS->Fill();
740         jcount++;
741       }
742     }//IsADCDataWord
743   }//rawStream.Next
744   // write the output tree
745   fLoader->WriteSDigits("OVERWRITE");
746   fLoader->UnloadSDigits();
747    
748   return kTRUE;
749 }
750
751 //_____________________________________________________________________________
752 Int_t AliZDC::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
753 {
754   // Returns a pedestal for detector det, PM quad, channel with res.
755   //
756   // Getting calibration object for ZDC set
757   AliCDBManager *man = AliCDBManager::Instance();
758   AliCDBEntry  *entry = man->Get("ZDC/Calib/Pedestals");
759   if(!entry) AliFatal("No calibration data loaded!");  
760   AliZDCPedestals *calibPed = (AliZDCPedestals*) entry->GetObject();
761   //
762   if(!calibPed){
763     printf("\t No calibration object found for ZDC!");
764     return -1;
765   }
766   //
767   Int_t index=0, kNch=24;
768   if(Quad!=5){
769     if(Det==1)        index = Quad+kNch*Res;     // ZN1
770     else if(Det==2)   index = Quad+5+kNch*Res;           // ZP1
771     else if(Det==3)   index = Quad+9+kNch*Res; // ZEM
772     else if(Det==4)   index = Quad+12+kNch*Res; // ZN2
773     else if(Det==5)   index = Quad+17+kNch*Res; // ZP2
774   }
775   else index = (Det-1)/3+22+kNch*Res; // Reference PMs
776   //
777   //
778   Float_t meanPed = calibPed->GetMeanPed(index);
779   Float_t pedWidth = calibPed->GetMeanPedWidth(index);
780   Float_t pedValue = gRandom->Gaus(meanPed,pedWidth);
781   //
782   //printf("\t AliZDC::Pedestal - det(%d, %d) - Ped[%d] = %d\n",Det, Quad, index,(Int_t) pedValue); // Chiara debugging!
783   
784   
785
786   return (Int_t) pedValue;
787 }
788
789
790 //_____________________________________________________________________________
791 Int_t AliZDC::ADCch2Phe(Int_t Det, Int_t Quad, Int_t ADCVal, Int_t Res) const
792 {
793   // Evaluation of the no. of phe produced
794   Float_t pmGain[6][5];
795   Float_t resADC[2];
796   for(Int_t j = 0; j < 5; j++){
797     pmGain[0][j] = 50000.;
798     pmGain[1][j] = 100000.;
799     pmGain[2][j] = 100000.;
800     pmGain[3][j] = 50000.;
801     pmGain[4][j] = 100000.;
802     pmGain[5][j] = 100000.;
803   }
804   // ADC Caen V965
805   resADC[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
806   resADC[1] = 0.0000064; // ADC Resolution low gain:  25  fC/adcCh
807   //
808   Int_t nPhe = (Int_t) (ADCVal / (pmGain[Det-1][Quad] * resADC[Res]));
809   //
810   //printf("\t AliZDC::ADCch2Phe -> det(%d, %d) - ADC %d  phe %d\n",Det,Quad,ADCVal,nPhe);
811
812   return nPhe;
813 }
814
815 //______________________________________________________________________
816 void AliZDC::SetTreeAddress(){
817
818   // Set branch address for the Trees.
819   if(fLoader->TreeH() && (fHits == 0x0))
820     fHits   = new TClonesArray("AliZDCHit",1000);
821       
822   AliDetector::SetTreeAddress();
823 }
824
825 //_____________________________________________________________________________
826 AliZDCChMap* AliZDC::GetChMap() const
827 {
828
829   // Getting calibration object for ZDC
830
831   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/ChMap");
832   if(!entry) AliFatal("No calibration data loaded!");  
833
834   AliZDCChMap *calibdata = dynamic_cast<AliZDCChMap*> (entry->GetObject());
835   if(!calibdata) AliFatal("Wrong calibration object in calibration  file!");
836
837   return calibdata;
838 }