CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
MuonPathConfirmator Class Reference

#include <MuonPathConfirmator.h>

Public Member Functions

void finish ()
 
void initialise (const edm::EventSetup &iEventSetup)
 
 MuonPathConfirmator (const edm::ParameterSet &pset, edm::ConsumesCollector &iC)
 
void run (edm::Event &iEvent, const edm::EventSetup &iEventSetup, std::vector< cmsdt::metaPrimitive > inMetaPrimitives, edm::Handle< DTDigiCollection > dtdigis, std::vector< cmsdt::metaPrimitive > &outMetaPrimitives)
 
 ~MuonPathConfirmator ()
 

Private Member Functions

void analyze (cmsdt::metaPrimitive mp, edm::Handle< DTDigiCollection > dtdigis, std::vector< cmsdt::metaPrimitive > &outMetaPrimitives)
 

Private Attributes

bool debug_
 
int LYRANDAHALF_RES_SHR = 4
 
int max_drift_tdc = -1
 
edm::FileInPath maxdrift_filename_
 
int maxdriftinfo_ [5][4][14]
 
double minx_match_2digis_
 
int PARTIALS_PRECISSION = 4
 
int SEMICHAMBER_H = int(SEMICHAMBER_H_REAL)
 
int SEMICHAMBER_H_PRECISSION = 13 + PARTIALS_PRECISSION
 
float SEMICHAMBER_H_REAL = ((235. / 2.) / (16. * 6.5)) * std::pow(2, SEMICHAMBER_H_PRECISSION)
 
int SEMICHAMBER_RES_SHR = SEMICHAMBER_H_PRECISSION
 
edm::FileInPath shift_filename_
 
std::map< int, float > shiftinfo_
 

Detailed Description

Definition at line 32 of file MuonPathConfirmator.h.

Constructor & Destructor Documentation

◆ MuonPathConfirmator()

MuonPathConfirmator::MuonPathConfirmator ( const edm::ParameterSet pset,
edm::ConsumesCollector iC 
)

Definition at line 12 of file MuonPathConfirmator.cc.

References debug_, Exception, edm::FileInPath::fullPath(), LogDebug, maxdrift_filename_, maxdriftinfo_, muonDTDigis_cfi::pset, nano_mu_digi_cff::rawId, edm::shift, shift_filename_, and shiftinfo_.

13  : debug_(pset.getUntrackedParameter<bool>("debug")),
14  minx_match_2digis_(pset.getParameter<double>("minx_match_2digis")) {
15  if (debug_)
16  LogDebug("MuonPathConfirmator") << "MuonPathConfirmator: constructor";
17 
18  //shift phi
19  int rawId;
20  shift_filename_ = pset.getParameter<edm::FileInPath>("shift_filename");
21  std::ifstream ifin3(shift_filename_.fullPath());
22  double shift;
23  if (ifin3.fail()) {
24  throw cms::Exception("Missing Input File")
25  << "MuonPathConfirmator::MuonPathConfirmator() - Cannot find " << shift_filename_.fullPath();
26  }
27  while (ifin3.good()) {
28  ifin3 >> rawId >> shift;
30  }
31 
32  int wh, st, se, maxdrift;
33  maxdrift_filename_ = pset.getParameter<edm::FileInPath>("maxdrift_filename");
34  std::ifstream ifind(maxdrift_filename_.fullPath());
35  if (ifind.fail()) {
36  throw cms::Exception("Missing Input File")
37  << "MPSLFilter::MPSLFilter() - Cannot find " << maxdrift_filename_.fullPath();
38  }
39  while (ifind.good()) {
40  ifind >> wh >> st >> se >> maxdrift;
41  maxdriftinfo_[wh][st][se] = maxdrift;
42  }
43 }
edm::FileInPath maxdrift_filename_
edm::FileInPath shift_filename_
int maxdriftinfo_[5][4][14]
std::map< int, float > shiftinfo_
static unsigned int const shift
const std::string & fullPath() const
Definition: FileInPath.cc:144
#define LogDebug(id)

◆ ~MuonPathConfirmator()

MuonPathConfirmator::~MuonPathConfirmator ( )

Definition at line 45 of file MuonPathConfirmator.cc.

References debug_, and LogDebug.

45  {
46  if (debug_)
47  LogDebug("MuonPathConfirmator") << "MuonPathAnalyzer: destructor";
48 }
#define LogDebug(id)

Member Function Documentation

◆ analyze()

void MuonPathConfirmator::analyze ( cmsdt::metaPrimitive  mp,
edm::Handle< DTDigiCollection dtdigis,
std::vector< cmsdt::metaPrimitive > &  outMetaPrimitives 
)
private

Definition at line 89 of file MuonPathConfirmator.cc.

References funct::abs(), cmsdt::CELL_SEMILENGTH, cmsdt::metaPrimitive::chi2, IntegrityClient_cfi::ChId, cmsdt::CHIGHQ, cmsdt::CLOWQ, createfilelist::int, cmsdt::metaPrimitive::lat1, cmsdt::metaPrimitive::lat2, cmsdt::metaPrimitive::lat3, cmsdt::metaPrimitive::lat4, cmsdt::metaPrimitive::lat5, cmsdt::metaPrimitive::lat6, cmsdt::metaPrimitive::lat7, cmsdt::metaPrimitive::lat8, DTLayerId::layer(), cmsdt::LHC_CLK_FREQ, cmsdt::LOWQ, LYRANDAHALF_RES_SHR, max_drift_tdc, minx_match_2digis_, PARTIALS_PRECISSION, cmsdt::metaPrimitive::phi, cmsdt::metaPrimitive::phi_cmssw, cmsdt::metaPrimitive::phiB, cmsdt::metaPrimitive::phiB_cmssw, cmsdt::metaPrimitive::quality, DetId::rawId(), cmsdt::metaPrimitive::rawId, edm::second(), DTChamberId::sector(), SEMICHAMBER_H, SEMICHAMBER_RES_SHR, shiftinfo_, cmsdt::SL1_CELLS_OFFSET, DTChamberId::station(), DTSuperLayerId::superLayer(), cmsdt::metaPrimitive::t0, cmsdt::metaPrimitive::tanPhi, cmsdt::metaPrimitive::tdc1, cmsdt::metaPrimitive::tdc2, cmsdt::metaPrimitive::tdc3, cmsdt::metaPrimitive::tdc4, cmsdt::metaPrimitive::tdc5, cmsdt::metaPrimitive::tdc6, cmsdt::metaPrimitive::tdc7, cmsdt::metaPrimitive::tdc8, cmsdt::TIME_TO_TDC_COUNTS, DTChamberId::wheel(), cmsdt::metaPrimitive::wi1, cmsdt::metaPrimitive::wi2, cmsdt::metaPrimitive::wi3, cmsdt::metaPrimitive::wi4, cmsdt::metaPrimitive::wi5, cmsdt::metaPrimitive::wi6, cmsdt::metaPrimitive::wi7, cmsdt::metaPrimitive::wi8, and cmsdt::metaPrimitive::x.

Referenced by run().

91  {
92  int dum_sl_rawid = mp.rawId;
93  DTSuperLayerId dumSlId(dum_sl_rawid);
94  DTChamberId ChId(dumSlId.wheel(), dumSlId.station(), dumSlId.sector());
95  DTSuperLayerId sl1Id(ChId.rawId(), 1);
96  DTSuperLayerId sl3Id(ChId.rawId(), 3);
97 
98  DTWireId wireIdSL1(sl1Id, 2, 1);
99  DTWireId wireIdSL3(sl3Id, 2, 1);
100  auto sl_shift_cm = shiftinfo_[wireIdSL1.rawId()] - shiftinfo_[wireIdSL3.rawId()];
101  bool isSL1 = (mp.rawId == sl1Id.rawId());
102  bool isSL3 = (mp.rawId == sl3Id.rawId());
103  if (!isSL1 && !isSL3)
104  outMetaPrimitives.emplace_back(mp);
105  else {
106  int best_tdc = -1;
107  int next_tdc = -1;
108  int best_wire = -1;
109  int next_wire = -1;
110  int best_layer = -1;
111  int next_layer = -1;
112  int best_lat = -1;
113  int next_lat = -1;
114  int lat = -1;
115  int matched_digis = 0;
116 
117  int position_prec = ((int)(mp.x)) << PARTIALS_PRECISSION;
118  int slope_prec = ((int)(mp.tanPhi)) << PARTIALS_PRECISSION;
119 
120  int slope_x_halfchamb = (((long int)slope_prec) * SEMICHAMBER_H) >> SEMICHAMBER_RES_SHR;
121  int slope_x_3semicells = (slope_prec * 3) >> LYRANDAHALF_RES_SHR;
122  int slope_x_1semicell = (slope_prec * 1) >> LYRANDAHALF_RES_SHR;
123 
124  for (const auto &dtLayerId_It : *dtdigis) {
125  const DTLayerId dtLId = dtLayerId_It.first;
126  // creating a new DTSuperLayerId object to compare with the required SL id
127  const DTSuperLayerId dtSLId(dtLId.wheel(), dtLId.station(), dtLId.sector(), dtLId.superLayer());
128  bool hitFromSL1 = (dtSLId.rawId() == sl1Id.rawId());
129  bool hitFromSL3 = (dtSLId.rawId() == sl3Id.rawId());
130  if (!(hitFromSL1 || hitFromSL3)) // checking hits are from one of the other SL of the same chamber
131  continue;
132  double minx = 10 * minx_match_2digis_ * ((double)max_drift_tdc / (double)CELL_SEMILENGTH);
133  double min2x = 10 * minx_match_2digis_ * ((double)max_drift_tdc / (double)CELL_SEMILENGTH);
134  if (isSL1 != hitFromSL1) { // checking hits have the opposite SL than the TP
135  for (auto digiIt = (dtLayerId_It.second).first; digiIt != (dtLayerId_It.second).second; ++digiIt) {
136  if ((*digiIt).time() < mp.t0)
137  continue;
138  int wp_semicells = ((*digiIt).wire() - 1 - SL1_CELLS_OFFSET) * 2 + 1;
139  int ly = dtLId.layer() - 1;
140  if (ly % 2 == 1)
141  wp_semicells -= 1;
142  if (hitFromSL3)
143  wp_semicells -= (int)round((sl_shift_cm * 10) / CELL_SEMILENGTH);
144  double hit_position = wp_semicells * max_drift_tdc +
145  ((*digiIt).time() - mp.t0) * (double)TIME_TO_TDC_COUNTS / (double)LHC_CLK_FREQ;
146  double hit_position_left = wp_semicells * max_drift_tdc -
147  ((*digiIt).time() - mp.t0) * (double)TIME_TO_TDC_COUNTS / (double)LHC_CLK_FREQ;
148  // extrapolating position to the layer of the hit
149  // mp.position is referred to the center between SLs, so one has to add half the distance between SLs
150  // + half a cell height to get to the first wire + ly * cell height to reach the desired ly
151  // 10 * VERT_PHI1_PHI3 / 2 + (CELL_HEIGHT / 2) + ly * CELL_HEIGHT = (10 * VERT_PHI1_PHI3 + (2 * ly + 1) * CELL_HEIGHT) / 2
152 
153  int position_in_layer = position_prec + (1 - 2 * (int)hitFromSL1) * slope_x_halfchamb;
154  if (ly == 0)
155  position_in_layer -= slope_x_3semicells;
156  if (ly == 1)
157  position_in_layer -= slope_x_1semicell;
158  if (ly == 2)
159  position_in_layer += slope_x_1semicell;
160  if (ly == 3)
161  position_in_layer += slope_x_3semicells;
162  position_in_layer = position_in_layer >> PARTIALS_PRECISSION;
163 
164  if (std::abs(position_in_layer - hit_position_left) < std::abs(position_in_layer - hit_position)) {
165  lat = 0;
166  hit_position = hit_position_left;
167  }
168  if (std::abs(position_in_layer - hit_position) < minx) {
169  // different layer than the stored in best, hit added, matched_digis++;. This approach in somewhat
170  // buggy, as we could have stored as best LayerX -> LayerY -> LayerX, and this should
171  // count only as 2 hits. However, as we confirm with at least 2 hits, having 2 or more
172  // makes no difference
173  if (dtLId.layer() != best_layer) {
174  minx = std::abs(position_in_layer - hit_position);
175  next_wire = best_wire;
176  next_tdc = best_tdc;
177  next_layer = best_layer;
178  next_lat = best_lat;
179  matched_digis++;
180  }
181  best_wire = (*digiIt).wire() - 1;
182  best_tdc = (*digiIt).time();
183  best_layer = dtLId.layer();
184  best_lat = lat;
185 
186  } else if ((std::abs(position_in_layer - hit_position) >= minx) &&
187  (std::abs(position_in_layer - hit_position) < min2x)) {
188  // same layer than the stored in best, no hit added
189  if (dtLId.layer() == best_layer)
190  continue;
191  // different layer than the stored in next, hit added. This approach in somewhat
192  // buggy, as we could have stored as next LayerX -> LayerY -> LayerX, and this should
193  // count only as 2 hits. However, as we confirm with at least 2 hits, having 2 or more
194  // makes no difference
195  matched_digis++;
196  // whether the layer is the same for this hit and the stored in next, we substitute
197  // the one stored and modify the min distance
198  min2x = std::abs(position_in_layer - hit_position);
199  next_wire = (*digiIt).wire() - 1;
200  next_tdc = (*digiIt).time();
201  next_layer = dtLId.layer();
202  next_lat = lat;
203  }
204  }
205  }
206  }
207  int new_quality = mp.quality;
208  std::vector<int> wi_c(4, -1), tdc_c(4, -1), lat_c(4, -1);
209  if (matched_digis >= 2 and best_layer != -1 and next_layer != -1) { // actually confirm
210  new_quality = CHIGHQ;
211  if (mp.quality == LOWQ)
212  new_quality = CLOWQ;
213 
214  wi_c[next_layer - 1] = next_wire;
215  tdc_c[next_layer - 1] = next_tdc;
216  lat_c[next_layer - 1] = next_lat;
217 
218  wi_c[best_layer - 1] = best_wire;
219  tdc_c[best_layer - 1] = best_tdc;
220  lat_c[best_layer - 1] = best_lat;
221  }
222  if (isSL1) {
223  outMetaPrimitives.emplace_back(metaPrimitive(
224  {mp.rawId, mp.t0, mp.x, mp.tanPhi, mp.phi, mp.phiB, mp.phi_cmssw, mp.phiB_cmssw, mp.chi2, new_quality,
225  mp.wi1, mp.tdc1, mp.lat1, mp.wi2, mp.tdc2, mp.lat2, mp.wi3, mp.tdc3, mp.lat3, mp.wi4,
226  mp.tdc4, mp.lat4, wi_c[0], tdc_c[0], lat_c[0], wi_c[1], tdc_c[1], lat_c[1], wi_c[2], tdc_c[2],
227  lat_c[2], wi_c[3], tdc_c[3], lat_c[3], -1}));
228  } else {
229  outMetaPrimitives.emplace_back(
230  metaPrimitive({mp.rawId, mp.t0, mp.x, mp.tanPhi, mp.phi, mp.phiB, mp.phi_cmssw,
231  mp.phiB_cmssw, mp.chi2, new_quality, wi_c[0], tdc_c[0], lat_c[0], wi_c[1],
232  tdc_c[1], lat_c[1], wi_c[2], tdc_c[2], lat_c[2], wi_c[3], tdc_c[3],
233  lat_c[3], mp.wi5, mp.tdc5, mp.lat5, mp.wi6, mp.tdc6, mp.lat6,
234  mp.wi7, mp.tdc7, mp.lat7, mp.wi8, mp.tdc8, mp.lat8, -1}));
235  }
236  } //SL2
237 }
int station() const
Return the station number.
Definition: DTChamberId.h:45
int superLayer() const
Return the superlayer number.
constexpr int CELL_SEMILENGTH
Definition: constants.h:314
U second(std::pair< T, U > const &p)
constexpr int SL1_CELLS_OFFSET
Definition: constants.h:247
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::map< int, float > shiftinfo_
constexpr int TIME_TO_TDC_COUNTS
Definition: constants.h:235
constexpr int LHC_CLK_FREQ
Definition: constants.h:222
int layer() const
Return the layer number.
Definition: DTLayerId.h:45
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:42
int sector() const
Definition: DTChamberId.h:52

◆ finish()

void MuonPathConfirmator::finish ( )

Definition at line 80 of file MuonPathConfirmator.cc.

References debug_, and LogDebug.

Referenced by progressbar.ProgressBar::__next__().

80  {
81  if (debug_)
82  LogDebug("MuonPathConfirmator") << "MuonPathConfirmator: finish";
83 };
#define LogDebug(id)

◆ initialise()

void MuonPathConfirmator::initialise ( const edm::EventSetup iEventSetup)

Definition at line 75 of file MuonPathConfirmator.cc.

References debug_, and LogDebug.

75  {
76  if (debug_)
77  LogDebug("MuonPathConfirmator") << "MuonPathConfirmator::initialiase";
78 }
#define LogDebug(id)

◆ run()

void MuonPathConfirmator::run ( edm::Event iEvent,
const edm::EventSetup iEventSetup,
std::vector< cmsdt::metaPrimitive inMetaPrimitives,
edm::Handle< DTDigiCollection dtdigis,
std::vector< cmsdt::metaPrimitive > &  outMetaPrimitives 
)

Definition at line 54 of file MuonPathConfirmator.cc.

References analyze(), IntegrityClient_cfi::ChId, debug_, LogDebug, max_drift_tdc, maxdriftinfo_, DTChamberId::sector(), DTChamberId::station(), and DTChamberId::wheel().

58  {
59  if (debug_)
60  LogDebug("MuonPathConfirmator") << "MuonPathConfirmator: run";
61 
62  // fit per SL (need to allow for multiple outputs for a single mpath)
63  if (!inMetaPrimitives.empty()) {
64  int dum_sl_rawid = inMetaPrimitives[0].rawId;
65  DTSuperLayerId dumSlId(dum_sl_rawid);
66  DTChamberId ChId(dumSlId.wheel(), dumSlId.station(), dumSlId.sector());
67  max_drift_tdc = maxdriftinfo_[dumSlId.wheel() + 2][dumSlId.station() - 1][dumSlId.sector() - 1];
68  }
69 
70  for (auto &mp : inMetaPrimitives) {
71  analyze(mp, dtdigis, outMetaPrimitives);
72  }
73 }
void analyze(cmsdt::metaPrimitive mp, edm::Handle< DTDigiCollection > dtdigis, std::vector< cmsdt::metaPrimitive > &outMetaPrimitives)
int maxdriftinfo_[5][4][14]
#define LogDebug(id)

Member Data Documentation

◆ debug_

bool MuonPathConfirmator::debug_
private

◆ LYRANDAHALF_RES_SHR

int MuonPathConfirmator::LYRANDAHALF_RES_SHR = 4
private

Definition at line 67 of file MuonPathConfirmator.h.

Referenced by analyze().

◆ max_drift_tdc

int MuonPathConfirmator::max_drift_tdc = -1
private

Definition at line 62 of file MuonPathConfirmator.h.

Referenced by analyze(), and run().

◆ maxdrift_filename_

edm::FileInPath MuonPathConfirmator::maxdrift_filename_
private

Definition at line 60 of file MuonPathConfirmator.h.

Referenced by MuonPathConfirmator().

◆ maxdriftinfo_

int MuonPathConfirmator::maxdriftinfo_[5][4][14]
private

Definition at line 61 of file MuonPathConfirmator.h.

Referenced by MuonPathConfirmator(), and run().

◆ minx_match_2digis_

double MuonPathConfirmator::minx_match_2digis_
private

Definition at line 56 of file MuonPathConfirmator.h.

Referenced by analyze().

◆ PARTIALS_PRECISSION

int MuonPathConfirmator::PARTIALS_PRECISSION = 4
private

Definition at line 64 of file MuonPathConfirmator.h.

Referenced by analyze().

◆ SEMICHAMBER_H

int MuonPathConfirmator::SEMICHAMBER_H = int(SEMICHAMBER_H_REAL)
private

Definition at line 69 of file MuonPathConfirmator.h.

Referenced by analyze().

◆ SEMICHAMBER_H_PRECISSION

int MuonPathConfirmator::SEMICHAMBER_H_PRECISSION = 13 + PARTIALS_PRECISSION
private

Definition at line 65 of file MuonPathConfirmator.h.

◆ SEMICHAMBER_H_REAL

float MuonPathConfirmator::SEMICHAMBER_H_REAL = ((235. / 2.) / (16. * 6.5)) * std::pow(2, SEMICHAMBER_H_PRECISSION)
private

Definition at line 68 of file MuonPathConfirmator.h.

◆ SEMICHAMBER_RES_SHR

int MuonPathConfirmator::SEMICHAMBER_RES_SHR = SEMICHAMBER_H_PRECISSION
private

Definition at line 66 of file MuonPathConfirmator.h.

Referenced by analyze().

◆ shift_filename_

edm::FileInPath MuonPathConfirmator::shift_filename_
private

Definition at line 58 of file MuonPathConfirmator.h.

Referenced by MuonPathConfirmator().

◆ shiftinfo_

std::map<int, float> MuonPathConfirmator::shiftinfo_
private

Definition at line 59 of file MuonPathConfirmator.h.

Referenced by analyze(), and MuonPathConfirmator().