1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
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 //
25 ///////////////////////////////////////////////////////////////////////////////
29 #include <TGeometry.h>
34 // --- AliRoot header files
35 #include "AliDetector.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"
44 #include "AliRawDataHeader.h"
45 #include "AliLoader.h"
54 //_____________________________________________________________________________
58 // Default constructor for the Zero Degree Calorimeter base class
72 //_____________________________________________________________________________
73 AliZDC::AliZDC(const char *name, const char *title)
74 : AliDetector(name,title)
77 // Standard constructor for the Zero Degree Calorimeter base class
83 // Allocate the hits array
84 fHits = new TClonesArray("AliZDCHit",1000);
85 gAlice->GetMCApp()->AddHitList(fHits);
87 char sensname[5],senstitle[25];
88 sprintf(sensname,"ZDC");
89 sprintf(senstitle,"ZDC dummy");
90 SetName(sensname); SetTitle(senstitle);
98 //____________________________________________________________________________
109 //_____________________________________________________________________________
110 void AliZDC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
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
120 static Float_t primKinEn, xImpact, yImpact, sFlag;
122 AliZDCHit *newquad, *curprimquad;
123 newquad = new AliZDCHit(fIshunt, track, vol, hits);
124 TClonesArray &lhits = *fHits;
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
132 else if(track == primary){
133 newquad->SetSFlag(0); // PRIMARY particle entering the ZDC
135 sFlag = newquad->GetSFlag();
136 primKinEn = newquad->GetPrimKinEn();
137 xImpact = newquad->GetXImpact();
138 yImpact = newquad->GetYImpact();
141 newquad->SetPrimKinEn(primKinEn);
142 newquad->SetXImpact(xImpact);
143 newquad->SetYImpact(yImpact);
144 newquad->SetSFlag(sFlag);
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;
158 //Otherwise create a new hit
159 new(lhits[fNhits]) AliZDCHit(newquad);
165 //_____________________________________________________________________________
166 void AliZDC::BuildGeometry()
169 // Build the ROOT TNode geometry for event display
170 // in the Zero Degree Calorimeter
171 // This routine is dummy for the moment
176 const int kColorZDC = kBlue;
179 top=gAlice->GetGeometry()->GetNode("alice");
182 brik = new TBRIK("S_ZDC","ZDC box","void",300,300,5);
184 node = new TNode("ZDC","ZDC","S_ZDC",0,0,600,"");
185 node->SetLineColor(kColorZDC);
189 //_____________________________________________________________________________
190 Int_t AliZDC::DistancetoPrimitive(Int_t , Int_t )
193 // Distance from the mouse to the Zero Degree Calorimeter
199 //____________________________________________________________________________
200 Float_t AliZDC::ZMin(void) const
202 // Minimum dimension of the ZDC module in z
206 //____________________________________________________________________________
207 Float_t AliZDC::ZMax(void) const
209 // Maximum dimension of the ZDC module in z
214 //_____________________________________________________________________________
215 void AliZDC::MakeBranch(Option_t *opt, const char * /*file*/)
218 // Create Tree branches for the ZDC
222 sprintf(branchname,"%s",GetName());
224 const char *cH = strstr(opt,"H");
226 if (cH && fLoader->TreeH())
227 fHits = new TClonesArray("AliZDCHit",1000);
229 AliDetector::MakeBranch(opt);
232 //_____________________________________________________________________________
233 void AliZDC::Hits2SDigits()
235 // Create summable digits from hits
237 if (GetDebug()) printf("\n Entering AliZDC::Hits2Digits() ");
239 fLoader->LoadHits("read");
240 fLoader->LoadSDigits("recreate");
241 AliRunLoader* runLoader = fLoader->GetRunLoader();
243 AliZDCSDigit* psdigit = &sdigit;
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;
250 runLoader->GetEvent(iEvent);
251 TTree* treeH = fLoader->TreeH();
252 Int_t ntracks = (Int_t) treeH->GetEntries();
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()) {
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]);
269 Float_t lightQ = zdcHit->GetLightPMQ();
270 Float_t lightC = zdcHit->GetLightPMC();
272 if (sector[0] == 1) { //ZN
274 pmQZN[sector[1]-1] += lightQ;
275 } else if (sector[0] == 2) { //ZP
277 pmQZP[sector[1]-1] += lightQ;
278 } else if (sector[0] == 3) { //ZEM
279 if (sector[1] == 1) pmZEM1 += lightC;
280 else pmZEM2 += lightQ;
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);
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();
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();
313 // Create sdigits for ZEM
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();
322 // write the output tree
323 fLoader->WriteSDigits("OVERWRITE");
326 fLoader->UnloadHits();
327 fLoader->UnloadSDigits();
330 //_____________________________________________________________________________
331 AliDigitizer* AliZDC::CreateDigitizer(AliRunDigitizer* manager) const
333 // Create the digitizer for ZDC
335 return new AliZDCDigitizer(manager);
338 //_____________________________________________________________________________
339 void AliZDC::Digits2Raw()
341 // Convert ZDC digits to raw data
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
350 fLoader->LoadDigits("read");
352 AliZDCDigit* pdigit = &digit;
353 TTree* treeD = fLoader->TreeD();
355 treeD->SetBranchAddress("ZDC", &pdigit);
359 UInt_t ADCHeaderGEO = 0;
360 UInt_t ADCHeaderCRATE = 0;
361 UInt_t ADCHeaderCNT = (UInt_t) treeD->GetEntries();
363 ADCHeader = ADCHeaderGEO << 27 | 0x1 << 25 | ADCHeaderCRATE << 16 |
366 //printf("ADCHeader = %d\n",ADCHeader);
369 UInt_t ADCDataGEO = ADCHeaderGEO;
370 UInt_t ADCDataValue[24];
371 UInt_t ADCDataOvFlw[24];
372 for(Int_t i = 0; i < 24; i++){
376 UInt_t ADCDataChannel = 0;
379 for (Int_t iDigit = 0; iDigit < treeD->GetEntries(); iDigit++) {
380 treeD->GetEntry(iDigit);
381 if (!pdigit) continue;
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);
390 index = 19 + digit.GetSector(1);
391 ADCDataChannel = 5 + digit.GetSector(1)*8;
394 if ((index < 0) || (index >= 22)) {
395 Error("Digits2Raw", "sector[0] = %d, sector[1] = %d",
396 digit.GetSector(0), digit.GetSector(1));
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;
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);
410 //for (Int_t i=0;i<24;i++)printf("ADCData[%d] = %d\n",i,ADCData[i]);
413 UInt_t ADCEndBlockGEO = ADCHeaderGEO;
414 UInt_t ADCEndBlockEvCount = gAlice->GetEventNrInRun();
416 ADCEndBlock = ADCEndBlockGEO << 27 | 0x1 << 26 | ADCEndBlockEvCount;
418 //printf("ADCEndBlock = %d\n",ADCEndBlock);
421 // open the output file
423 sprintf(fileName, "ZDC_%d.ddl", AliZDCRawStream::kDDLOffset);
425 ofstream file(fileName, ios::binary);
427 ofstream file(fileName);
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));
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));
446 fLoader->UnloadDigits();
449 //______________________________________________________________________
450 void AliZDC::SetTreeAddress(){
451 // Set branch address for the Trees.
458 if (fLoader->TreeH() && (fHits == 0x0))
459 fHits = new TClonesArray("AliZDCHit",1000);
461 AliDetector::SetTreeAddress();
465 //Calibration methods (by Alberto Colla)
468 //________________________________________________________________
469 void AliZDC::CreateCalibData()
472 //if (fCalibData) delete fCalibData; // delete previous version
473 fCalibData = new AliZDCCalibData(GetName());
475 //________________________________________________________________
476 void AliZDC::WriteCalibData(Int_t option)
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);
486 TFile* cdfile = TFile::Open(fnam,"UPDATE","",kCompressLevel);
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);
499 //________________________________________________________________
500 void AliZDC::LoadCalibData()
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);
509 TFile* cdfile = TFile::Open(fnam);
511 // Loads Calibration Data from current directory.
512 // User MUST take care of corresponding file opening and ->cd()...!!!
514 if (fCalibData) delete fCalibData; // delete previous version
515 TString dtname = "Calib_";
517 fCalibData = (AliZDCCalibData*) gDirectory->Get(dtname.Data());
519 Error("LoadCalibData","No Calibration data found for %s",GetName());
528 //Calibration methods (by Alberto Colla)