CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
L1TMuon::DTCollector Class Reference

#include <DTCollector.h>

Inheritance diagram for L1TMuon::DTCollector:
L1TMuon::SubsystemCollector

Public Member Functions

 DTCollector (const edm::ParameterSet &)
 
void extractPrimitives (const edm::Event &, const edm::EventSetup &, std::vector< TriggerPrimitive > &) const override
 
 ~DTCollector () override
 
- Public Member Functions inherited from L1TMuon::SubsystemCollector
 SubsystemCollector (const edm::ParameterSet &)
 
virtual ~SubsystemCollector ()
 

Private Member Functions

int findBTIGroupForThetaDigi (const L1MuDTChambThDigi &, const int position) const
 
TriggerPrimitive processDigis (const L1MuDTChambPhDigi &, const int &segment_number) const
 
TriggerPrimitive processDigis (const L1MuDTChambThDigi &, const int bti_group) const
 
TriggerPrimitive processDigis (const L1MuDTChambPhDigi &, const L1MuDTChambThDigi &, const int bti_group) const
 

Private Attributes

std::unique_ptr< DTBunchCrossingCleaner_bxc
 
const int bx_max
 
const int bx_min
 

Additional Inherited Members

- Protected Attributes inherited from L1TMuon::SubsystemCollector
edm::InputTag _src
 

Detailed Description

Definition at line 22 of file DTCollector.h.

Constructor & Destructor Documentation

DTCollector::DTCollector ( const edm::ParameterSet ps)

Definition at line 13 of file DTCollector.cc.

References _bxc, edm::ParameterSet::getParameter(), and edm::ParameterSet::getParameterSet().

13  :
15  bx_min(ps.getParameter<int>("BX_min")),
16  bx_max(ps.getParameter<int>("BX_max")) {
17  if( ps.getParameter<bool>("runBunchCrossingCleaner") ) {
18  const edm::ParameterSet& bxccfg = ps.getParameterSet("bxCleanerCfg");
19  _bxc.reset(new DTBunchCrossingCleaner(bxccfg));
20  } else {
21  _bxc.reset(nullptr);
22  }
23 }
T getParameter(std::string const &) const
std::unique_ptr< DTBunchCrossingCleaner > _bxc
Definition: DTCollector.h:40
ParameterSet const & getParameterSet(std::string const &) const
SubsystemCollector(const edm::ParameterSet &)
L1TMuon::DTCollector::~DTCollector ( )
inlineoverride

Definition at line 25 of file DTCollector.h.

References extractPrimitives(), findBTIGroupForThetaDigi(), position, and processDigis().

25 {}

Member Function Documentation

void DTCollector::extractPrimitives ( const edm::Event ,
const edm::EventSetup ,
std::vector< TriggerPrimitive > &   
) const
overridevirtual

Implements L1TMuon::SubsystemCollector.

Definition at line 25 of file DTCollector.cc.

References _bxc, L1TMuon::SubsystemCollector::_src, bx_max, bx_min, L1MuDTChambPhContainer::chPhiSegm1(), L1MuDTChambPhContainer::chPhiSegm2(), L1MuDTChambThContainer::chThetaSegm(), findBTIGroupForThetaDigi(), edm::Event::getByLabel(), processDigis(), relativeConstraints::station, groupFilesInBlocks::temp, and makeMuonMisalignmentScenario::wheel.

Referenced by ~DTCollector().

27  {
28  TriggerPrimitiveCollection cleaned, temp, chamb_list;
31  ev.getByLabel(_src,phiDigis);
32  ev.getByLabel(_src,thetaDigis);
33  for( int wheel = -2; wheel <= 2 ; ++wheel ) {
34  for( int station = 1; station <= 4; ++station ) {
35  for( int sector = 0; sector <= 11; ++sector ) {
36  chamb_list.clear();
37  for( int bx = bx_min; bx <= bx_max; ++bx) {
38  std::unique_ptr<const L1MuDTChambPhDigi> phi_segm_1(
39  phiDigis->chPhiSegm1(wheel,station,sector,bx)
40  );
41  std::unique_ptr<const L1MuDTChambPhDigi> phi_segm_2(
42  phiDigis->chPhiSegm2(wheel,station,sector,bx)
43  );
44  std::unique_ptr<const L1MuDTChambThDigi> theta_segm(
45  thetaDigis->chThetaSegm(wheel,station,sector,bx)
46  );
47 
48  int bti_group_1=-1, bti_group_2=-1;
49 
50  if( theta_segm ) {
51  bti_group_1 = findBTIGroupForThetaDigi(*theta_segm,1);
52  bti_group_2 = findBTIGroupForThetaDigi(*theta_segm,2);
53  }
54 
55  if( phi_segm_1 && bti_group_1 != -1 ) {
56  chamb_list.push_back(processDigis(*phi_segm_1,
57  *theta_segm,
58  bti_group_1));
59  } else if ( phi_segm_1 && bti_group_1 == -1 ) {
60  chamb_list.push_back(processDigis(*phi_segm_1,1));
61  } else if ( !phi_segm_1 && bti_group_1 != -1 ) {
62  chamb_list.push_back(processDigis(*theta_segm,
63  bti_group_1));
64  }
65 
66  if( phi_segm_2 && bti_group_2 != -1) {
67  chamb_list.push_back(processDigis(*phi_segm_2,
68  *theta_segm,
69  bti_group_2));
70  } else if ( phi_segm_2 && bti_group_2 == -1 ) {
71  chamb_list.push_back(processDigis(*phi_segm_2,2));
72  } else if ( !phi_segm_2 && bti_group_2 != -1 ) {
73  chamb_list.push_back(processDigis(*phi_segm_2,bti_group_2));
74  }
75 
76  phi_segm_1.release();
77  phi_segm_2.release();
78  theta_segm.release();
79  }
80  if( _bxc ) {
81  temp = _bxc->clean(chamb_list);
82  cleaned.insert(cleaned.end(),temp.begin(),temp.end());
83  } else {
84  cleaned.insert(cleaned.end(),chamb_list.begin(),chamb_list.end());
85  }
86  }
87  }
88  }
89  out.insert(out.end(),cleaned.begin(),cleaned.end());
90 }
L1MuDTChambPhDigi const * chPhiSegm1(int wheel, int stat, int sect, int bx) const
L1MuDTChambPhDigi const * chPhiSegm2(int wheel, int stat, int sect, int bx) const
bool ev
std::unique_ptr< DTBunchCrossingCleaner > _bxc
Definition: DTCollector.h:40
int findBTIGroupForThetaDigi(const L1MuDTChambThDigi &, const int position) const
Definition: DTCollector.cc:112
TriggerPrimitive processDigis(const L1MuDTChambPhDigi &, const int &segment_number) const
Definition: DTCollector.cc:92
std::vector< TriggerPrimitive > TriggerPrimitiveCollection
L1MuDTChambThDigi const * chThetaSegm(int wheel, int stat, int sect, int bx) const
int DTCollector::findBTIGroupForThetaDigi ( const L1MuDTChambThDigi digi,
const int  position 
) const
private

Definition at line 112 of file DTCollector.cc.

References DEFINE_EDM_PLUGIN, mps_fire::i, L1MuDTChambThDigi::position(), and mps_fire::result.

Referenced by extractPrimitives(), processDigis(), and ~DTCollector().

113  {
114  //if( digi.stNum() == 4 ) return -1; // there is no theta layer there
115  int result = -1;
116  for( int i = 0; i < 7; ++i ) {
117  if( digi.position(i) == pos ) result = i;
118  }
119  return result;
120 }
int position(const int i) const
TriggerPrimitive DTCollector::processDigis ( const L1MuDTChambPhDigi digi,
const int &  segment_number 
) const
private

Definition at line 92 of file DTCollector.cc.

References L1MuDTChambPhDigi::scNum(), L1MuDTChambPhDigi::stNum(), and L1MuDTChambPhDigi::whNum().

Referenced by extractPrimitives(), and ~DTCollector().

93  {
94  DTChamberId detid(digi.whNum(),digi.stNum(),digi.scNum()+1);
95  return TriggerPrimitive(detid,digi,segment_number);
96 }
L1TMuon::TriggerPrimitive TriggerPrimitive
Definition: Common.h:33
TriggerPrimitive DTCollector::processDigis ( const L1MuDTChambThDigi digi_th,
const int  bti_group 
) const
private

Definition at line 98 of file DTCollector.cc.

References L1MuDTChambThDigi::scNum(), L1MuDTChambThDigi::stNum(), and L1MuDTChambThDigi::whNum().

99  {
100  DTChamberId detid(digi_th.whNum(),digi_th.stNum(),digi_th.scNum()+1);
101  return TriggerPrimitive(detid,digi_th,bti_group);
102 }
L1TMuon::TriggerPrimitive TriggerPrimitive
Definition: Common.h:33
TriggerPrimitive DTCollector::processDigis ( const L1MuDTChambPhDigi digi_phi,
const L1MuDTChambThDigi digi_theta,
const int  bti_group 
) const
private

Definition at line 104 of file DTCollector.cc.

References findBTIGroupForThetaDigi(), L1MuDTChambPhDigi::scNum(), L1MuDTChambPhDigi::stNum(), and L1MuDTChambPhDigi::whNum().

106  {
107  DTChamberId detid(digi_phi.whNum(),digi_phi.stNum(),digi_phi.scNum()+1);
108  return TriggerPrimitive(detid,digi_phi,digi_theta,bti_group);
109 }
L1TMuon::TriggerPrimitive TriggerPrimitive
Definition: Common.h:33

Member Data Documentation

std::unique_ptr<DTBunchCrossingCleaner> L1TMuon::DTCollector::_bxc
private

Definition at line 40 of file DTCollector.h.

Referenced by DTCollector(), and extractPrimitives().

const int L1TMuon::DTCollector::bx_max
private

Definition at line 39 of file DTCollector.h.

Referenced by extractPrimitives().

const int L1TMuon::DTCollector::bx_min
private

Definition at line 39 of file DTCollector.h.

Referenced by extractPrimitives().