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