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 // TRD raw data conversion class //
22 ///////////////////////////////////////////////////////////////////////////////
29 #include "AliRawDataHeaderSim.h"
30 #include "AliRawReader.h"
32 #include "AliFstream.h"
34 #include "AliTRDrawData.h"
35 #include "AliTRDdigitsManager.h"
36 #include "AliTRDgeometry.h"
37 #include "AliTRDarrayDictionary.h"
38 #include "AliTRDarrayADC.h"
39 #include "AliTRDrawStreamBase.h"
40 #include "AliTRDrawOldStream.h"
41 #include "AliTRDcalibDB.h"
42 #include "AliTRDSignalIndex.h"
43 #include "AliTRDfeeParam.h"
44 #include "AliTRDmcmSim.h"
46 ClassImp(AliTRDrawData)
48 Int_t AliTRDrawData::fgRawFormatVersion = AliTRDrawData::kRawNewFormat;
49 Int_t AliTRDrawData::fgDataSuppressionLevel = 1;
51 //_____________________________________________________________________________
52 AliTRDrawData::AliTRDrawData()
59 ,fTrackletContainer(NULL)
65 // Default constructor
68 fFee = AliTRDfeeParam::Instance();
69 fNumberOfDDLs = AliDAQ::NumberOfDdls("TRD");
73 //_____________________________________________________________________________
74 AliTRDrawData::AliTRDrawData(const AliTRDrawData &r)
81 ,fTrackletContainer(NULL)
90 fFee = AliTRDfeeParam::Instance();
91 fNumberOfDDLs = AliDAQ::NumberOfDdls("TRD");
95 //_____________________________________________________________________________
96 AliTRDrawData::~AliTRDrawData()
102 if (fTrackletContainer){
103 delete fTrackletContainer;
104 fTrackletContainer = NULL;
109 //_____________________________________________________________________________
110 Bool_t AliTRDrawData::Digits2Raw(TTree *digitsTree, TTree *tracks )
113 // Initialize necessary parameters and call one
114 // of the raw data simulator selected by SetRawVersion.
116 // Currently tracklet output is not spported yet and it
117 // will be supported in higher version simulator.
120 AliTRDdigitsManager* digitsManager = new AliTRDdigitsManager();
122 if (!digitsManager->ReadDigits(digitsTree)) {
123 delete digitsManager;
127 if (tracks != NULL) {
128 delete digitsManager;
129 AliError("Tracklet input is not supported yet.");
133 fGeo = new AliTRDgeometry();
135 if (!AliTRDcalibDB::Instance()) {
136 AliError("Could not get calibration object");
138 delete digitsManager;
142 Int_t retval = kTRUE;
143 Int_t rv = fFee->GetRAWversion();
145 // Call appropriate Raw Simulator
146 if ( rv > 0 && rv <= 3 ) retval = Digits2Raw(digitsManager);
149 AliWarning(Form("Unsupported raw version (%d).", rv));
154 delete digitsManager;
160 //_____________________________________________________________________________
161 Bool_t AliTRDrawData::Digits2Raw(AliTRDdigitsManager *digitsManager)
164 // Raw data simulator for all versions > 0. This is prepared for real data.
165 // This version simulate only raw data with ADC data and not with tracklet.
168 const Int_t kMaxHcWords = (fGeo->TBmax()/3)
174 // Buffer to temporary store half chamber data
175 UInt_t *hcBuffer = new UInt_t[kMaxHcWords];
177 Bool_t newEvent = kFALSE; // only for correct readout tree
178 Bool_t newSM = kFALSE; // new SM flag, for writing SM index words
179 Bool_t newStack = kFALSE; // new stack flag, for writing stack index words
181 // sect is same as iDDL, so I use only sect here.
182 for (Int_t sect = 0; sect < fGeo->Nsector(); sect++) {
185 sprintf(name,"TRD_%d.ddl",sect + AliTRDrawOldStream::kDDLOffset);
187 AliFstream* of = new AliFstream(name);
189 // Write a dummy data header
190 AliRawDataHeaderSim header; // the event header
191 UInt_t hpos = of->Tellp();
192 of->WriteBuffer((char *) (& header), sizeof(header));
194 // Reset payload byte size (payload does not include header).
195 Int_t npayloadbyte = 0;
198 if ( fgRawFormatVersion == 0 ){
199 // GTU common data header (5x4 bytes per super module, shows link mask)
200 for( Int_t stack = 0; stack < fGeo->Nstack(); stack++ ) {
201 UInt_t gtuCdh = (UInt_t)(0xe << 28);
202 for( Int_t layer = 0; layer < fGeo->Nlayer(); layer++) {
203 Int_t iDet = fGeo->GetDetector(layer, stack, sect);
205 // If chamber status is ok, we assume that the optical link is also OK.
206 // This is shown in the GTU link mask.
207 if ( AliTRDcalibDB::Instance()->GetChamberStatus(iDet) )
208 gtuCdh = gtuCdh | (3 << (2*layer));
210 of->WriteBuffer((char *) (& gtuCdh), sizeof(gtuCdh));
216 // check the existance of the data
217 // SM index word and Stack index word
218 if ( fgRawFormatVersion == 1 ){
219 UInt_t *iwbuffer = new UInt_t[42]; // index word buffer; max 42 = 2 SM headers + 5*8 stack headers
221 UInt_t StackMask = 0x0;
222 Bool_t StackHasData = kFALSE;
223 Bool_t SMHasData = kFALSE;
224 iwbuffer[nheader++] = 0x0001a020; // SM index words
225 iwbuffer[nheader++] = 0x10404071; // SM header
227 for ( Int_t stack= 0; stack < fGeo->Nstack(); stack++) {
228 UInt_t LinkMask = 0x0;
229 for( Int_t layer = 0; layer < fGeo->Nlayer(); layer++) {
230 Int_t iDet = fGeo->GetDetector(layer,stack,sect);
231 AliTRDarrayADC *digits = (AliTRDarrayADC *) digitsManager->GetDigits(iDet);
232 if ( digits->HasData() ) {
233 StackMask = StackMask | ( 1 << stack ); // active stack mask for new stack
234 LinkMask = LinkMask | ( 3 << (2*layer) ); // 3 = 0011
235 StackHasData = kTRUE;
240 if ( fgDataSuppressionLevel==0 || StackHasData ){
241 //if ( StackHasData ){
242 iwbuffer[nheader++] = 0x0007a000 | LinkMask; // stack index word + link masks
243 //if (fgDataSuppressionLevel==0) iwbuffer[nheader-1] = 0x0007afff; // no suppression
244 iwbuffer[nheader++] = 0x04045b01; // stack header
245 for (Int_t i=0;i<6;i++) iwbuffer[nheader++] = 0x00000000; // 6 dummy words
246 StackHasData = kFALSE;
250 if ( fgDataSuppressionLevel==0 || SMHasData ){
251 iwbuffer[0] = iwbuffer[0] | StackMask; // add stack masks to SM index word
252 if (fgDataSuppressionLevel==0) iwbuffer[0] = 0x0001a03f; // no suppression
253 of->WriteBuffer((char *) iwbuffer, nheader*4);
254 AliDebug(11, Form("SM %d index word: %08x", iwbuffer[0]));
255 AliDebug(11, Form("SM %d header: %08x", iwbuffer[1]));
258 // end of SM & stack header ------------------------------------------------------------------------
259 // -------------------------------------------------------------------------------------------------
261 // Prepare chamber data
262 for( Int_t stack = 0; stack < fGeo->Nstack(); stack++) {
263 for( Int_t layer = 0; layer < fGeo->Nlayer(); layer++) {
265 Int_t iDet = fGeo->GetDetector(layer,stack,sect);
267 newEvent = kTRUE; // it is expected that each event has at least one tracklet;
268 // this is only needed for correct readout tree
270 AliDebug(11, Form("New event!! Event counter: %d",fEventCounter));
273 if ( stack==0 && layer==0 ) newSM = kTRUE; // new SM flag
274 if ( layer==0 ) newStack = kTRUE; // new stack flag
275 AliDebug(15, Form("stack : %d, layer : %d, iDec : %d\n",stack,layer,iDet));
276 // Get the digits array
277 AliTRDarrayADC *digits = (AliTRDarrayADC *) digitsManager->GetDigits(iDet);
278 if (fgDataSuppressionLevel==0 || digits->HasData() ) { // second part is new!! and is for indicating a new event
280 if (digits->HasData()) digits->Expand();
283 Int_t rv = fFee->GetRAWversion();
286 if ( fgRawFormatVersion == 0 ){
287 // Process A side of the chamber
288 if ( rv >= 1 && rv <= 2 ) {
289 hcwords = ProduceHcDataV1andV2(digits,0,iDet,hcBuffer,kMaxHcWords);
293 hcwords = ProduceHcDataV3 (digits,0,iDet,hcBuffer,kMaxHcWords,newEvent);
294 //hcwords = ProduceHcDataV3 (digits,0,iDet,hcBuffer,kMaxHcWords);
295 if(newEvent == kTRUE) newEvent = kFALSE;
298 of->WriteBuffer((char *) hcBuffer, hcwords*4);
299 npayloadbyte += hcwords*4;
301 // Process B side of the chamber
302 if ( rv >= 1 && rv <= 2 ) {
303 hcwords = ProduceHcDataV1andV2(digits,1,iDet,hcBuffer,kMaxHcWords);
307 hcwords = ProduceHcDataV3 (digits,1,iDet,hcBuffer,kMaxHcWords,newEvent);
308 //hcwords = ProduceHcDataV3 (digits,1,iDet,hcBuffer,kMaxHcWords);
311 of->WriteBuffer((char *) hcBuffer, hcwords*4);
312 npayloadbyte += hcwords*4;
314 } else { // real data format
316 if (digits->HasData()){
317 // Process A side of the chamber
318 hcwords = ProduceHcData(digits,0,iDet,hcBuffer,kMaxHcWords,newEvent,newSM);
320 // AssignStackMask(hcBuffer, stack); // active stack mask for this stack
321 // hcwords += AddStackIndexWords(hcBuffer, stack, hcwords);
322 // newStack = kFALSE;
324 //if ( newSM ) newSM = kFALSE;
325 if ( newEvent ) newEvent = kFALSE;
326 //AssignLinkMask(hcBuffer, layer); // active link mask for this layer(2*HC)
327 of->WriteBuffer((char *) hcBuffer, hcwords*4);
328 npayloadbyte += hcwords*4;
329 //for ( Int_t i=0; i<hcwords; i++ ) AliInfo(Form("Buf : %X",hcBuffer[i]));
331 // Process B side of the chamber
332 hcwords = ProduceHcData(digits,1,iDet,hcBuffer,kMaxHcWords,newEvent,newSM);
333 of->WriteBuffer((char *) hcBuffer, hcwords*4);
334 npayloadbyte += hcwords*4;
336 hcBuffer[hcwords++] = fgkEndOfTrackletMarker;
337 hcBuffer[hcwords++] = fgkEndOfTrackletMarker;
338 hcBuffer[hcwords++] = (1<<31) | (0<<24) | (0<<17) | (1<<14) | (sect<<9) | (layer<<6) | (stack<<3) | (0<<2) | 1;
339 hcBuffer[hcwords++] = (24<<26) | (99<<10) | (15<<6) | (11<<2) | 1;
340 hcBuffer[hcwords++] = kEndofrawdatamarker;
341 hcBuffer[hcwords++] = kEndofrawdatamarker;
342 hcBuffer[hcwords++] = kEndofrawdatamarker;
343 hcBuffer[hcwords++] = kEndofrawdatamarker;
344 npayloadbyte += hcwords*4;
346 hcBuffer[hcwords++] = fgkEndOfTrackletMarker;
347 hcBuffer[hcwords++] = fgkEndOfTrackletMarker;
348 hcBuffer[hcwords++] = (1<<31) | (0<<24) | (0<<17) | (1<<14) | (sect<<9) | (layer<<6) | (stack<<3) | (1<<2) | 1;
349 hcBuffer[hcwords++] = (24<<26) | (99<<10) | (15<<6) | (11<<2) | 1;
350 hcBuffer[hcwords++] = kEndofrawdatamarker;
351 hcBuffer[hcwords++] = kEndofrawdatamarker;
352 hcBuffer[hcwords++] = kEndofrawdatamarker;
353 hcBuffer[hcwords++] = kEndofrawdatamarker;
354 npayloadbyte += hcwords*4;
356 of->WriteBuffer((char *) hcBuffer, hcwords*4);
366 header.fSize = UInt_t(of->Tellp()) - hpos;
367 header.SetAttribute(0); // Valid data
368 of->Seekp(hpos); // Rewind to header position
369 of->WriteBuffer((char *) (& header), sizeof(header));
379 //_____________________________________________________________________________
380 void AliTRDrawData::ProduceSMIndexData(UInt_t *buf, Int_t& nw){
382 // This function generates
383 // 1) SM index words : ssssssss ssssssss vvvv rrrr r d t mmmmm
384 // - s : size of SM header (number of header, default = 0x0001)
385 // - v : SM header version (default = 0xa)
386 // - r : reserved for future use (default = 00000)
387 // - d : track data enabled bit (default = 0)
388 // - t : tracklet data enabled bit (default = 1)
389 // - m : stack mask (each bit corresponds a stack, default = 11111)
391 // 2) SM header : rrr c vvvv vvvvvvvv vvvv rrrr bbbbbbbb
392 // - r : reserved for future use (default = 000)
393 // - c : clean check out flag (default = 1)
394 // - v : hardware design revision (default = 0x0404)
395 // - r : reserved for future use (default = 0x0)
396 // - b : physical board ID (default = 0x71)
398 // 3) stack index words : ssssssss ssssssss vvvv mmmm mmmmmmmm
399 // - s : size of stack header (number of header, (default = 0x0007)
400 // - v : header version (default = 0xa)
401 // - m : link mask (default = 0xfff)
403 // 4) stack header : vvvvvvvv vvvvvvvv bbbbbbbb rrrr rrr c
404 // - v : hardware design revision (default = 0x0404)
405 // - b : physical board ID (default = 0x5b)
406 // - r : reserved for future use (default = 0000 000)
407 // - c : clean checkout flag (default = 1)
409 // and 6 dummy words(0x00000000)
412 //buf[nw++] = 0x0001a03f; // SM index words
413 fSMindexPos = nw; // memorize position of the SM index word for re-allocating stack mask
414 buf[nw++] = 0x0001a020; // SM index words
415 buf[nw++] = 0x10404071; // SM header
417 fStackindexPos = nw; // memorize position of the stack index word for future adding
419 for (Int_t istack=0; istack<5; istack++){
420 buf[nw++] = 0x0007afff; // stack index words
421 buf[nw++] = 0x04045b01; // stack header
422 for (Int_t i=0;i<6;i++) buf[nw++] = 0x00000000; // 6 dummy words
423 } // loop over 5 stacks
427 //_____________________________________________________________________________
428 void AliTRDrawData::AssignStackMask(UInt_t *buf, Int_t nStack){
430 // This function re-assign stack mask active(from 0 to 1) in the SM index word
432 buf[fSMindexPos] = buf[fSMindexPos] | ( 1 << nStack );
435 //_____________________________________________________________________________
436 Int_t AliTRDrawData::AddStackIndexWords(UInt_t *buf, Int_t /*nStack*/, Int_t nMax){
438 // This function add stack index words and stack header when there is data for the stack
440 // 1) stack index words : ssssssss ssssssss vvvv mmmm mmmmmmmm
441 // - s : size of stack header (number of header, (default = 0x0007)
442 // - v : header version (default = 0xa)
443 // - m : link mask (default = 0xfff)
444 // - m : link mask (starting value = 0x000)
446 // 2) stack header : vvvvvvvv vvvvvvvv bbbbbbbb rrrr rrr c
447 // - v : hardware design revision (default = 0x0404)
448 // - b : physical board ID (default = 0x5b)
449 // - r : reserved for future use (default = 0000 000)
450 // - c : clean checkout flag (default = 1)
452 // and 6 dummy words(0x00000000)
455 Int_t nAddedWords = 0; // Number of added words
456 if ( ShiftWords(buf, fStackindexPos, 8, nMax)== kFALSE ){
457 AliError("Adding stack header failed.");
461 buf[fStackindexPos++] = 0x0007a000; // stack index words
462 buf[fStackindexPos++] = 0x04045b01; // stack header
463 for (Int_t i=0;i<6;i++) buf[fStackindexPos++] = 0x00000000; // 6 dummy words
469 //_____________________________________________________________________________
470 void AliTRDrawData::AssignLinkMask(UInt_t *buf, Int_t nLayer){
472 // This function re-assign link mask active(from 0 to 1) in the stack index word
474 buf[fStackindexPos-8] = buf[fStackindexPos-8] | ( 3 << (2*nLayer) ); // 3 = 0011
477 //_____________________________________________________________________________
478 Bool_t AliTRDrawData::ShiftWords(UInt_t *buf, Int_t nStart, Int_t nWords, Int_t nMax){
480 // This function shifts n words
482 //if ( nStart+nWords > sizeof(buf)/sizeof(UInt_t) ){
483 // AliError("Words shift failed. No more buffer space.");
487 for ( Int_t iw=nMax; iw>nStart-1; iw--){
488 buf[iw+nWords] = buf[iw];
493 //_____________________________________________________________________________
494 Int_t AliTRDrawData::ProduceHcData(AliTRDarrayADC *digits, Int_t side, Int_t det, UInt_t *buf, Int_t maxSize, Bool_t /*newEvent = kFALSE*/, Bool_t /*newSM = kFALSE*/){
496 // This function can be used for both ZS and NZS data
499 Int_t nw = 0; // Number of written words
500 Int_t of = 0; // Number of overflowed words
501 Int_t *tempnw = 0x0; // Number of written words for temp. buffer
502 Int_t *tempof = 0x0; // Number of overflowed words for temp. buffer
503 Int_t layer = fGeo->GetLayer( det ); // Layer
504 Int_t stack = fGeo->GetStack( det ); // Stack
505 Int_t sect = fGeo->GetSector( det ); // Sector (=iDDL)
506 const Int_t kCtype = fGeo->GetStack(det) == 2 ? 0 : 1; // Chamber type (0:C0, 1:C1)
508 Bool_t tracklet_on = fFee->GetTracklet(); // tracklet simulation active?
510 AliDebug(1,Form("Producing raw data for sect=%d layer=%d stack=%d side=%d",sect,layer,stack,side));
512 AliTRDmcmSim* mcm = new AliTRDmcmSim();
514 UInt_t *tempBuffer = buf; // tempBuffer used to write ADC data
515 // different in case of tracklet writing
518 tempBuffer = new UInt_t[maxSize];
519 tempnw = new Int_t(0);
520 tempof = new Int_t(0);
527 WriteIntermediateWordsV2(tempBuffer,*tempnw,*tempof,maxSize,det,side); //??? no tracklet or NZS
529 // scanning direction such, that tracklet-words are sorted in ascending z and then in ascending y order
530 // ROB numbering on chamber and MCM numbering on ROB increase with decreasing z and increasing y
531 for (Int_t iRobRow = 0; iRobRow <= (kCtype + 3)-1; iRobRow++ ) { // ROB number should be increasing
532 Int_t iRob = iRobRow * 2 + side;
534 for (Int_t iMcmRB = 0; iMcmRB < fGeo->MCMmax(); iMcmRB++ ) {
535 Int_t iMcm = 16 - 4*(iMcmRB/4 + 1) + (iMcmRB%4);
537 mcm->Init(det, iRob, iMcm);
538 mcm->SetData(digits); // no filtering done here (already done in digitizer)
541 Int_t tempNw = mcm->ProduceTrackletStream(&buf[nw], maxSize - nw);
545 AliError(Form("Buffer overflow detected. Please increase the buffer size and recompile."));
550 mcm->ZSMapping(); // Calculate zero suppression mapping
551 // at the moment it has to be rerun here
552 // Write MCM data to temp. buffer
553 Int_t tempNw = mcm->ProduceRawStream( &tempBuffer[*tempnw], maxSize - *tempnw, fEventCounter );
556 *tempnw += maxSize - nw;
557 AliError(Form("Buffer overflow detected. Please increase the buffer size and recompile."));
567 // in case of tracklet writing copy temp data to final buffer
569 if (nw + *tempnw < maxSize) {
570 memcpy(&buf[nw], tempBuffer, *tempnw * sizeof(UInt_t));
574 AliError("Buffer overflow detected");
578 // Write end of raw data marker
579 if (nw+3 < maxSize) {
580 buf[nw++] = 0x00000000; // fFee->GetRawDataEndmarker();
581 buf[nw++] = 0x00000000; // fFee->GetRawDataEndmarker();
582 buf[nw++] = 0x00000000; // fFee->GetRawDataEndmarker();
583 buf[nw++] = 0x00000000; // fFee->GetRawDataEndmarker();
589 delete [] tempBuffer;
595 AliError("Buffer overflow. Data is truncated. Please increase buffer size and recompile.");
601 //_____________________________________________________________________________
602 Int_t AliTRDrawData::ProduceHcDataV1andV2(AliTRDarrayADC *digits, Int_t side
603 , Int_t det, UInt_t *buf, Int_t maxSize)
606 // This function simulates: 1) SM-I commissiong data Oct. 06 (Raw Version == 1).
607 // 2) Full Raw Production Version (Raw Version == 2)
609 // Produce half chamber data (= an ORI data) for the given chamber (det) and side (side)
612 // side=0 means A side with ROB positions 0, 2, 4, 6.
613 // side=1 means B side with ROB positions 1, 3, 5, 7.
615 // Chamber type (C0 orC1) is determined by "det" automatically.
616 // Appropriate size of buffer (*buf) must be prepared prior to calling this function.
617 // Pointer to the buffer and its size must be given to "buf" and "maxSize".
618 // Return value is the number of valid data filled in the buffer in unit of 32 bits
620 // If buffer size if too small, the data is truncated with the buffer size however
621 // the function will finish without crash (this behaviour is similar to the MCM).
624 Int_t nw = 0; // Number of written words
625 Int_t of = 0; // Number of overflowed words
626 Int_t layer = fGeo->GetLayer( det ); // Layer
627 Int_t stack = fGeo->GetStack( det ); // Stack
628 Int_t sect = fGeo->GetSector( det ); // Sector (=iDDL)
629 Int_t nRow = fGeo->GetRowMax( layer, stack, sect );
630 Int_t nCol = fGeo->GetColMax( layer );
631 const Int_t kNTBin = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
632 Int_t kCtype = 0; // Chamber type (0:C0, 1:C1)
633 Int_t iEv = 0xA; // Event ID. Now fixed to 10, how do I get event id?
634 UInt_t x = 0; // General used number
635 Int_t rv = fFee->GetRAWversion();
637 // Check the nCol and nRow.
639 (nRow == 16 || nRow == 12)) {
640 kCtype = (nRow-12) / 4;
643 AliError(Form("This type of chamber is not supported (nRow=%d, nCol=%d)."
648 AliDebug(1,Form("Producing raw data for sect=%d layer=%d stack=%d side=%d"
649 ,sect,layer,stack,side));
651 // Tracklet should be processed here but not implemented yet
653 // Write end of tracklet marker
655 buf[nw++] = kEndoftrackletmarker;
661 // Half Chamber header
663 // Now it is the same version as used in SM-I commissioning.
664 Int_t dcs = det+100; // DCS Serial (in simulation, it is meaningless
665 x = (dcs<<20) | (sect<<15) | (layer<<12) | (stack<<9) | (side<<8) | 1;
673 else if ( rv == 2 ) {
674 // h[0] (there are 3 HC header)
675 Int_t minorv = 0; // The minor version number
676 Int_t add = 2; // The number of additional header words to follow
677 x = (1<<31) | (rv<<24) | (minorv<<17) | (add<<14) | (sect<<9) | (layer<<6) | (stack<<3) | (side<<2) | 1;
685 Int_t bcCtr = 99; // bunch crossing counter. Here it is set to 99 always for no reason
686 Int_t ptCtr = 15; // pretrigger counter. Here it is set to 15 always for no reason
687 Int_t ptPhase = 11; // pretrigger phase. Here it is set to 11 always for no reason
688 x = (bcCtr<<16) | (ptCtr<<12) | (ptPhase<<8) | ((kNTBin-1)<<2) | 1;
696 Int_t pedSetup = 1; // Pedestal filter setup (0:1). Here it is always 1 for no reason
697 Int_t gainSetup = 1; // Gain filter setup (0:1). Here it is always 1 for no reason
698 Int_t tailSetup = 1; // Tail filter setup (0:1). Here it is always 1 for no reason
699 Int_t xtSetup = 0; // Cross talk filter setup (0:1). Here it is always 0 for no reason
700 Int_t nonlinSetup = 0; // Nonlinearity filter setup (0:1). Here it is always 0 for no reason
701 Int_t bypassSetup = 0; // Filter bypass (for raw data) setup (0:1). Here it is always 0 for no reason
702 Int_t commonAdditive = 10; // Digital filter common additive (0:63). Here it is always 10 for no reason
703 x = (pedSetup<<31) | (gainSetup<<30) | (tailSetup<<29) | (xtSetup<<28) | (nonlinSetup<<27)
704 | (bypassSetup<<26) | (commonAdditive<<20) | 1;
713 // Scan for ROB and MCM
714 for (Int_t iRobRow = 0; iRobRow < (kCtype + 3); iRobRow++ ) {
715 Int_t iRob = iRobRow * 2 + side;
716 for (Int_t iMcm = 0; iMcm < fGeo->MCMmax(); iMcm++ ) {
717 Int_t padrow = iRobRow * 4 + iMcm / 4;
720 x = ((iRob * fGeo->MCMmax() + iMcm) << 24) | ((iEv % 0x100000) << 4) | 0xC;
729 for (Int_t iAdc = 0; iAdc < 21; iAdc++ ) {
730 Int_t padcol = fFee->GetPadColFromADC(iRob, iMcm, iAdc);
731 UInt_t aa = !(iAdc & 1) + 2; // 3 for the even ADC channel , 2 for the odd ADC channel
732 UInt_t *a = new UInt_t[kNTBin+2];
733 // 3 timebins are packed into one 32 bits word
734 for (Int_t iT = 0; iT < kNTBin; iT+=3) {
735 if ((padcol >= 0) && (padcol < nCol)) {
736 a[iT ] = ((iT ) < kNTBin ) ? digits->GetData(padrow,padcol,iT ) : 0;
737 a[iT+1] = ((iT + 1) < kNTBin ) ? digits->GetData(padrow,padcol,iT + 1) : 0;
738 a[iT+2] = ((iT + 2) < kNTBin ) ? digits->GetData(padrow,padcol,iT + 2) : 0;
741 a[iT] = a[iT+1] = a[iT+2] = 0; // This happenes at the edge of chamber (should be pedestal! How?)
743 x = (a[iT+2] << 22) | (a[iT+1] << 12) | (a[iT] << 2) | aa;
754 for (Int_t iT = 0; iT < kNTBin; iT++) {
755 avg += (Float_t) (a[iT]);
757 avg /= (Float_t) kNTBin;
758 for (Int_t iT = 0; iT < kNTBin; iT++) {
759 rms += ((Float_t) (a[iT]) - avg) * ((Float_t) (a[iT]) - avg);
761 rms = TMath::Sqrt(rms / (Float_t) kNTBin);
763 AliDebug(2,Form("Large RMS (>1.7) (ROB,MCM,ADC)=(%02d,%02d,%02d), avg=%03.1f, rms=%03.1f"
764 ,iRob,iMcm,iAdc,avg,rms));
771 // Write end of raw data marker
773 buf[nw++] = kEndofrawdatamarker;
779 AliWarning("Buffer overflow. Data is truncated. Please increase buffer size and recompile.");
786 //_____________________________________________________________________________
788 //Int_t AliTRDrawData::ProduceHcDataV3(AliTRDarrayADC *digits, Int_t side , Int_t det, UInt_t *buf, Int_t maxSize)
789 Int_t AliTRDrawData::ProduceHcDataV3(AliTRDarrayADC *digits, Int_t side , Int_t det, UInt_t *buf, Int_t maxSize, Bool_t newEvent = kFALSE)
792 // This function simulates: Raw Version == 3 (Zero Suppression Prototype)
795 Int_t nw = 0; // Number of written words
796 Int_t of = 0; // Number of overflowed words
797 Int_t layer = fGeo->GetLayer( det ); // Layer
798 Int_t stack = fGeo->GetStack( det ); // Stack
799 Int_t sect = fGeo->GetSector( det ); // Sector (=iDDL)
800 Int_t nRow = fGeo->GetRowMax( layer, stack, sect );
801 Int_t nCol = fGeo->GetColMax( layer );
802 const Int_t kNTBin = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
803 Int_t kCtype = 0; // Chamber type (0:C0, 1:C1)
804 //Int_t iEv = 0xA; // Event ID. Now fixed to 10, how do I get event id?
808 Bool_t tracklet_on = fFee->GetTracklet(); // **new**
810 // Check the nCol and nRow.
812 (nRow == 16 || nRow == 12)) {
813 kCtype = (nRow-12) / 4;
816 AliError(Form("This type of chamber is not supported (nRow=%d, nCol=%d)."
821 AliDebug(1,Form("Producing raw data for sect=%d layer=%d stack=%d side=%d"
822 ,sect,layer,stack,side));
824 AliTRDmcmSim** mcm = new AliTRDmcmSim*[(kCtype + 3)*(fGeo->MCMmax())];
826 // in case no tracklet-words are processed: write the tracklet-endmarker as well as all additional words immediately and write
827 // raw-data in one go; if tracklet-processing is enabled, first all tracklet-words of a half-chamber have to be processed before the
828 // additional words (tracklet-endmarker,headers,...)are written. Raw-data is written in a second loop;
831 WriteIntermediateWords(buf,nw,of,maxSize,det,side);
834 // Scan for ROB and MCM
835 // scanning direction such, that tracklet-words are sorted in ascending z and then in ascending y order
836 // ROB numbering on chamber and MCM numbering on ROB increase with decreasing z and increasing y
837 for (Int_t iRobRow = (kCtype + 3)-1; iRobRow >= 0; iRobRow-- ) {
838 Int_t iRob = iRobRow * 2 + side;
840 for (Int_t iMcmRB = 0; iMcmRB < fGeo->MCMmax(); iMcmRB++ ) {
841 Int_t iMcm = 16 - 4*(iMcmRB/4 + 1) + (iMcmRB%4);
842 Int_t entry = iRobRow*(fGeo->MCMmax()) + iMcm;
844 mcm[entry] = new AliTRDmcmSim();
845 mcm[entry]->Init( det, iRob, iMcm , newEvent);
846 //mcm[entry]->Init( det, iRob, iMcm);
847 if (newEvent == kTRUE) newEvent = kFALSE; // only one mcm is concerned with new event
848 Int_t padrow = mcm[entry]->GetRow();
850 // Copy ADC data to MCM simulator
851 for (Int_t iAdc = 0; iAdc < 21; iAdc++ ) {
852 Int_t padcol = mcm[entry]->GetCol( iAdc );
853 if ((padcol >= 0) && (padcol < nCol)) {
854 for (Int_t iT = 0; iT < kNTBin; iT++) {
855 mcm[entry]->SetData( iAdc, iT, digits->GetData( padrow, padcol, iT) );
858 else { // this means it is out of chamber, and masked ADC
859 mcm[entry]->SetDataPedestal( iAdc );
863 // Simulate process in MCM
864 mcm[entry]->Filter(); // Apply filter
865 mcm[entry]->ZSMapping(); // Calculate zero suppression mapping
866 //jkl mcm[entry]->CopyArrays();
867 //jkl mcm[entry]->GeneratefZSM1Dim();
868 //jkl mcm[entry]->RestoreZeros();
871 mcm[entry]->Tracklet();
872 Int_t tempNw = mcm[entry]->ProduceTrackletStream( &buf[nw], maxSize - nw );
877 AliError(Form("Buffer overflow detected. Please increase the buffer size and recompile."));
882 // no tracklets: write raw-data already in this loop
884 // Write MCM data to buffer
885 Int_t tempNw = mcm[entry]->ProduceRawStream( &buf[nw], maxSize - nw );
889 AliError(Form("Buffer overflow detected. Please increase the buffer size and recompile."));
900 //mcm->DumpData( "trdmcmdata.txt", "RFZS" ); // debugging purpose
904 // if tracklets are switched on, raw-data can be written only after all tracklets
906 WriteIntermediateWords(buf,nw,of,maxSize,det,side);
909 // Scan for ROB and MCM
910 for (Int_t iRobRow = (kCtype + 3)-1; iRobRow >= 0; iRobRow-- ) {
911 //Int_t iRob = iRobRow * 2 + side;
913 for (Int_t iMcmRB = 0; iMcmRB < fGeo->MCMmax(); iMcmRB++ ) {
914 Int_t iMcm = 16 - 4*(iMcmRB/4 + 1) + (iMcmRB%4);
916 Int_t entry = iRobRow*(fGeo->MCMmax()) + iMcm;
918 // Write MCM data to buffer
919 Int_t tempNw = mcm[entry]->ProduceRawStream( &buf[nw], maxSize - nw );
923 AliError(Form("Buffer overflow detected. Please increase the buffer size and recompile."));
936 // Write end of raw data marker
938 buf[nw++] = kEndofrawdatamarker;
944 AliError("Buffer overflow. Data is truncated. Please increase buffer size and recompile.");
952 //_____________________________________________________________________________
953 AliTRDdigitsManager *AliTRDrawData::Raw2Digits(AliRawReader *rawReader)
956 // Vx of the raw data reading
959 rawReader->Select("TRD"); //[mj]
961 AliTRDarrayADC *digits = 0;
962 AliTRDarrayDictionary *track0 = 0;
963 AliTRDarrayDictionary *track1 = 0;
964 AliTRDarrayDictionary *track2 = 0;
966 //AliTRDSignalIndex *indexes = 0;
967 // Create the digits manager
968 AliTRDdigitsManager* digitsManager = new AliTRDdigitsManager();
969 digitsManager->CreateArrays();
971 if (!fTrackletContainer) {
972 //if (!fTrackletContainer && ( fReconstructor->IsWritingTracklets() || fReconstructor->IsProcessingTracklets() )) {
973 // maximum tracklets for one HC
974 const Int_t kTrackletChmb=256;
975 fTrackletContainer = new UInt_t *[2];
976 fTrackletContainer[0] = new UInt_t[kTrackletChmb];
977 fTrackletContainer[1] = new UInt_t[kTrackletChmb];
980 AliTRDrawStreamBase *pinput = AliTRDrawStreamBase::GetRawStream(rawReader);
981 AliTRDrawStreamBase &input = *pinput;
982 input.SetRawVersion( fFee->GetRAWversion() ); //<= ADDED by MinJung
984 AliInfo(Form("Stream version: %s", input.IsA()->GetName()));
986 // Loop through the digits
991 //det = input.NextChamber(digitsManager);
992 det = input.NextChamber(digitsManager,fTrackletContainer);
994 //if (!fReconstructor->IsWritingTracklets()) continue;
995 if (*(fTrackletContainer[0]) > 0 || *(fTrackletContainer[1]) > 0) WriteTracklets(det);
1000 digits = (AliTRDarrayADC *) digitsManager->GetDigits(det);
1001 track0 = (AliTRDarrayDictionary *) digitsManager->GetDictionary(det,0);
1002 track1 = (AliTRDarrayDictionary *) digitsManager->GetDictionary(det,1);
1003 track2 = (AliTRDarrayDictionary *) digitsManager->GetDictionary(det,2);
1005 if (digits) digits->Compress();
1006 if (track0) track0->Compress();
1007 if (track1) track1->Compress();
1008 if (track2) track2->Compress();
1012 if (fTrackletContainer){
1013 delete [] fTrackletContainer[0];
1014 delete [] fTrackletContainer[1];
1015 delete [] fTrackletContainer;
1016 fTrackletContainer = NULL;
1022 return digitsManager;
1025 //_____________________________________________________________________________
1026 void AliTRDrawData::WriteIntermediateWords(UInt_t* buf, Int_t& nw, Int_t& of, const Int_t& maxSize, const Int_t& det, const Int_t& side) {
1028 Int_t layer = fGeo->GetLayer( det ); // Layer
1029 Int_t stack = fGeo->GetStack( det ); // Stack
1030 Int_t sect = fGeo->GetSector( det ); // Sector (=iDDL)
1031 Int_t rv = fFee->GetRAWversion();
1032 const Int_t kNTBin = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
1035 // Write end of tracklet marker
1037 buf[nw++] = kEndoftrackletmarker;
1043 // Half Chamber header
1044 // h[0] (there are 3 HC header)
1045 Int_t minorv = 0; // The minor version number
1046 Int_t add = 2; // The number of additional header words to follow
1047 x = (1<<31) | (rv<<24) | (minorv<<17) | (add<<14) | (sect<<9) | (layer<<6) | (stack<<3) | (side<<2) | 1;
1055 Int_t bcCtr = 99; // bunch crossing counter. Here it is set to 99 always for no reason
1056 Int_t ptCtr = 15; // pretrigger counter. Here it is set to 15 always for no reason
1057 Int_t ptPhase = 11; // pretrigger phase. Here it is set to 11 always for no reason
1058 x = (bcCtr<<16) | (ptCtr<<12) | (ptPhase<<8) | ((kNTBin-1)<<2) | 1;
1066 Int_t pedSetup = 1; // Pedestal filter setup (0:1). Here it is always 1 for no reason
1067 Int_t gainSetup = 1; // Gain filter setup (0:1). Here it is always 1 for no reason
1068 Int_t tailSetup = 1; // Tail filter setup (0:1). Here it is always 1 for no reason
1069 Int_t xtSetup = 0; // Cross talk filter setup (0:1). Here it is always 0 for no reason
1070 Int_t nonlinSetup = 0; // Nonlinearity filter setup (0:1). Here it is always 0 for no reason
1071 Int_t bypassSetup = 0; // Filter bypass (for raw data) setup (0:1). Here it is always 0 for no reason
1072 Int_t commonAdditive = 10; // Digital filter common additive (0:63). Here it is always 10 for no reason
1073 x = (pedSetup<<31) | (gainSetup<<30) | (tailSetup<<29) | (xtSetup<<28) | (nonlinSetup<<27)
1074 | (bypassSetup<<26) | (commonAdditive<<20) | 1;
1083 //_____________________________________________________________________________
1084 void AliTRDrawData::WriteIntermediateWordsV2(UInt_t* buf, Int_t& nw, Int_t& of, const Int_t& maxSize, const Int_t& det, const Int_t& side) {
1086 Int_t layer = fGeo->GetLayer( det ); // Layer
1087 Int_t stack = fGeo->GetStack( det ); // Stack
1088 Int_t sect = fGeo->GetSector( det ); // Sector (=iDDL)
1089 Int_t rv = fFee->GetRAWversion();
1090 const Int_t kNTBin = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
1091 Bool_t tracklet_on = fFee->GetTracklet();
1094 // Write end of tracklet marker
1096 buf[nw++] = fgkEndOfTrackletMarker;
1097 buf[nw++] = fgkEndOfTrackletMarker; // the number of tracklet end marker should be more than 2
1104 // Half Chamber header
1105 // h[0] (there are 2 HC headers) xmmm mmmm nnnn nnnq qqss sssp ppcc ci01
1106 // , where x : Raw version speacial number (=1)
1107 // m : Raw version major number (test pattern, ZS, disable tracklet, 0, options)
1108 // n : Raw version minor number
1109 // q : number of addtional header words (default = 1)
1110 // s : SM sector number (ALICE numbering)
1111 // p : plane(layer) number
1112 // c : chamber(stack) number
1113 // i : side number (0:A, 1:B)
1114 Int_t majorv = 0; // The major version number
1115 Int_t minorv = 0; // The minor version number
1116 Int_t add = 1; // The number of additional header words to follow : now 1, previous 2
1117 Int_t TP = 0; // test pattern (default=0)
1118 Int_t ZS = (rv==3) ? 1 : 0; // zero suppression
1119 Int_t DT = (tracklet_on) ? 0 : 1; // disable tracklet
1121 majorv = (TP<<6) | (ZS<<5) | (DT<<4) | 1; // major version
1123 x = (1<<31) | (majorv<<24) | (minorv<<17) | (add<<14) | (sect<<9) | (layer<<6) | (stack<<3) | (side<<2) | 1;
1124 if (nw < maxSize) buf[nw++] = x; else of++;
1126 // h[1] tttt ttbb bbbb bbbb bbbb bbpp pphh hh01
1127 // , where t : number of time bins
1128 // b : bunch crossing number
1129 // p : pretrigger counter
1130 // h : pretrigger phase
1131 Int_t bcCtr = 99; // bunch crossing counter. Here it is set to 99 always for no reason
1132 Int_t ptCtr = 15; // pretrigger counter. Here it is set to 15 always for no reason
1133 Int_t ptPhase = 11; // pretrigger phase. Here it is set to 11 always for no reason
1134 //x = (bcCtr<<16) | (ptCtr<<12) | (ptPhase<<8) | ((kNTBin-1)<<2) | 1; // old format
1135 x = ((kNTBin)<<26) | (bcCtr<<10) | (ptCtr<<6) | (ptPhase<<2) | 1;
1136 if (nw < maxSize) buf[nw++] = x; else of++;
1140 //_____________________________________________________________________________
1141 AliTRDdigitsManager *AliTRDrawData::Raw2DigitsOLD(AliRawReader *rawReader)
1144 // Vx of the raw data reading
1147 AliTRDarrayADC *digits = 0;
1148 AliTRDarrayDictionary *track0 = 0;
1149 AliTRDarrayDictionary *track1 = 0;
1150 AliTRDarrayDictionary *track2 = 0;
1152 AliTRDSignalIndex *indexes = 0;
1153 // Create the digits manager
1154 AliTRDdigitsManager* digitsManager = new AliTRDdigitsManager();
1155 digitsManager->CreateArrays();
1157 AliTRDrawOldStream input(rawReader);
1158 input.SetRawVersion( fFee->GetRAWversion() );
1161 AliInfo(Form("Stream version: %s", input.IsA()->GetName()));
1163 // Loop through the digits
1167 while (input.Next()) {
1169 det = input.GetDet();
1171 if (det != lastdet) { // If new detector found
1175 if (digits) digits->Compress();
1176 if (track0) track0->Compress();
1177 if (track1) track1->Compress();
1178 if (track2) track2->Compress();
1180 // Add a container for the digits of this detector
1181 digits = (AliTRDarrayADC *) digitsManager->GetDigits(det);
1182 track0 = (AliTRDarrayDictionary *) digitsManager->GetDictionary(det,0);
1183 track1 = (AliTRDarrayDictionary *) digitsManager->GetDictionary(det,1);
1184 track2 = (AliTRDarrayDictionary *) digitsManager->GetDictionary(det,2);
1186 // Allocate memory space for the digits buffer
1187 if (digits->GetNtime() == 0)
1189 digits->Allocate(input.GetMaxRow(),input.GetMaxCol(), input.GetNumberOfTimeBins());
1190 track0->Allocate(input.GetMaxRow(),input.GetMaxCol(), input.GetNumberOfTimeBins());
1191 track1->Allocate(input.GetMaxRow(),input.GetMaxCol(), input.GetNumberOfTimeBins());
1192 track2->Allocate(input.GetMaxRow(),input.GetMaxCol(), input.GetNumberOfTimeBins());
1195 indexes = digitsManager->GetIndexes(det);
1196 indexes->SetSM(input.GetSM());
1197 indexes->SetStack(input.GetStack());
1198 indexes->SetLayer(input.GetLayer());
1199 indexes->SetDetNumber(det);
1200 if (indexes->IsAllocated() == kFALSE)
1201 indexes->Allocate(input.GetMaxRow(), input.GetMaxCol(), input.GetNumberOfTimeBins());
1204 // 3 timebin data are stored per word
1205 for (it = 0; it < 3; it++)
1207 if ( input.GetTimeBin() + it < input.GetNumberOfTimeBins() )
1209 if (input.GetSignals()[it] > 0)
1211 digits->SetData(input.GetRow(), input.GetCol(),input.GetTimeBin() + it, input.GetSignals()[it]);
1213 indexes->AddIndexRC(input.GetRow(), input.GetCol());
1214 track0->SetData(input.GetRow(), input.GetCol(), input.GetTimeBin() + it, 0);
1215 track1->SetData(input.GetRow(), input.GetCol(), input.GetTimeBin() + it, 0);
1216 track2->SetData(input.GetRow(), input.GetCol(), input.GetTimeBin() + it, 0);
1222 if (digits) digits->Compress();
1223 if (track0) track0->Compress();
1224 if (track1) track1->Compress();
1225 if (track2) track2->Compress();
1227 return digitsManager;
1231 //_____________________________________________________________________________
1232 Bool_t AliTRDrawData::WriteTracklets(Int_t det)
1235 // Write the raw data tracklets into seperate file
1238 UInt_t **leaves = new UInt_t *[2];
1239 for (Int_t i=0; i<2 ;i++){
1240 leaves[i] = new UInt_t[258];
1241 leaves[i][0] = det; // det
1242 leaves[i][1] = i; // side
1243 memcpy(leaves[i]+2, fTrackletContainer[i], sizeof(UInt_t) * 256);
1246 if (!fTrackletTree){
1247 AliDataLoader *dl = fRunLoader->GetLoader("TRDLoader")->GetDataLoader("tracklets");
1249 fTrackletTree = dl->Tree();
1252 TBranch *trkbranch = fTrackletTree->GetBranch("trkbranch");
1254 trkbranch = fTrackletTree->Branch("trkbranch",leaves[0],"det/i:side/i:tracklets[256]/i");
1257 for (Int_t i=0; i<2; i++){
1258 if (leaves[i][2]>0) {
1259 trkbranch->SetAddress(leaves[i]);
1260 fTrackletTree->Fill();
1264 AliDataLoader *dl = fRunLoader->GetLoader("TRDLoader")->GetDataLoader("tracklets");
1265 dl->WriteData("OVERWRITE");
1273 //_____________________________________________________________________________
1274 Bool_t AliTRDrawData::OpenOutput()
1277 // Connect the output tree
1282 //if (fReconstructor->IsWritingTracklets()){
1283 TString evfoldname = AliConfig::GetDefaultEventFolderName();
1284 fRunLoader = AliRunLoader::GetRunLoader(evfoldname);
1287 fRunLoader = AliRunLoader::Open("galice.root");
1290 AliError(Form("Can not open session for file galice.root."));
1294 UInt_t **leaves = new UInt_t *[2];
1295 AliDataLoader *dl = fRunLoader->GetLoader("TRDLoader")->GetDataLoader("tracklets");
1297 AliError("Could not get the tracklets data loader!");
1298 dl = new AliDataLoader("TRD.Tracklets.root","tracklets", "tracklets");
1299 fRunLoader->GetLoader("TRDLoader")->AddDataLoader(dl);
1301 fTrackletTree = dl->Tree();
1305 fTrackletTree = dl->Tree();
1307 TBranch *trkbranch = fTrackletTree->GetBranch("trkbranch");
1309 fTrackletTree->Branch("trkbranch",leaves[0],"det/i:side/i:tracklets[256]/i");