CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
DTTFFEDReader Class Reference

#include <DTTFFEDReader.h>

Inheritance diagram for DTTFFEDReader:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

Public Member Functions

 DTTFFEDReader (const edm::ParameterSet &pset)
 Constructor. More...
 
bool fillRawData (edm::Event &e, L1MuDTChambPhContainer::Phi_Container &phi_data, L1MuDTChambThContainer::The_Container &the_data, L1MuDTTrackContainer::TrackContainer &tra_data)
 Generate and fill FED raw data for a full event. More...
 
void produce (edm::Event &e, const edm::EventSetup &c)
 Produce digis out of raw data. More...
 
virtual ~DTTFFEDReader ()
 Destructor. More...
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

Private Member Functions

void analyse (edm::Event &e)
 
int bxNr (int channel)
 
void calcCRC (int myD1, int myD2, int &myC)
 
int channel (int wheel, int sector, int bx)
 
void clear ()
 
edm::InputTag getDTTFInputTag ()
 
const
L1MuDTTrackContainer::TrackContainer
k_data ()
 
void match ()
 
const
L1MuDTChambPhContainer::Phi_Container
p_data ()
 
void process (edm::Event &e)
 
int sector (int channel)
 
const
L1MuDTChambThContainer::The_Container
t_data ()
 
int wheel (int channel)
 

Private Attributes

edm::InputTag DTTFInputTag
 
L1MuDTTrackContainer::TrackContainer dtTracks
 
unsigned int efTrack [3][12][6][2]
 
unsigned int etTrack [3][12][6][2]
 
L1MuDTChambPhContainer::Phi_Container phiSegments
 
L1MuDTChambThContainer::The_Container theSegments
 
bool verbose_
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 

Detailed Description

L1 DT Track Finder Raw-to-Digi

Date:
2009/11/18 13:27:12
Revision:
1.7

J. Troconiz UAM Madrid E. Delmeire UAM Madrid

Definition at line 29 of file DTTFFEDReader.h.

Constructor & Destructor Documentation

DTTFFEDReader::DTTFFEDReader ( const edm::ParameterSet pset)

Constructor.

Definition at line 30 of file DTTFFEDReader.cc.

References edm::ParameterSet::getParameter(), and edm::ParameterSet::getUntrackedParameter().

30  {
31 
32  produces<L1MuDTChambPhContainer>();
33  produces<L1MuDTChambThContainer>();
34  produces<L1MuDTTrackContainer>("DATA");
35 
36  DTTFInputTag = pset.getParameter<edm::InputTag>("DTTF_FED_Source");
37 
38  verbose_ = pset.getUntrackedParameter<bool>("verbose",false);
39 
40 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::InputTag DTTFInputTag
Definition: DTTFFEDReader.h:50
DTTFFEDReader::~DTTFFEDReader ( )
virtual

Destructor.

Definition at line 42 of file DTTFFEDReader.cc.

42 {}

Member Function Documentation

void DTTFFEDReader::analyse ( edm::Event e)
private

Definition at line 83 of file DTTFFEDReader.cc.

References match(), and align_tpl::process.

83  {
84  clear();
85  process(e);
86  match();
87  return;
88 }
void process(edm::Event &e)
int DTTFFEDReader::bxNr ( int  channel)
private

Definition at line 490 of file DTTFFEDReader.cc.

490  {
491 
492  int myChannel = channel;
493 
494  if (myChannel > 127) myChannel -= 2;
495 
496  if (myChannel < 0 || myChannel > 251 ){ return -999; }
497 
498  int myBx = 1-(myChannel%3);
499 
500  return myBx;
501 }
int channel(int wheel, int sector, int bx)
void DTTFFEDReader::calcCRC ( int  myD1,
int  myD2,
int &  myC 
)
private

Definition at line 527 of file DTTFFEDReader.cc.

References funct::C, and i.

527  {
528 
529  int myCRC[16],D[64],C[16];
530 
531  for( int i=0; i < 32; i++ ){ D[i]=(myD2>>i)&0x1; }
532  for( int i=0; i < 32; i++ ){ D[i+32]=(myD1>>i)&0x1; }
533  for( int i=0; i < 16; i++ ){ C[i]=(myC>>i)&0x1; }
534 
535  myCRC[0] = ( D[63] + D[62] + D[61] + D[60] + D[55] + D[54] +
536  D[53] + D[52] + D[51] + D[50] + D[49] + D[48] +
537  D[47] + D[46] + D[45] + D[43] + D[41] + D[40] +
538  D[39] + D[38] + D[37] + D[36] + D[35] + D[34] +
539  D[33] + D[32] + D[31] + D[30] + D[27] + D[26] +
540  D[25] + D[24] + D[23] + D[22] + D[21] + D[20] +
541  D[19] + D[18] + D[17] + D[16] + D[15] + D[13] +
542  D[12] + D[11] + D[10] + D[9] + D[8] + D[7] +
543  D[6] + D[5] + D[4] + D[3] + D[2] + D[1] +
544  D[0] + C[0] + C[1] + C[2] + C[3] + C[4] +
545  C[5] + C[6] + C[7] + C[12] + C[13] + C[14] +
546  C[15] )%2;
547 
548  myCRC[1] = ( D[63] + D[62] + D[61] + D[56] + D[55] + D[54] +
549  D[53] + D[52] + D[51] + D[50] + D[49] + D[48] +
550  D[47] + D[46] + D[44] + D[42] + D[41] + D[40] +
551  D[39] + D[38] + D[37] + D[36] + D[35] + D[34] +
552  D[33] + D[32] + D[31] + D[28] + D[27] + D[26] +
553  D[25] + D[24] + D[23] + D[22] + D[21] + D[20] +
554  D[19] + D[18] + D[17] + D[16] + D[14] + D[13] +
555  D[12] + D[11] + D[10] + D[9] + D[8] + D[7] +
556  D[6] + D[5] + D[4] + D[3] + D[2] + D[1] +
557  C[0] + C[1] + C[2] + C[3] + C[4] + C[5] +
558  C[6] + C[7] + C[8] + C[13] + C[14] + C[15] )%2;
559 
560  myCRC[2] = ( D[61] + D[60] + D[57] + D[56] + D[46] + D[42] +
561  D[31] + D[30] + D[29] + D[28] + D[16] + D[14] +
562  D[1] + D[0] + C[8] + C[9] + C[12] + C[13] )%2;
563 
564  myCRC[3] = ( D[62] + D[61] + D[58] + D[57] + D[47] + D[43] +
565  D[32] + D[31] + D[30] + D[29] + D[17] + D[15] +
566  D[2] + D[1] + C[9] + C[10] + C[13] + C[14] )%2;
567 
568  myCRC[4] = ( D[63] + D[62] + D[59] + D[58] + D[48] + D[44] +
569  D[33] + D[32] + D[31] + D[30] + D[18] + D[16] +
570  D[3] + D[2] + C[0] + C[10] + C[11] + C[14] +
571  C[15] )%2;
572 
573  myCRC[5] = ( D[63] + D[60] + D[59] + D[49] + D[45] + D[34] +
574  D[33] + D[32] + D[31] + D[19] + D[17] + D[4] +
575  D[3] + C[1] + C[11] + C[12] + C[15] )%2;
576 
577  myCRC[6] = ( D[61] + D[60] + D[50] + D[46] + D[35] + D[34] +
578  D[33] + D[32] + D[20] + D[18] + D[5] + D[4] +
579  C[2] + C[12] + C[13] )%2;
580 
581  myCRC[7] = ( D[62] + D[61] + D[51] + D[47] + D[36] + D[35] +
582  D[34] + D[33] + D[21] + D[19] + D[6] + D[5] +
583  C[3] + C[13] + C[14] )%2;
584 
585  myCRC[8] = ( D[63] + D[62] + D[52] + D[48] + D[37] + D[36] +
586  D[35] + D[34] + D[22] + D[20] + D[7] + D[6] +
587  C[0] + C[4] + C[14] + C[15] )%2;
588 
589  myCRC[9] = ( D[63] + D[53] + D[49] + D[38] + D[37] + D[36] +
590  D[35] + D[23] + D[21] + D[8] + D[7] + C[1] +
591  C[5] + C[15] )%2;
592 
593  myCRC[10] = ( D[54] + D[50] + D[39] + D[38] + D[37] + D[36] +
594  D[24] + D[22] + D[9] + D[8] + C[2] + C[6] )%2;
595 
596  myCRC[11] = ( D[55] + D[51] + D[40] + D[39] + D[38] + D[37] +
597  D[25] + D[23] + D[10] + D[9] + C[3] + C[7] )%2;
598 
599  myCRC[12] = ( D[56] + D[52] + D[41] + D[40] + D[39] + D[38] +
600  D[26] + D[24] + D[11] + D[10] + C[4] + C[8] )%2;
601 
602  myCRC[13] = ( D[57] + D[53] + D[42] + D[41] + D[40] + D[39] +
603  D[27] + D[25] + D[12] + D[11] + C[5] + C[9] )%2;
604 
605  myCRC[14] = ( D[58] + D[54] + D[43] + D[42] + D[41] + D[40] +
606  D[28] + D[26] + D[13] + D[12] + C[6] + C[10] )%2;
607 
608  myCRC[15] = ( D[63] + D[62] + D[61] + D[60] + D[59] + D[54] +
609  D[53] + D[52] + D[51] + D[50] + D[49] + D[48] +
610  D[47] + D[46] + D[45] + D[44] + D[42] + D[40] +
611  D[39] + D[38] + D[37] + D[36] + D[35] + D[34] +
612  D[33] + D[32] + D[31] + D[30] + D[29] + D[26] +
613  D[25] + D[24] + D[23] + D[22] + D[21] + D[20] +
614  D[19] + D[18] + D[17] + D[16] + D[15] + D[14] +
615  D[12] + D[11] + D[10] + D[9] + D[8] + D[7] +
616  D[6] + D[5] + D[4] + D[3] + D[2] + D[1] +
617  D[0] + C[0] + C[1] + C[2] + C[3] + C[4] +
618  C[5] + C[6] + C[11] + C[12] + C[13] + C[14] +
619  C[15] )%2;
620 
621  int tempC = 0x0;
622  for(int i=0; i<16 ; i++){ tempC = tempC + (myCRC[i]<<i); }
623  myC = tempC ;
624  return;
625 }
int i
Definition: DBlmapReader.cc:9
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:150
int DTTFFEDReader::channel ( int  wheel,
int  sector,
int  bx 
)
private

Definition at line 470 of file DTTFFEDReader.cc.

References abs.

470  {
471 
472  // wheel : -3 -2 -1 +1 +2 +3 <=> PHTF's : N2, N1, N0, P0, P1, P2
473  // 0 <=> ETTF
474  // sector : 0 -> 11
475  // bx : -1 -> +1
476 
477  int myChannel = 255;
478 
479  if ( abs(bx) > 1) { return myChannel; }
480  if ( sector < 0 || sector > 11) { return myChannel; }
481  if ( abs(wheel) > 3) { return myChannel; }
482 
483  myChannel = sector*21 + wheel*3 - bx + 10 ;
484 
485  if (myChannel > 125) myChannel += 2;
486 
487  return myChannel;
488 }
#define abs(x)
Definition: mlp_lapack.h:159
int wheel(int channel)
int sector(int channel)
void DTTFFEDReader::clear ( void  )
private

Definition at line 451 of file DTTFFEDReader.cc.

References i, j, gen::k, and prof2calltree::l.

Referenced by python.Vispa.Views.WidgetView.WidgetView::closeEvent(), python.Vispa.Views.BoxDecayView.BoxDecayView::closeEvent(), python.Vispa.Share.FindAlgorithm.FindAlgorithm::findUsingFindDialog(), python.Vispa.Views.LineDecayView.LineDecayView::setDataObjects(), python.Vispa.Views.WidgetView.WidgetView::setDataObjects(), python.Vispa.Views.TreeView.TreeView::updateContent(), python.Vispa.Views.TableView.TableView::updateContent(), python.Vispa.Views.BoxDecayView.BoxDecayView::updateContent(), and python.Vispa.Views.PropertyView.PropertyView::updateContent().

451  {
452  phiSegments.clear();
453  theSegments.clear();
454  dtTracks.clear();
455 
456  for(int i=0; i<3; i++){
457  for(int j=0; j<12; j++){
458  for(int k=0; k<6; k++){
459  for(int l=0; l<2; l++){
460  etTrack[i][j][k][l] = 0;
461  efTrack[i][j][k][l] = 0;
462  }
463  }
464  }
465  }
466 
467  return;
468 }
int i
Definition: DBlmapReader.cc:9
L1MuDTChambPhContainer::Phi_Container phiSegments
Definition: DTTFFEDReader.h:76
unsigned int etTrack[3][12][6][2]
Definition: DTTFFEDReader.h:82
int j
Definition: DBlmapReader.cc:9
L1MuDTChambThContainer::The_Container theSegments
Definition: DTTFFEDReader.h:78
L1MuDTTrackContainer::TrackContainer dtTracks
Definition: DTTFFEDReader.h:80
int k[5][pyjets_maxn]
unsigned int efTrack[3][12][6][2]
Definition: DTTFFEDReader.h:84
bool DTTFFEDReader::fillRawData ( edm::Event e,
L1MuDTChambPhContainer::Phi_Container phi_data,
L1MuDTChambThContainer::The_Container the_data,
L1MuDTTrackContainer::TrackContainer tra_data 
)

Generate and fill FED raw data for a full event.

Definition at line 66 of file DTTFFEDReader.cc.

69  {
70 
71  analyse(e);
72 
73  phi_data = p_data();
74  the_data = t_data();
75  tra_data = k_data();
76 
77  return true;
78 }
const L1MuDTChambThContainer::The_Container & t_data()
const L1MuDTTrackContainer::TrackContainer & k_data()
void analyse(edm::Event &e)
const L1MuDTChambPhContainer::Phi_Container & p_data()
edm::InputTag DTTFFEDReader::getDTTFInputTag ( )
inlineprivate

Definition at line 97 of file DTTFFEDReader.h.

References DTTFInputTag.

97 { return DTTFInputTag; }
edm::InputTag DTTFInputTag
Definition: DTTFFEDReader.h:50
const L1MuDTTrackContainer::TrackContainer & DTTFFEDReader::k_data ( )
private

Definition at line 447 of file DTTFFEDReader.cc.

447  {
448  return dtTracks;
449 }
L1MuDTTrackContainer::TrackContainer dtTracks
Definition: DTTFFEDReader.h:80
void DTTFFEDReader::match ( )
private

Definition at line 420 of file DTTFFEDReader.cc.

References i.

420  {
421 
423  i != dtTracks.end();
424  i++ ) {
425  int bxTh = i->bx()+1;
426  int sectorTh = i->scNum();
427  int wheelTh = i->whNum()+3;
428  if(wheelTh > 3) wheelTh-=1;
429  int muonTh = i->TrkTag();
430 
431  i->setEtaPacked(etTrack[bxTh][sectorTh][wheelTh][muonTh]);
432  i->setFineHaloPacked(efTrack[bxTh][sectorTh][wheelTh][muonTh]);
433  }
434 
435  return;
436 }
int i
Definition: DBlmapReader.cc:9
unsigned int etTrack[3][12][6][2]
Definition: DTTFFEDReader.h:82
L1MuDTTrackContainer::TrackContainer dtTracks
Definition: DTTFFEDReader.h:80
unsigned int efTrack[3][12][6][2]
Definition: DTTFFEDReader.h:84
TrackContainer::iterator TrackIterator
const L1MuDTChambPhContainer::Phi_Container & DTTFFEDReader::p_data ( )
private

Definition at line 439 of file DTTFFEDReader.cc.

439  {
440  return phiSegments;
441 }
L1MuDTChambPhContainer::Phi_Container phiSegments
Definition: DTTFFEDReader.h:76
void DTTFFEDReader::process ( edm::Event e)
private

Definition at line 91 of file DTTFFEDReader.cc.

References abs, FEDRawData::data(), runTheMatrix::data, edm::Event::getByLabel(), beamvalidation::lines, TopDecayID::muonID, phi, ExpressReco_HICollisions_FallBack::pt, and FEDRawData::size().

Referenced by python.Vispa.Plugins.ConfigEditor.ConfigDataAccessor.ConfigDataAccessor::dumpPython(), python.Vispa.Plugins.ConfigEditor.ConfigDataAccessor.ConfigDataAccessor::open(), python.Vispa.Plugins.ConfigEditor.ConfigDataAccessor.ConfigDataAccessor::outputEventContent(), python.Vispa.Plugins.ConfigEditor.ConfigDataAccessor.ConfigDataAccessor::setProcess(), and python.Vispa.Plugins.ConfigEditor.ConfigDataAccessor.ConfigDataAccessor::setProperty().

91  {
92 
93  // Container
94  vector<int> DTTFWordContainer;
95  vector<int>::iterator DTTFiterator;
96 
97  // Header constituents
98  int BOEevTy, DTTFId;
99 
100  // DTTF Payload constituents
101  int DTTFWord;
102  int DTTFChan, bitsID;
103  int addr1[2] = { 3, 3};
104  int addr2[2] = {15, 15};
105  int addr3[2] = {15, 15};
106  int addr4[2] = {15, 15};
107 
108  // Trailer constituents
109  int evtLgth , CRC;
110 
111  //--> Header
112 
114  e.getByLabel(getDTTFInputTag(),data);
115  FEDRawData dttfdata = data->FEDData(0x030C);
116  if ( dttfdata.size() == 0 ) return;
117 
118  int* dataWord1 = new int;
119  int* dataWord2 = new int;
120  unsigned char* LineFED=dttfdata.data();
121  *dataWord2=*((int*)LineFED);
122  LineFED+=4;
123  *dataWord1=*((int*)LineFED);
124  int lines = 1; // already counting header
125 
126  BOEevTy = ((*dataWord1)&0xFF000000)>>24; // positions 57 ->64
127  DTTFId = ((*dataWord2)&0x000FFF00)>>8; // positions 9 ->20
128 
129  if( (BOEevTy != 0x50) || ( DTTFId != 0x030C) ){
130  if ( verbose_ ) edm::LogWarning("dttf_unpacker")
131  << "Not a DTTF header " << hex << *dataWord1;
132  delete dataWord1;
133  delete dataWord2;
134  return;
135  }
136 
137  int newCRC = 0xFFFF;
138  calcCRC(*dataWord1, *dataWord2, newCRC);
139 
140 
141  //--> DTTF data
142 
143  LineFED+=4;
144  *dataWord2=*((int*)LineFED);
145  LineFED+=4;
146  *dataWord1=*((int*)LineFED);
147  int chkEOE = ((*dataWord1)&0xFFF00000)>>20;
148  lines++;
149 
150  while(chkEOE != 0xA00){
151 
152  calcCRC(*dataWord1, *dataWord2, newCRC);
153 
154  DTTFWord = *dataWord1;
155  DTTFWordContainer.push_back(DTTFWord);
156  DTTFWord = *dataWord2;
157  DTTFWordContainer.push_back(DTTFWord);
158 
159  LineFED+=4;
160  *dataWord2=*((int*)LineFED);
161  LineFED+=4;
162  *dataWord1=*((int*)LineFED);
163  chkEOE = ((*dataWord1)&0xFFF00000)>>20;
164  lines++;
165 
166  if(lines > 3026){
167  if ( verbose_ ) edm::LogWarning("dttf_unpacker")
168  << "Warning : number of DTTF lines > 3026 "; // 3026 = 1(header) + 3024(max # PHTF-ETTF 64 bits words) + 1(trailer)
169  delete dataWord1;
170  delete dataWord2;
171  return;
172  }
173 
174  } // end while-Data loop
175 
176 
177  //--> Trailer
178 
179  evtLgth = ((*dataWord1)&0x00FFFFFF); // positions 33 ->56
180  CRC = ((*dataWord2)&0xFFFF0000)>>16; // positions 17 ->32
181 
182  calcCRC(*dataWord1, (*dataWord2)&0xFFFF, newCRC);
183  if( newCRC != CRC){
184  if ( verbose_ ) edm::LogWarning("dttf_unpacker")
185  << "Calculated CRC " << hex << newCRC << " differs from CRC in trailer " << hex << CRC;
186  delete dataWord1;
187  delete dataWord2;
188  return;
189  }
190 
191  if( lines != evtLgth){
192  if ( verbose_ ) edm::LogWarning("dttf_unpacker")
193  << "Number of words read != event length " << dec << lines << " " << evtLgth;
194  delete dataWord1;
195  delete dataWord2;
196  return;
197  }
198 
199 
200  // --> analyse event
201 
202  for( DTTFiterator = DTTFWordContainer.begin();
203  DTTFiterator != DTTFWordContainer.end();
204  DTTFiterator++ ){
205 
206  DTTFChan = ((*DTTFiterator)&0xFF000000)>>24;
207  DTTFiterator++;
208  bitsID = ((*DTTFiterator)&0xF0000000)>>28;
209 
210  int bxID = bxNr(DTTFChan);
211  if(bxID == -999) continue;
212  int wheelID = wheel(DTTFChan);
213  if(wheelID == -999) continue;
214  int sectorID = sector(DTTFChan);
215  if(sectorID == -999) continue;
216 
217  //Input
218  if(wheelID!=0 && bitsID<=0x9){
219 
220  int wheelPh = (abs(wheelID)-1)*wheelID/abs(wheelID);
221  int stationID = 0;
222  int ra = 0;
223  int ba = 0;
224  int tsqual = 0;
225  int ts2tag = 0;
226 
227  if ( ( bitsID >> 1 ) == 0 ){ stationID = 1;}
228  if ( ( bitsID >> 1 ) == 1 ){ stationID = 2;}
229  if ( ( bitsID >> 1 ) == 4 ){ stationID = 3;}
230  if ( ( bitsID >> 1 ) == 2 ){ stationID = 4;}
231 
232  if(stationID != 3){
233 
234  ts2tag = (bitsID)&0x1;
235  tsqual = (~(*DTTFiterator)&0x07)-1;
236  ba = (~(*DTTFiterator)&0x1FF8)>>3;
237  if( ba>0x1FF) ba-=0x400;
238  ra = (~(*DTTFiterator)&0x1FFE000)>>13;
239  if( ra>0x7FF) ra-=0x1000;
240  }
241  else{
242 
243  ts2tag = (bitsID)&0x1;
244  tsqual = (~(*DTTFiterator)&0x07)-1;
245  ra = (~(*DTTFiterator)&0x7FF8)>>3;
246  if( ra>0x7FF) ra-=0x1000;
247  }
248 
249  if(tsqual!=7 && wheelID!=-1){
250  phiSegments.push_back(
251  L1MuDTChambPhDigi( bxID+ts2tag, wheelPh, sectorID, stationID,
252  ra, ba, tsqual, ts2tag, 0) );
253  }
254  }
255  //Input
256 
257  //Input
258  if(wheelID==0 && bitsID<=0x4){
259 
260  int wheelTh = bitsID-2;
261 
262  int posALL, posBTI[7];
263 
264  if ( wheelTh == -2 || wheelTh == -1 ||
265  ( wheelTh == 0 && (sectorID == 0 || sectorID == 3 || sectorID == 4 || sectorID == 7 || sectorID == 8 || sectorID == 11) ) ) {
266 
267  posALL = ~(*DTTFiterator)&0x7F;
268  posBTI[0] = ~(*DTTFiterator)&0x01;
269  posBTI[1] = (~(*DTTFiterator)&0x02)>>1;
270  posBTI[2] = (~(*DTTFiterator)&0x04)>>2;
271  posBTI[3] = (~(*DTTFiterator)&0x08)>>3;
272  posBTI[4] = (~(*DTTFiterator)&0x10)>>4;
273  posBTI[5] = (~(*DTTFiterator)&0x20)>>5;
274  posBTI[6] = (~(*DTTFiterator)&0x40)>>6;
275 
276  if(posALL){
277  theSegments.push_back(
278  L1MuDTChambThDigi( bxID, wheelTh, sectorID, 1, posBTI) );
279  }
280 
281  posALL = ~(*DTTFiterator)&0x3F80;
282  posBTI[0] = (~(*DTTFiterator)&0x0080)>>7;
283  posBTI[1] = (~(*DTTFiterator)&0x0100)>>8;
284  posBTI[2] = (~(*DTTFiterator)&0x0200)>>9;
285  posBTI[3] = (~(*DTTFiterator)&0x0400)>>10;
286  posBTI[4] = (~(*DTTFiterator)&0x0800)>>11;
287  posBTI[5] = (~(*DTTFiterator)&0x1000)>>12;
288  posBTI[6] = (~(*DTTFiterator)&0x2000)>>13;
289 
290  if(posALL){
291  theSegments.push_back(
292  L1MuDTChambThDigi( bxID, wheelTh, sectorID, 2, posBTI) );
293  }
294 
295  posALL = ~(*DTTFiterator)&0x1FC000;
296  posBTI[0] = (~(*DTTFiterator)&0x004000)>>14;
297  posBTI[1] = (~(*DTTFiterator)&0x008000)>>15;
298  posBTI[2] = (~(*DTTFiterator)&0x010000)>>16;
299  posBTI[3] = (~(*DTTFiterator)&0x020000)>>17;
300  posBTI[4] = (~(*DTTFiterator)&0x040000)>>18;
301  posBTI[5] = (~(*DTTFiterator)&0x080000)>>19;
302  posBTI[6] = (~(*DTTFiterator)&0x100000)>>20;
303 
304  if(posALL){
305  theSegments.push_back(
306  L1MuDTChambThDigi( bxID, wheelTh, sectorID, 3, posBTI) );
307  }
308  }
309 
310  else {
311 
312  posALL = ~(*DTTFiterator)&0x7F;
313  posBTI[6] = ~(*DTTFiterator)&0x01;
314  posBTI[5] = (~(*DTTFiterator)&0x02)>>1;
315  posBTI[4] = (~(*DTTFiterator)&0x04)>>2;
316  posBTI[3] = (~(*DTTFiterator)&0x08)>>3;
317  posBTI[2] = (~(*DTTFiterator)&0x10)>>4;
318  posBTI[1] = (~(*DTTFiterator)&0x20)>>5;
319  posBTI[0] = (~(*DTTFiterator)&0x40)>>6;
320 
321  if(posALL){
322  theSegments.push_back(
323  L1MuDTChambThDigi( bxID, wheelTh, sectorID, 1, posBTI) );
324  }
325 
326  posALL = ~(*DTTFiterator)&0x3F80;
327  posBTI[6] = (~(*DTTFiterator)&0x0080)>>7;
328  posBTI[5] = (~(*DTTFiterator)&0x0100)>>8;
329  posBTI[4] = (~(*DTTFiterator)&0x0200)>>9;
330  posBTI[3] = (~(*DTTFiterator)&0x0400)>>10;
331  posBTI[2] = (~(*DTTFiterator)&0x0800)>>11;
332  posBTI[1] = (~(*DTTFiterator)&0x1000)>>12;
333  posBTI[0] = (~(*DTTFiterator)&0x2000)>>13;
334 
335  if(posALL){
336  theSegments.push_back(
337  L1MuDTChambThDigi( bxID, wheelTh, sectorID, 2, posBTI) );
338  }
339 
340  posALL = ~(*DTTFiterator)&0x1FC000;
341  posBTI[6] = (~(*DTTFiterator)&0x004000)>>14;
342  posBTI[5] = (~(*DTTFiterator)&0x008000)>>15;
343  posBTI[4] = (~(*DTTFiterator)&0x010000)>>16;
344  posBTI[3] = (~(*DTTFiterator)&0x020000)>>17;
345  posBTI[2] = (~(*DTTFiterator)&0x040000)>>18;
346  posBTI[1] = (~(*DTTFiterator)&0x080000)>>19;
347  posBTI[0] = (~(*DTTFiterator)&0x100000)>>20;
348 
349  if(posALL){
350  theSegments.push_back(
351  L1MuDTChambThDigi( bxID, wheelTh, sectorID, 3, posBTI) );
352  }
353  }
354  }
355  //Input
356 
357  //Addresses
358  if(wheelID!=0 && bitsID>=0xA && bitsID<=0xB){
359 
360  int candID = bitsID - 0xA;
361 
362  addr4[candID] = ((*DTTFiterator)&0x0F);
363  addr3[candID] = ((*DTTFiterator)&0xF0)>>4;
364  addr2[candID] = ((*DTTFiterator)&0xF00)>>8;
365  addr1[candID] = ((*DTTFiterator)&0x3000)>>12;
366  }
367  //Addresses
368 
369  //Output
370  if(wheelID!=0 && bitsID>=0xC){
371 
372  int muonID = 0;
373  int pt = 0;
374  int ch = 0;
375  int phi = 0;
376  int qual = 0;
377 
378  muonID = (bitsID&0x1);
379  qual = (~(*DTTFiterator)&0x07);
380  phi = ((*DTTFiterator)&0x7F8)>>3;
381  ch = (~(*DTTFiterator)&0x800)>>11;
382  pt = (~(*DTTFiterator)&0x1F000)>>12;
383 
384  if(qual!=0){
385  dtTracks.push_back(
386  L1MuDTTrackCand( 0, phi, 0, pt, ch, 1, 0, qual, bxID, wheelID, sectorID,
387  muonID, addr1[muonID], addr2[muonID], addr3[muonID], addr4[muonID] ) );
388  }
389  }
390  //Output
391 
392  //Output
393  if(wheelID==0 && bitsID>=0x8){
394 
395  int wheelTh = bitsID&0x7;
396 
397  int etaALL;
398 
399  etaALL = ~(*DTTFiterator)&0x007F;
400  if (etaALL) {
401  etTrack[bxID+1][sectorID][wheelTh][0] = (*DTTFiterator)&0x003F;
402  efTrack[bxID+1][sectorID][wheelTh][0] = (~(*DTTFiterator)&0x0040)>>6;
403  }
404 
405  etaALL = (~(*DTTFiterator)&0x3F80)>>7;
406  if (etaALL) {
407  etTrack[bxID+1][sectorID][wheelTh][1]= ( (*DTTFiterator)&0x1F80)>>7;
408  efTrack[bxID+1][sectorID][wheelTh][1]= (~(*DTTFiterator)&0x2000)>>13;
409  }
410  }
411  //Output
412 
413  } // end for-loop container content
414 
415  delete dataWord1;
416  delete dataWord2;
417  return;
418 }
L1MuDTChambPhContainer::Phi_Container phiSegments
Definition: DTTFFEDReader.h:76
#define abs(x)
Definition: mlp_lapack.h:159
unsigned int etTrack[3][12][6][2]
Definition: DTTFFEDReader.h:82
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:47
L1MuDTChambThContainer::The_Container theSegments
Definition: DTTFFEDReader.h:78
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
L1MuDTTrackContainer::TrackContainer dtTracks
Definition: DTTFFEDReader.h:80
unsigned int efTrack[3][12][6][2]
Definition: DTTFFEDReader.h:84
static const int muonID
Definition: TopGenEvent.h:20
void calcCRC(int myD1, int myD2, int &myC)
edm::InputTag getDTTFInputTag()
Definition: DTTFFEDReader.h:97
const unsigned char * data() const
Return a const pointer to the beginning of the data buffer.
Definition: FEDRawData.cc:29
int wheel(int channel)
int bxNr(int channel)
int sector(int channel)
Definition: DDAxes.h:10
void DTTFFEDReader::produce ( edm::Event e,
const edm::EventSetup c 
)
virtual

Produce digis out of raw data.

Implements edm::EDProducer.

Definition at line 44 of file DTTFFEDReader.cc.

References edm::Event::put().

Referenced by python.JSONExport.JsonExport::export(), and python.HTMLExport.HTMLExport::export().

44  {
45 
46  auto_ptr<L1MuDTChambPhContainer> phi_product(new L1MuDTChambPhContainer);
47  auto_ptr<L1MuDTChambThContainer> the_product(new L1MuDTChambThContainer);
48  auto_ptr<L1MuDTTrackContainer> tra_product(new L1MuDTTrackContainer);
49 
53 
54  if (!fillRawData(e, phi_data, the_data, tra_data)) return;
55 
56  phi_product->setContainer(phi_data);
57  the_product->setContainer(the_data);
58  tra_product->setContainer(tra_data);
59 
60  e.put(phi_product);
61  e.put(the_product);
62  e.put(tra_product,"DATA");
63 
64 }
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
std::vector< L1MuDTTrackCand > TrackContainer
std::vector< L1MuDTChambPhDigi > Phi_Container
std::vector< L1MuDTChambThDigi > The_Container
bool fillRawData(edm::Event &e, L1MuDTChambPhContainer::Phi_Container &phi_data, L1MuDTChambThContainer::The_Container &the_data, L1MuDTTrackContainer::TrackContainer &tra_data)
Generate and fill FED raw data for a full event.
int DTTFFEDReader::sector ( int  channel)
private

Definition at line 503 of file DTTFFEDReader.cc.

Referenced by geometryXMLparser.DTAlignable::index().

503  {
504 
505  int myChannel = channel;
506 
507  if (myChannel > 127) myChannel -= 2;
508 
509  if (myChannel < 0 || myChannel > 251 ){ return -999; }
510 
511  return myChannel/21;
512 }
int channel(int wheel, int sector, int bx)
const L1MuDTChambThContainer::The_Container & DTTFFEDReader::t_data ( )
private

Definition at line 443 of file DTTFFEDReader.cc.

443  {
444  return theSegments;
445 }
L1MuDTChambThContainer::The_Container theSegments
Definition: DTTFFEDReader.h:78
int DTTFFEDReader::wheel ( int  channel)
private

Definition at line 514 of file DTTFFEDReader.cc.

Referenced by geometryXMLparser.DTAlignable::index().

514  {
515 
516  int myChannel = channel;
517 
518  if (myChannel > 127) myChannel -= 2;
519 
520  if (myChannel < 0 || myChannel > 251 ){ return -999; }
521 
522  int myWheel = ((myChannel%21)/3)-3;
523 
524  return myWheel;
525 }
int channel(int wheel, int sector, int bx)

Member Data Documentation

edm::InputTag DTTFFEDReader::DTTFInputTag
private

Definition at line 50 of file DTTFFEDReader.h.

Referenced by getDTTFInputTag().

L1MuDTTrackContainer::TrackContainer DTTFFEDReader::dtTracks
private

Definition at line 80 of file DTTFFEDReader.h.

unsigned int DTTFFEDReader::efTrack[3][12][6][2]
private

Definition at line 84 of file DTTFFEDReader.h.

unsigned int DTTFFEDReader::etTrack[3][12][6][2]
private

Definition at line 82 of file DTTFFEDReader.h.

L1MuDTChambPhContainer::Phi_Container DTTFFEDReader::phiSegments
private

Definition at line 76 of file DTTFFEDReader.h.

L1MuDTChambThContainer::The_Container DTTFFEDReader::theSegments
private

Definition at line 78 of file DTTFFEDReader.h.

bool DTTFFEDReader::verbose_
private

Definition at line 52 of file DTTFFEDReader.h.