]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0RawData.cxx
DQM shifter histos based only on physics data
[u/mrichter/AliRoot.git] / T0 / AliT0RawData.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 //  T0 raw data conversion class                                            //
21 //                                                                           //
22 ///////////////////////////////////////////////////////////////////////////////
23
24 //#include <Riostream.h>
25 //#include <TTree.h>
26 #include <TMap.h>
27 #include "AliT0.h"
28 #include "AliT0RawData.h"
29 #include "AliT0digit.h"
30 #include "AliBitPacking.h"
31 #include "AliRawDataHeader.h"
32 #include "AliRawDataHeaderSim.h"
33 #include "AliBitPacking.h"
34 #include "AliFstream.h"
35 #include "AliRunLoader.h"
36 #include "AliDAQ.h"
37 #include "AliT0LookUpValue.h"
38 #include "AliT0LookUpKey.h"
39
40 using std::cout;
41 using std::endl;
42 ClassImp(AliT0RawData)
43
44 //_____________________________________________________________________________
45   AliT0RawData::AliT0RawData():TObject(),
46                                fVerbose(0),      
47                                fIndex(-1) ,     
48                                fEventNumber(0), 
49                                fTimeCFD(new TArrayI(24)),    
50                                fADC1( new TArrayI(24)),     
51                                fTimeLED( new TArrayI(24)), 
52                                fADC0( new TArrayI(24)),     
53                                fFile(0x0),   
54                                fDataHeaderPos(0),
55                                fDRMDataHeaderPos(0),
56                                fTRMDataHeaderPos(0),
57                                fParam(0),
58                                fLookUp(0)
59   
60
61 {
62   /*
63 -  48 channels (2 words each as in TOF DDL) for :
64 word 1 :0-5bit number of PMT; word 2: 0-7 error sign, 8-31 TDC
65 and the same but for amplified signal. Now I wrote the same time because
66 CDF are not ready and differences didn't measured yet.
67
68 -  48 channel for amplitude: very preliminary, QTC features are not
69 known now, preliminary i put as T1 time signal for this PMT in first
70 channel and T1+A in second, where A=Log(Amplitude);
71 and the same for amplified but A=Log(10*Amplitude).
72
73 - T0-A and T0-C 2 channels
74 - T0A-T0C vertex information
75 - Time Meaner where T0C TOF increase to the T0A TOF distance
76 - 6 multiplicity signals the same way as amplitude and with the same
77 uncertances
78   */
79
80   //open the output file
81   // char fileName[15];
82   // strcpy(fileName,AliDAQ::DdlFileName("T0",0)); //The name of the output file
83   TString fileName = Form("%s",AliDAQ::DdlFileName("T0",0));
84   fFile = new AliFstream(fileName.Data());
85   memset(fBuffer,0,512*sizeof(UInt_t));
86
87   //get event number 
88   AliRunLoader *runloader = AliRunLoader::Instance();
89   if (runloader) {
90     fEventNumber = runloader->GetEventNumber();
91   }
92
93   // Inverse lookup table for simulation
94
95   fParam = AliT0Parameters::Instance();
96   fParam->Init();
97   AliT0LookUpKey* lookkey;//= new AliT0LookUpKey();
98   AliT0LookUpValue*  lookvalue;//= new AliT0LookUpValue();
99   TMap *lookup = fParam->GetMapLookup();
100   TMapIter iter(lookup);
101
102   for( Int_t iline=0; iline<106; iline++)
103     {
104       lookvalue = ( AliT0LookUpValue*) iter.Next();
105       lookkey = (AliT0LookUpKey*) lookup->GetValue(lookvalue);
106       fLookUp.Add(lookkey, lookvalue);
107       //lookkey= new AliT0LookUpKey();
108       //lookvalue= new AliT0LookUpValue();
109     }
110     
111 }
112
113 //_____________________________________________________________________________
114
115 AliT0RawData::AliT0RawData(const AliT0RawData &r):TObject(),
116                                                   fVerbose(0),      
117                                                   fIndex(-1) ,     
118                                                   fEventNumber(0), 
119                                                   fTimeCFD(new TArrayI(24)),    
120                                                   fADC1( new TArrayI(24)),     
121                                                   fTimeLED( new TArrayI(24)), 
122                                                   fADC0( new TArrayI(24)),     
123                                                   fFile(0x0),   
124                                                   fDataHeaderPos(0),
125                                                   fDRMDataHeaderPos(0),
126                                                   fTRMDataHeaderPos(0),
127                                                   fParam(0),
128                                                   fLookUp(0)
129
130 {
131   //
132   // AliT0rawData copy constructor
133   //
134
135   ((AliT0RawData &) r).Copy(*this);
136
137 }
138
139 //_____________________________________________________________________________
140 AliT0RawData::~AliT0RawData()
141 {
142   //
143   // Destructor
144   //
145   delete fTimeCFD;
146   delete fADC1;
147   delete fTimeLED;
148   delete fADC0;
149   delete fFile;
150 }
151
152 //_____________________________________________________________________________
153 AliT0RawData &AliT0RawData::operator=(const AliT0RawData &r)
154 {
155   //
156   // Assignment operator
157   //
158
159   if (this != &r) ((AliT0RawData &) r).Copy(*this);
160   return *this;
161
162 }
163
164 //_____________________________________________________________________________
165 void AliT0RawData::GetDigits(AliT0digit *fDigits)
166 //void AliT0RawData::GetDigits(fDigits)
167 {
168  
169   //This method packs the T0 digits in a proper 32 bits structure
170
171   //read T0 digits and fill TDC and ADC arrays
172
173
174   //  Int_t error=0;
175   Int_t time,  positionOfTRMHeader;
176   
177   // Get the digits array
178   
179   fDigits->GetTimeCFD(*fTimeCFD);
180   fDigits->GetQT0(*fADC1);
181   fDigits->GetTimeLED(*fTimeLED);
182   fDigits->GetQT1(*fADC0);
183   Int_t meantime = fDigits->MeanTime(); 
184   Int_t timediff = fDigits->TimeDiff(); 
185   Int_t mult0=fDigits->SumMult();
186   Int_t mult1=fDigits->SumMult();
187   Int_t timeA = fDigits->BestTimeC();
188   Int_t timeC = fDigits->BestTimeA();
189   
190   
191   //  TArrayI  *allData = new TArrayI(110);
192   Int_t allData[110][1];
193   for (Int_t i=0; i<110; i++) allData[i][0] = 0;
194
195   allData[0][0]=0;
196   for (Int_t i=1; i<13; i++) {
197     allData[i][0]    = fTimeCFD->At(i-1);
198     allData[i+12][0] = fTimeLED->At(i-1);
199     allData[i+56][0] = fTimeCFD->At(i-1+12);
200     allData[i+68][0] = fTimeLED->At(i-1+12);
201   }
202   
203   for (Int_t iii=0; iii<12; iii++) {
204     allData[2*iii+25][0] = fADC1->At(iii);
205     allData[2*iii+26][0] = fADC0->At(iii);
206   }
207   for (Int_t ii=12; ii<24; ii++) {
208     allData[2*ii+57][0] = fADC1->At(ii);
209     allData[2*ii+58][0] = fADC0->At(ii);
210   }
211   
212   allData[49][0] = meantime;
213   allData[50][0] = timediff;
214   allData[51][0] = timeA;
215   allData[52][0] = timeC;
216   allData[53][0] = mult0;
217   allData[54][0] = mult1;
218   allData[55][0] = mult0;
219   allData[56][0] = mult1;
220
221   //    cout.setf( ios_base::hex, ios_base::basefield );
222   //space for DRM header
223   fIndex += 6;
224
225
226   Int_t startTRM=fIndex;
227   //space for 1st TRM header
228   fIndex ++;
229   positionOfTRMHeader= fIndex;
230   //space for chain  header
231   fIndex ++;
232   WriteChainDataHeader(1, 1); // 
233
234   //  fIndex++;
235   // Loop through all PMT
236   Int_t chain=0; 
237   Int_t iTDC = 0;
238   Int_t channel=0;
239   Int_t trm1words=0;
240   Int_t itrm=7;
241   Int_t inside =0;
242   Int_t isData = 0;
243   AliT0LookUpKey * lookkey  = new AliT0LookUpKey();
244   AliT0LookUpValue * lookvalue ;//= new AliT0LookUpValue(trm,tdc,chain,channel);
245   for (Int_t det = 0; det < 105; det++) {
246     time = allData[det][0];
247     if (time >0 && time !=99999) {
248       lookkey->SetKey(det);
249       lookvalue = (AliT0LookUpValue*) fLookUp.GetValue((TObject*)lookkey);     
250       if (lookvalue ) 
251         {
252           isData++;
253           itrm= lookvalue->GetTRM();
254           if (det >56 &&inside == 0)  {
255             WriteChainDataTrailer(1); // 1st chain trailer
256             fIndex++;
257             WriteChainDataHeader(2, 1);
258             //      fIndex++;
259             inside++;
260           }         
261           chain = lookvalue->GetChain();
262           iTDC = lookvalue->GetTDC();
263           channel = lookvalue->GetChannel();
264           FillTime(channel,  iTDC,  time);
265           AliDebug(1,Form("look %i  itrm %i ,  chain %i , iTDC %i, channel %i time %i", det,itrm,chain,iTDC,channel, time));
266         }
267       else
268         {
269           cout<<" no lookup value for key "<<det<<endl;
270           //  break;
271         }
272     }
273     
274   }
275   if (inside==0) {
276     WriteChainDataTrailer(1); // 1st chain trailer
277     fIndex++;
278     WriteChainDataHeader(2, 1);
279   }
280     //  WriteChainDataHeader(2, 1); // 
281   WriteChainDataTrailer(3); // 2st chain trailer
282   WriteTrailer(15,0,fEventNumber,5); // 1st TRM trailer
283   
284   
285   trm1words = fIndex - startTRM;
286   //space for 2st TRM header
287   
288   WriteTRMDataHeader(itrm, trm1words , positionOfTRMHeader);
289   
290   //DRM trailer
291   WriteTrailer(1,0,fEventNumber,5);
292     
293     WriteDRMDataHeader();
294     
295 }
296
297 //_____________________________________________________________________________
298
299 void  AliT0RawData::WriteDRMDataHeader()
300 {
301 //Write a (dummy or real) DDL DRM  data header, 
302 //set the compression bit if compressed
303 //  UInt_t drmheader[4];  
304   UInt_t word;
305   UInt_t baseWord=0;
306   //fill DRM headers
307   //DRM Global Header
308   word = 1;
309   AliBitPacking::PackWord(word,baseWord,0, 3); // 0001 
310   word = fIndex ;
311
312   AliBitPacking::PackWord(word,baseWord,4, 20); // event words 
313   word=124;
314   AliBitPacking::PackWord(word,baseWord, 21, 27); // event words 
315   word=4;
316   AliBitPacking::PackWord(word,baseWord,28, 31);// 0100 marks header
317   fBuffer[0]=  baseWord;
318
319
320   //DRM status header 1
321   word = 1;
322   AliBitPacking::PackWord(word,baseWord,0, 3); // 0001 
323   word = 1;
324   AliBitPacking::PackWord(word,baseWord,4, 14); // slotID now 0000000001
325   word = 1;
326   AliBitPacking::PackWord(word,baseWord,15, 15); //if 1  LHC clock is coorectly recieved from CPDM 
327   word=0;
328   AliBitPacking::PackWord(word,baseWord,16,27); // reserve for future use
329   word=4;
330   AliBitPacking::PackWord(word,baseWord,28,31); // 0100 marks header
331   fBuffer[1] = baseWord;
332   
333    word=0;
334    baseWord=0;
335    
336    //DRM status header 2
337    word = 1;
338    AliBitPacking::PackWord(word,baseWord, 0, 3); // 0001 
339    word = 3;
340    AliBitPacking::PackWord(word,baseWord, 4, 14); //enable slotID now 00000000011
341    word = 0;
342    AliBitPacking::PackWord(word,baseWord, 15, 15); // something
343    word=0;
344    AliBitPacking::PackWord(word,baseWord, 16, 27); // fault ID for simulation 0
345    word=4;
346    AliBitPacking::PackWord(word,baseWord,28,31); // 0100 marks header
347    fBuffer[2]=  baseWord;
348         
349    word=0;
350    baseWord=0;
351    //DRM status header 3
352    word = 1;
353     AliBitPacking::PackWord(word,baseWord,0, 3); // 0001 
354     word = 0;
355     AliBitPacking::PackWord(word,baseWord,4, 27); // TTC event counter
356     word=4;
357     AliBitPacking::PackWord(word,baseWord,28,31); // 0100 marks header
358     fBuffer[3]=  baseWord;
359
360     // new DRM format
361     fBuffer[4]=  baseWord;
362     fBuffer[5]=  baseWord;
363    
364     word=0;
365     baseWord=0;
366     
367 }
368   
369 //_____________________________________________________________________________
370
371 void  AliT0RawData::WriteTRMDataHeader(UInt_t slotID, Int_t nWordsInTRM,
372                                           Int_t  positionOfTRMHeader)
373 {
374 //Write a (dummy or real) DDL TRM  data header, 
375 //set the compression bit if compressed
376 //  UInt_t trmheader;  
377   UInt_t word;
378   UInt_t baseWord=0;
379   //fill TRM headers
380   //TRM Global Header
381   word = slotID;
382   AliBitPacking::PackWord(word,baseWord,0, 3); // slotID
383   word = nWordsInTRM;
384  //+this word - DRM header 
385
386   AliBitPacking::PackWord(word,baseWord,4, 16); // event words 
387   word=0;
388   AliBitPacking::PackWord(word,baseWord,17,18); // ACQ
389   word=0;
390   AliBitPacking::PackWord(word,baseWord,19,19); //  L SEY inside LUT
391   word=0;
392   AliBitPacking::PackWord(word,baseWord,20,27); //  MBZ
393   word=4;
394   AliBitPacking::PackWord(word,baseWord,28,31); // 0100 marks header
395   fBuffer[positionOfTRMHeader] =  baseWord;
396
397   word=0; 
398   baseWord=0;
399      
400 }
401
402 //_____________________________________________________________________________
403
404 void  AliT0RawData::WriteChainDataHeader(UInt_t chainNumber,UInt_t slotID)
405 {
406 //Write a (dummy or real) DDL Chain  data header, 
407 //set the compression bit if compressed
408 //  chainNumber 00 or 10
409   UInt_t word;
410   UInt_t baseWord=0;
411   //fill TRM headers
412   //TRM Global Header
413   word = slotID; // ask Tatiana 7 or 9 
414   AliBitPacking::PackWord(word,baseWord,0, 3); // slotID
415   word = 0;
416   AliBitPacking::PackWord(word,baseWord,4, 15); // bunchID
417   word=0;
418   AliBitPacking::PackWord(word,baseWord,16,23); // PB24 temperature
419   word=0;
420   AliBitPacking::PackWord(word,baseWord,24,26); //  PB24 ID
421   word=0;
422   AliBitPacking::PackWord(word,baseWord,27,27); //  TS
423   word=chainNumber;
424   AliBitPacking::PackWord(word,baseWord,28,31); // 0100 marks header
425   fBuffer[fIndex] =  baseWord;
426   //cout<<" chain header "<<baseWord<<" number "<<chainNumber<<endl;
427   word=0;
428   baseWord=0;     
429   
430 }
431 //_____________________________________________________________________________
432
433 void  AliT0RawData::WriteChainDataTrailer(UInt_t chainNumber )
434 {
435 //Write a (dummy or real) DDL Chain  data trailer 
436 //set the compression bit if compressed
437 //  chainNumber 00 or 10
438   UInt_t word;
439   UInt_t baseWord=0;
440   word = 0; // ask Tatiana 7 or 9 
441   AliBitPacking::PackWord(word,baseWord,0, 3); // status
442   word = 0;
443   AliBitPacking::PackWord(word,baseWord,4, 15); // MBZ
444   word=fEventNumber;
445   AliBitPacking::PackWord(word,baseWord,16,27); // event counter
446   word=chainNumber;
447   AliBitPacking::PackWord(word,baseWord,28,31); // chain number
448   fIndex++;
449   fBuffer[fIndex] =  baseWord;
450
451   word=0;
452   baseWord=0;     
453   
454 }
455 //_____________________________________________________________________________
456
457 void  AliT0RawData::WriteDataHeader(Bool_t dummy, Bool_t compressed)
458 {
459 //Write a (dummy or real) DDL data header, 
460 //set the compression bit if compressed
461
462   AliRawDataHeaderSim header;
463   
464   if (dummy) {
465     //if size=0 it means that this data header is a dummy data header
466     fDataHeaderPos = fFile->Tellp();
467     fFile->WriteBuffer((char*)(&header), sizeof(header));
468   } else {
469     UInt_t currentFilePos = fFile->Tellp();
470     fFile->Seekp(fDataHeaderPos);
471     header.fSize = currentFilePos-fDataHeaderPos;
472     header.SetAttribute(0);  // valid data
473     if (compressed) header.SetAttribute(1); 
474     fFile->WriteBuffer((char*)(&header), sizeof(header));
475     fFile->Seekp(currentFilePos);
476   }
477   
478 }
479
480 //___ __________________________________________________________________________
481
482
483 void  AliT0RawData::WriteTrailer(UInt_t slot, Int_t word1, UInt_t word2, UInt_t word3)
484 {
485 //Write a (dummy or real) DDL Chain  data trailer 
486
487   UInt_t word;
488   UInt_t baseWord=0;
489   word = slot;
490   AliBitPacking::PackWord(word,baseWord,0, 3); // 0001 
491   word=word1;
492   AliBitPacking::PackWord(word,baseWord,4, 15); // CRC ?
493   word = word2;
494   AliBitPacking::PackWord(word,baseWord,16,27); // event counter
495   word=word3;
496   AliBitPacking::PackWord(word,baseWord,28,31); //  marks trailer
497   fIndex++;
498   fBuffer[fIndex] =  baseWord;
499
500   word=0;
501   baseWord=0;
502
503 }
504
505 //---------------------------------------------------------------------------------------
506 //---------------------------------------------------------------------------------------
507 void  AliT0RawData::FillTime(Int_t ch, Int_t iTDC, Int_t time)
508 {
509   //  put all time fields thother in 1 word
510
511   UInt_t word;
512   UInt_t baseWord=0;
513
514   word=time;
515   AliBitPacking::PackWord(word,baseWord, 0, 20); // Time 
516
517   word=ch;
518   AliBitPacking::PackWord(word,baseWord, 21, 23); // number of channel 
519   word=iTDC;
520   AliBitPacking::PackWord(word,baseWord, 24, 27); // TDC ID
521
522   word=0;
523   AliBitPacking::PackWord(word,baseWord, 28, 28); // E = 0 in simulation
524   word=0;
525   AliBitPacking::PackWord(word,baseWord, 29, 30); // PS bit data 00
526   word=1;
527   AliBitPacking::PackWord(word,baseWord, 31, 31); // 1
528   fIndex++;
529   fBuffer[fIndex]=baseWord;
530
531   word=0;
532   baseWord=0;
533 }
534 //---------------------------------------------------------------------------------------
535
536 Int_t AliT0RawData::RawDataT0(AliT0digit *fDigits)
537   //Int_t AliT0RawData::RawDataT0(*fDigits)
538 {
539    //This method creates the Raw data files for T0 detector
540
541
542   // const Int_t kSize=512; //2*AliTOFGeometry::NpadXSector() 
543                           //max number of digits per DDL file times 2
544   //  UInt_t fBuffer[kSize];
545   //  UInt_t baseWord;
546   // UInt_t word;
547
548   fIndex=-1;
549  
550
551    WriteDataHeader(kTRUE, kFALSE);
552   GetDigits(fDigits);
553   //write packing digits
554   
555   
556   fFile->WriteBuffer((char*) fBuffer,((fIndex+1)*sizeof(UInt_t)));
557   //write real data header on its place
558    WriteDataHeader(kFALSE, kFALSE);
559   
560   
561   //end for
562   
563   return 0;  
564   
565 }