]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDC.cxx
3be63336b227daf3f74b21947e225841e6ad8159
[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 <TGeometry.h>
30 #include <TNode.h>
31 #include <TTree.h>
32 #include <TSystem.h>
33
34 // --- AliRoot header files
35 #include "AliDetector.h"
36 #include "AliZDC.h"
37 #include "AliZDCHit.h"
38 #include "AliZDCSDigit.h"
39 #include "AliZDCDigit.h"
40 #include "AliZDCDigitizer.h"
41 #include "AliZDCRawStream.h"
42 #include "AliZDCCalibData.h"
43
44 #include "AliRawDataHeader.h"
45 #include "AliLoader.h"
46 #include "AliRun.h"
47 #include "AliMC.h"
48
49  
50 ClassImp(AliZDC)
51
52 AliZDC *gZDC;
53  
54 //_____________________________________________________________________________
55 AliZDC::AliZDC()
56 {
57   //
58   // Default constructor for the Zero Degree Calorimeter base class
59   //
60   
61   fIshunt     = 1;
62   fNoShower   = 0;
63
64   fHits       = 0;
65   fNhits      = 0;
66
67   fDigits     = 0;
68   fNdigits    = 0;
69
70 }
71  
72 //_____________________________________________________________________________
73 AliZDC::AliZDC(const char *name, const char *title)
74   : AliDetector(name,title)
75 {
76   //
77   // Standard constructor for the Zero Degree Calorimeter base class
78   //
79
80   fIshunt   = 1;
81   fNoShower = 0;
82
83   // Allocate the hits array  
84   fHits   = new TClonesArray("AliZDCHit",1000);
85   gAlice->GetMCApp()->AddHitList(fHits);
86   
87   char sensname[5],senstitle[25];
88   sprintf(sensname,"ZDC");
89   sprintf(senstitle,"ZDC dummy");
90   SetName(sensname); SetTitle(senstitle);
91
92   fDigits     = 0;
93   fNdigits    = 0;
94   
95   gZDC=this;
96
97 }
98 //____________________________________________________________________________ 
99 AliZDC::~AliZDC()
100 {
101   //
102   // ZDC destructor
103   //
104
105   fIshunt   = 0;
106   gZDC=0;
107
108 }
109 //_____________________________________________________________________________
110 void AliZDC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
111 {
112   //
113   //            Add a ZDC hit to the hit list.
114   // -> We make use of 2 array of hits:
115   // [1]  fHits (the usual one) that contains hits for each PRIMARY
116   // [2]  fStHits that contains hits for each EVENT and is used to
117   //      obtain digits at the end of each event
118   //
119   
120   static Float_t primKinEn, xImpact, yImpact, sFlag;
121
122   AliZDCHit *newquad, *curprimquad;
123   newquad = new AliZDCHit(fIshunt, track, vol, hits);
124   TClonesArray &lhits = *fHits;
125   
126   if(fNhits==0){
127       // First hit -> setting flag for primary or secondary particle
128       Int_t primary = gAlice->GetMCApp()->GetPrimary(track);     
129       if(track != primary){
130         newquad->SetSFlag(1);  // SECONDARY particle entering the ZDC
131       }
132       else if(track == primary){
133         newquad->SetSFlag(0);  // PRIMARY particle entering the ZDC
134       }  
135       sFlag     = newquad->GetSFlag();
136       primKinEn = newquad->GetPrimKinEn();
137       xImpact   = newquad->GetXImpact();
138       yImpact   = newquad->GetYImpact();
139    }
140    else{       
141       newquad->SetPrimKinEn(primKinEn);
142       newquad->SetXImpact(xImpact);
143       newquad->SetYImpact(yImpact);
144       newquad->SetSFlag(sFlag);
145    }
146  
147   Int_t j;
148   for(j=0; j<fNhits; j++){
149     // If hits are equal (same track, same volume), sum them.
150      curprimquad = (AliZDCHit*) lhits[j];
151      if(*curprimquad == *newquad){
152         *curprimquad = *curprimquad+*newquad;
153         delete newquad;
154         return;
155      } 
156   }
157
158     //Otherwise create a new hit
159     new(lhits[fNhits]) AliZDCHit(newquad);
160     fNhits++;
161     
162     delete newquad;
163 }
164
165 //_____________________________________________________________________________
166 void AliZDC::BuildGeometry()
167 {
168   //
169   // Build the ROOT TNode geometry for event display 
170   // in the Zero Degree Calorimeter
171   // This routine is dummy for the moment
172   //
173
174   TNode *node, *top;
175   TBRIK *brik;
176   const int kColorZDC  = kBlue;
177   
178   //
179   top=gAlice->GetGeometry()->GetNode("alice");
180   
181   // ZDC
182     brik = new TBRIK("S_ZDC","ZDC box","void",300,300,5);
183     top->cd();
184     node = new TNode("ZDC","ZDC","S_ZDC",0,0,600,"");
185     node->SetLineColor(kColorZDC);
186     fNodes->Add(node);
187 }
188
189 //_____________________________________________________________________________
190 Int_t AliZDC::DistancetoPrimitive(Int_t , Int_t )
191 {
192   //
193   // Distance from the mouse to the Zero Degree Calorimeter
194   // Dummy routine
195   //
196   return 9999;
197 }
198
199 //____________________________________________________________________________
200 Float_t AliZDC::ZMin(void) const
201 {
202   // Minimum dimension of the ZDC module in z
203   return -11600.;
204 }
205
206 //____________________________________________________________________________
207 Float_t AliZDC::ZMax(void) const
208 {
209   // Maximum dimension of the ZDC module in z
210   return -11750.;
211 }
212   
213
214 //_____________________________________________________________________________
215 void AliZDC::MakeBranch(Option_t *opt, const char * /*file*/)
216 {
217   //
218   // Create Tree branches for the ZDC
219   //
220
221   char branchname[10];
222   sprintf(branchname,"%s",GetName());
223
224   const char *cH = strstr(opt,"H");
225   
226   if (cH && fLoader->TreeH())
227    fHits   = new TClonesArray("AliZDCHit",1000); 
228   
229   AliDetector::MakeBranch(opt);
230 }
231
232 //_____________________________________________________________________________
233 void AliZDC::Hits2SDigits()
234 {
235   // Create summable digits from hits
236   
237   if (GetDebug()) printf("\n    Entering AliZDC::Hits2Digits() ");
238   
239   fLoader->LoadHits("read");
240   fLoader->LoadSDigits("recreate");
241   AliRunLoader* runLoader = fLoader->GetRunLoader();
242   AliZDCSDigit sdigit;
243   AliZDCSDigit* psdigit = &sdigit;
244
245   // Event loop
246   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
247     Float_t pmCZN = 0, pmCZP = 0, pmQZN[4], pmQZP[4], pmZEM1 = 0, pmZEM2 = 0;
248     for (Int_t i = 0; i < 4; i++) pmQZN[i] = pmQZP[i] = 0;
249
250     runLoader->GetEvent(iEvent);
251     TTree* treeH = fLoader->TreeH();
252     Int_t ntracks = (Int_t) treeH->GetEntries();
253     ResetHits();
254
255     // Tracks loop
256     Int_t sector[2];
257     for (Int_t itrack = 0; itrack < ntracks; itrack++) {
258       treeH->GetEntry(itrack);
259       for (AliZDCHit* zdcHit = (AliZDCHit*)FirstHit(-1); zdcHit;
260                       zdcHit = (AliZDCHit*)NextHit()) { 
261                       
262         sector[0] = zdcHit->GetVolume(0);
263         sector[1] = zdcHit->GetVolume(1);
264         if ((sector[1] < 1) || (sector[1] > 4)) {
265           Error("Hits2SDigits", "sector[0] = %d, sector[1] = %d", 
266                 sector[0], sector[1]);
267           continue;
268         }
269         Float_t lightQ = zdcHit->GetLightPMQ();
270         Float_t lightC = zdcHit->GetLightPMC();
271      
272         if (sector[0] == 1) {          //ZN 
273           pmCZN += lightC;
274           pmQZN[sector[1]-1] += lightQ;
275         } else if (sector[0] == 2) {   //ZP 
276           pmCZP += lightC;
277           pmQZP[sector[1]-1] += lightQ;
278         } else if (sector[0] == 3) {   //ZEM 
279           if (sector[1] == 1) pmZEM1 += lightC;
280           else                pmZEM2 += lightQ;
281         }
282       }//Hits loop
283     }
284
285     // create the output tree
286     fLoader->MakeTree("S");
287     TTree* treeS = fLoader->TreeS();
288     const Int_t kBufferSize = 4000;
289     treeS->Branch(GetName(), "AliZDCSDigit", &psdigit, kBufferSize);
290
291     // Create sdigits for ZN
292     sector[0] = 1; // Detector = ZN
293     sector[1] = 0; // Common PM ADC
294     new(psdigit) AliZDCSDigit(sector, pmCZN);
295     if (pmCZN > 0) treeS->Fill();
296     for (Int_t j = 0; j < 4; j++) {
297       sector[1] = j+1; // Towers PM ADCs
298       new(psdigit) AliZDCSDigit(sector, pmQZN[j]);
299       if (pmQZN[j] > 0) treeS->Fill();
300     }
301   
302     // Create sdigits for ZP
303     sector[0] = 2; // Detector = ZP
304     sector[1] = 0; // Common PM ADC
305     new(psdigit) AliZDCSDigit(sector, pmCZP);
306     if (pmCZP > 0) treeS->Fill();
307     for (Int_t j = 0; j < 4; j++) {
308       sector[1] = j+1; // Towers PM ADCs
309       new(psdigit) AliZDCSDigit(sector, pmQZP[j]);
310       if (pmQZP[j] > 0) treeS->Fill();
311     }
312
313     // Create sdigits for ZEM
314     sector[0] = 3; 
315     sector[1] = 1; // Detector = ZEM1
316     new(psdigit) AliZDCSDigit(sector, pmZEM1);
317     if (pmZEM1 > 0) treeS->Fill();
318     sector[1] = 2; // Detector = ZEM2
319     new(psdigit) AliZDCSDigit(sector, pmZEM2);
320     if (pmZEM2 > 0) treeS->Fill();
321
322     // write the output tree
323     fLoader->WriteSDigits("OVERWRITE");
324   }
325
326   fLoader->UnloadHits();
327   fLoader->UnloadSDigits();
328 }
329
330 //_____________________________________________________________________________
331 AliDigitizer* AliZDC::CreateDigitizer(AliRunDigitizer* manager) const
332 {
333   // Create the digitizer for ZDC
334
335   return new AliZDCDigitizer(manager);
336 }
337
338 //_____________________________________________________________________________
339 void AliZDC::Digits2Raw()
340 {
341   // Convert ZDC digits to raw data
342
343   // preliminary format: 12 interger values (ZNC, ZNQ1-4, ZPC, ZPQ1-4, ZEM1,2)
344   // For the CAEN module V965 we have an header, the Data Words and an End Of Block
345   UInt_t ADCHeader; 
346   UInt_t ADCData[24];
347   UInt_t ADCEndBlock;
348
349   // load the digits
350   fLoader->LoadDigits("read");
351   AliZDCDigit digit;
352   AliZDCDigit* pdigit = &digit;
353   TTree* treeD = fLoader->TreeD();
354   if (!treeD) return;
355   treeD->SetBranchAddress("ZDC", &pdigit);
356
357   // Fill data array
358   // ADC header
359   UInt_t ADCHeaderGEO = 0;
360   UInt_t ADCHeaderCRATE = 0;
361   UInt_t ADCHeaderCNT = (UInt_t) treeD->GetEntries();
362     
363   ADCHeader = ADCHeaderGEO << 27 | 0x1 << 25 | ADCHeaderCRATE << 16 |
364               ADCHeaderCNT << 8 ;
365
366   //printf("ADCHeader = %d\n",ADCHeader);
367       
368   // ADC data word
369   UInt_t ADCDataGEO = ADCHeaderGEO;
370   UInt_t ADCDataValue[24];
371   UInt_t ADCDataOvFlw[24];
372   for(Int_t i = 0; i < 24; i++){
373     ADCDataValue[i] = 0;
374     ADCDataOvFlw[i] = 0;
375   }
376   UInt_t ADCDataChannel = 0;
377   
378   // loop over digits
379   for (Int_t iDigit = 0; iDigit < treeD->GetEntries(); iDigit++) {
380     treeD->GetEntry(iDigit);
381     if (!pdigit) continue;
382
383     //ADC data
384     Int_t index = 0;
385     if(digit.GetSector(0)!=3){
386       index = (digit.GetSector(0)-1) + digit.GetSector(1)*4;
387       ADCDataChannel = (digit.GetSector(0)-1)*8 + digit.GetSector(1);
388     }
389     else {
390       index = 19 + digit.GetSector(1);
391       ADCDataChannel = 5 + digit.GetSector(1)*8;
392     }
393      
394     if ((index < 0) || (index >= 22)) {
395       Error("Digits2Raw", "sector[0] = %d, sector[1] = %d", 
396             digit.GetSector(0), digit.GetSector(1));
397       continue;
398     }
399     
400     ADCDataValue[index] = digit.GetADCValue(0);
401     if (ADCDataValue[index] > 2047) ADCDataOvFlw[index] = 1; 
402     ADCDataValue[index+2] = digit.GetADCValue(1);
403     if (ADCDataValue[index+2] > 2047) ADCDataOvFlw[index+2] = 1; 
404     
405     ADCData[index] =   ADCDataGEO << 27 | ADCDataChannel << 17 | 
406                        ADCDataOvFlw[index] << 12 | (ADCDataValue[index] & 0xfff); 
407     ADCData[index+2] = ADCDataGEO << 27 | ADCDataChannel << 17 | 0x1 << 16 |
408                        ADCDataOvFlw[index+2] << 12 | (ADCDataValue[index+2] & 0xfff);                    
409   }
410   //for (Int_t i=0;i<24;i++)printf("ADCData[%d] = %d\n",i,ADCData[i]);
411   
412   // End of Block
413   UInt_t ADCEndBlockGEO = ADCHeaderGEO;
414   UInt_t ADCEndBlockEvCount = gAlice->GetEventNrInRun();
415   
416   ADCEndBlock = ADCEndBlockGEO << 27 | 0x1 << 26 | ADCEndBlockEvCount;
417   
418   //printf("ADCEndBlock = %d\n",ADCEndBlock);
419
420
421   // open the output file
422   char fileName[30];
423   sprintf(fileName, "ZDC_%d.ddl", AliZDCRawStream::kDDLOffset);
424 #ifndef __DECCXX
425   ofstream file(fileName, ios::binary);
426 #else
427   ofstream file(fileName);
428 #endif
429
430   // write the DDL data header
431   AliRawDataHeader header;
432   header.fSize = sizeof(header) + sizeof(ADCHeader) + 
433                  sizeof(ADCData) + sizeof(ADCEndBlock);
434   //printf("sizeof header = %d, ADCHeader = %d, ADCData = %d, ADCEndBlock = %d\n",
435   //        sizeof(header),sizeof(ADCHeader),sizeof(ADCData),sizeof(ADCEndBlock));
436   header.SetAttribute(0);  // valid data
437   file.write((char*)(&header), sizeof(header));
438
439   // write the raw data and close the file
440   file.write((char*) &ADCHeader, sizeof (ADCHeader));
441   file.write((char*)(ADCData), sizeof(ADCData));
442   file.write((char*) &ADCEndBlock, sizeof(ADCEndBlock));
443   file.close();
444
445   // unload the digits
446   fLoader->UnloadDigits();
447 }
448
449 //______________________________________________________________________
450 void AliZDC::SetTreeAddress(){
451   // Set branch address for the Trees.
452   // Inputs:
453   //      none.
454   // Outputs:
455   //      none.
456   // Return:
457   //      none.
458   if (fLoader->TreeH() && (fHits == 0x0))
459     fHits   = new TClonesArray("AliZDCHit",1000);
460       
461   AliDetector::SetTreeAddress();
462 }
463  
464  
465 //Calibration methods (by Alberto Colla)
466  
467  
468 //________________________________________________________________
469 void AliZDC::CreateCalibData()
470 {
471   // 
472   //if (fCalibData) delete fCalibData; // delete previous version
473   fCalibData = new AliZDCCalibData(GetName());
474 }
475 //________________________________________________________________
476 void AliZDC::WriteCalibData(Int_t option)
477 {
478   //
479   const int kCompressLevel = 9;
480   const char defname[] = "$(ALICE)/AliRoot/data/AliZDCCalib.root";
481   char* fnam = GetZDCCalibFName();
482   if (!fnam || fnam[0]=='\0') {
483     fnam = gSystem->ExpandPathName(defname);
484     Warning("WriteCalibData","No File Name is provided, using default %s",fnam);
485   }
486   TFile* cdfile = TFile::Open(fnam,"UPDATE","",kCompressLevel);
487
488   // Writes Calibration Data to current directory. 
489   // User MUST take care of corresponding file opening and ->cd()... !!!
490   // By default, the object is overwritten. Use 0 option for opposite.
491   if (option) option = TObject::kOverwrite;
492   if (fCalibData) fCalibData->Write(0,option);
493   else if (fCalibData) fCalibData->Write(0,option);
494
495   cdfile->Close();
496   delete cdfile;
497 }
498
499 //________________________________________________________________
500 void AliZDC::LoadCalibData()
501 {
502   //
503   char* fnam = GetZDCCalibFName();
504   if (!fnam || fnam[0]=='\0') return; 
505   if (!gAlice->IsFileAccessible(fnam)) {
506     Error("LoadCalibData","ZDC Calibration Data file is not accessible, %s",fnam);
507     exit(1);
508   }
509   TFile* cdfile = TFile::Open(fnam);
510
511   // Loads Calibration Data from current directory. 
512   // User MUST take care of corresponding file opening and ->cd()...!!!
513   //
514   if (fCalibData) delete fCalibData; // delete previous version
515   TString dtname = "Calib_";
516   dtname += GetName();
517   fCalibData = (AliZDCCalibData*) gDirectory->Get(dtname.Data());
518   if (!fCalibData) { 
519     Error("LoadCalibData","No Calibration data found for %s",GetName());
520     exit(1);
521   }
522
523   cdfile->Close();
524   delete cdfile;
525 }
526
527
528 //Calibration methods (by Alberto Colla)