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