CMS 3D CMS Logo

ME0Motherboard.cc
Go to the documentation of this file.
5 
7  const edm::ParameterSet& conf) :
8  theEndcap(endcap), theChamber(chamber)
9 {
10  edm::ParameterSet tmbParams = conf.getParameter<edm::ParameterSet>("tmbParam");
11  infoV = tmbParams.getParameter<int>("verbosity");
12 }
13 
15  theEndcap(1), theChamber(1)
16 {
17  infoV = 2;
18 }
19 
21 }
22 
24 {
25  for (int bx = 0; bx < MAX_TRIGGER_BINS; bx++) {
26  for (int i = 0; i < MAX_TRIGGERS; i++) {
27  Triggers[bx][i].clear();
28  }
29  }
30 }
31 
33 {
34  std::unique_ptr<ME0PadDigiCollection> me0Pads(new ME0PadDigiCollection());
35  declusterize(me0Clusters, *me0Pads);
36  run(me0Pads.get());
37 }
38 
40 {
41  clear();
42 }
43 
44 // Returns vector of read-out correlated Triggers, if any. Starts with
45 // the vector of all found Triggers and selects the ones in the read-out
46 // time window.
47 std::vector<ME0TriggerDigi> ME0Motherboard::readoutTriggers()
48 {
49  std::vector<ME0TriggerDigi> tmpV;
50 
51  std::vector<ME0TriggerDigi> all_trigs = getTriggers();
52  for (const auto& ptrig : all_trigs) {
53  // in the future, add a selection on the BX
54  tmpV.push_back(ptrig);
55  }
56  return tmpV;
57 }
58 
59 // Returns vector of all found correlated Triggers, if any.
60 std::vector<ME0TriggerDigi> ME0Motherboard::getTriggers()
61 {
62  std::vector<ME0TriggerDigi> tmpV;
63 
64  for (int bx = 0; bx < MAX_TRIGGER_BINS; bx++) {
65  for (int i = 0; i < MAX_TRIGGERS; i++) {
66  tmpV.push_back(Triggers[bx][i]);
67  }
68  }
69  return tmpV;
70 }
71 
72 // compare Triggers by quality
74 {
75  return trig1.getQuality() > trig2.getQuality();
76 }
77 
78 // compare Triggers by GEM bending angle
80 {
81  // todo: In the future I plan a member to be added to ME0TriggerDigi, getME0Dphi().
82  // That function will derive the bending angle from the pattern.
83  // The ME0TriggerDigi pattterns are at this point not defined yet.
84  return true;
85 }
86 
87 
88 
90  ME0PadDigiCollection& out_pads)
91 {
92  for (auto detUnitIt = in_clusters->begin();detUnitIt != in_clusters->end(); ++detUnitIt) {
93  const ME0DetId& id = (*detUnitIt).first;
94  const auto& range = (*detUnitIt).second;
95  for (auto digiIt = range.first; digiIt!=range.second; ++digiIt) {
96  for (const auto& p: digiIt->pads()){
97  out_pads.insertDigi(id, ME0PadDigi(p, digiIt->bx()));
98  }
99  }
100  }
101 }
T getParameter(std::string const &) const
static bool sortByME0Dphi(const ME0TriggerDigi &, const ME0TriggerDigi &)
static bool sortByQuality(const ME0TriggerDigi &, const ME0TriggerDigi &)
int getQuality() const
return the Quality
MuonDigiCollection< ME0DetId, ME0PadDigi > ME0PadDigiCollection
void declusterize(const ME0PadDigiClusterCollection *, ME0PadDigiCollection &)
std::vector< ME0TriggerDigi > getTriggers()
const unsigned theEndcap
std::vector< ME0TriggerDigi > readoutTriggers()
void clear()
clear this Trigger
const unsigned theChamber
ME0TriggerDigi Triggers[MAX_TRIGGER_BINS][MAX_TRIGGERS]
void run(const ME0PadDigiCollection *)