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