test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuonSeedPtExtractor.cc
Go to the documentation of this file.
5 
6 #include "TMath.h"
7 #include <sstream>
8 
10 : theBeamSpot(0.,0.,0.),
11  scaleDT_( par.getParameter<bool>("scaleDT") )
12 {
13  init(par);
14 }
15 
16 
18 {
19  // load pT seed parameters
20  // DT combinations
21  fillParametersForCombo("DT_12", par);
22  fillParametersForCombo("DT_13", par);
23  fillParametersForCombo("DT_14", par);
24  fillParametersForCombo("DT_23", par);
25  fillParametersForCombo("DT_24", par);
26  fillParametersForCombo("DT_34", par);
27  // CSC combinations
28  fillParametersForCombo("CSC_01", par);
29  fillParametersForCombo("CSC_12", par);
30  fillParametersForCombo("CSC_02", par);
31  fillParametersForCombo("CSC_13", par);
32  fillParametersForCombo("CSC_03", par);
33  fillParametersForCombo("CSC_14", par);
34  fillParametersForCombo("CSC_23", par);
35  fillParametersForCombo("CSC_24", par);
36  fillParametersForCombo("CSC_34", par);
37 
38  // Overlap combinations
39  fillParametersForCombo("OL_1213", par);
40  fillParametersForCombo("OL_1222", par);
41  fillParametersForCombo("OL_1232", par);
42  fillParametersForCombo("OL_2213", par);
43  fillParametersForCombo("OL_2222", par);
44 
45  // Single segments (CSC)
46  fillParametersForCombo("SME_11", par);
47  fillParametersForCombo("SME_12", par);
48  fillParametersForCombo("SME_13", par);
49  fillParametersForCombo("SME_21", par);
50  fillParametersForCombo("SME_22", par);
51  fillParametersForCombo("SME_31", par);
52  fillParametersForCombo("SME_32", par);
53  fillParametersForCombo("SME_41", par);
54  fillParametersForCombo("SME_42", par);
55 
56  // Single segments (DT)
57  fillParametersForCombo("SMB_10", par);
58  fillParametersForCombo("SMB_11", par);
59  fillParametersForCombo("SMB_12", par);
60  fillParametersForCombo("SMB_20", par);
61  fillParametersForCombo("SMB_21", par);
62  fillParametersForCombo("SMB_22", par);
63  fillParametersForCombo("SMB_30", par);
64  fillParametersForCombo("SMB_31", par);
65  fillParametersForCombo("SMB_32", par);
66 
67  fillScalesForCombo("CSC_01_1_scale", par);
68  fillScalesForCombo("CSC_12_1_scale", par);
69  fillScalesForCombo("CSC_12_2_scale", par);
70  fillScalesForCombo("CSC_12_3_scale", par);
71  fillScalesForCombo("CSC_13_2_scale", par);
72  fillScalesForCombo("CSC_13_3_scale", par);
73  fillScalesForCombo("CSC_14_3_scale", par);
74  fillScalesForCombo("CSC_23_1_scale", par);
75  fillScalesForCombo("CSC_23_2_scale", par);
76  fillScalesForCombo("CSC_24_1_scale", par);
77  fillScalesForCombo("CSC_34_1_scale", par);
78  fillScalesForCombo("DT_12_1_scale", par);
79  fillScalesForCombo("DT_12_2_scale", par);
80  fillScalesForCombo("DT_13_1_scale", par);
81  fillScalesForCombo("DT_13_2_scale", par);
82  fillScalesForCombo("DT_14_1_scale", par);
83  fillScalesForCombo("DT_14_2_scale", par);
84  fillScalesForCombo("DT_23_1_scale", par);
85  fillScalesForCombo("DT_23_2_scale", par);
86  fillScalesForCombo("DT_24_1_scale", par);
87  fillScalesForCombo("DT_24_2_scale", par);
88  fillScalesForCombo("DT_34_1_scale", par);
89  fillScalesForCombo("DT_34_2_scale", par);
90  fillScalesForCombo("OL_1213_0_scale", par);
91  fillScalesForCombo("OL_1222_0_scale", par);
92  fillScalesForCombo("OL_1232_0_scale", par);
93  fillScalesForCombo("OL_2213_0_scale", par);
94  fillScalesForCombo("OL_2222_0_scale", par);
95  fillScalesForCombo("SMB_10_0_scale", par);
96  fillScalesForCombo("SMB_11_0_scale", par);
97  fillScalesForCombo("SMB_12_0_scale", par);
98  fillScalesForCombo("SMB_20_0_scale", par);
99  fillScalesForCombo("SMB_21_0_scale", par);
100  fillScalesForCombo("SMB_22_0_scale", par);
101  fillScalesForCombo("SMB_30_0_scale", par);
102  fillScalesForCombo("SMB_31_0_scale", par);
103  fillScalesForCombo("SMB_32_0_scale", par);
104  fillScalesForCombo("SME_11_0_scale", par);
105  fillScalesForCombo("SME_12_0_scale", par);
106  fillScalesForCombo("SME_13_0_scale", par);
107  fillScalesForCombo("SME_21_0_scale", par);
108  fillScalesForCombo("SME_22_0_scale", par);
109 
110 }
111 
112 
114 }
115 
116 
118 {
119  theParametersForCombo[name] = pset.getParameter<std::vector<double> >(name);
120 }
121 
122 
124 {
125  theScalesForCombo[name] = pset.getParameter<std::vector<double> >(name);
126 }
127 
128 
131 {
132  GlobalPoint innerPoint = firstHit->globalPosition() - theBeamSpot;
133  GlobalPoint outerPoint = secondHit->globalPosition() - theBeamSpot;
136 
137  // ways in which the hits could be in the wrong order
138  if( (outerHit->isDT() && innerHit->isCSC())
139  || (outerHit->isDT() && innerHit->isDT() && (innerPoint.perp() > outerPoint.perp()))
140  || (outerHit->isCSC() && innerHit->isCSC() && (fabs(innerPoint.z()) > fabs(outerPoint.z()))) )
141  {
142  innerHit = secondHit;
143  outerHit = firstHit;
144  innerPoint = innerHit->globalPosition() - theBeamSpot;
145  outerPoint = outerHit->globalPosition() - theBeamSpot;
146  }
147 
148  double phiInner = innerPoint.phi();
149  double phiOuter = outerPoint.phi();
150 
151  double etaInner = innerPoint.eta();
152  double etaOuter = outerPoint.eta();
153  //std::cout<<" inner pos = "<< innerPoint << " phi eta " << phiInner << " " << etaInner << std::endl;
154  //std::cout<<" outer pos = "<< outerPoint << " phi eta " << phiOuter << " " << etaOuter << std::endl;
155  //double thetaInner = firstHit->globalPosition().theta();
156  // if some of the segments is missing r-phi measurement then we should
157  // use only the 4D phi estimate (also use 4D eta estimate only)
158  // the direction is not so important (it will be corrected)
159  /*
160  bool firstOK = (4==allValidSets[iSet][firstMeasurement]->hit->dimension());
161  bool lastOK = (4==allValidSets[iSet][lastMeasurement]->hit->dimension());
162  if(!(firstOK * lastOK)){
163  if(!firstOK){
164  }
165  if(!firstOK){
166  }
167  }
168  */
169  double dPhi = phiInner - phiOuter;
170  if(dPhi < -TMath::Pi()){
171  dPhi += 2*TMath::Pi();
172  }
173  else if(dPhi > TMath::Pi()){
174  dPhi -= 2*TMath::Pi();
175  }
176  int sign = 1;
177  if ( dPhi< 0.) {
178  dPhi = -dPhi;
179  sign = -1;
180  }
181 
182  if (dPhi < 1.0e-6){
183  dPhi = 1.0e-6;
184  }
185  double eta = fabs(etaOuter);// what if it is 2D segment? use inner?
186 
187 
188  std::vector <int> stationCoded(2);
189  stationCoded[0] = stationCoded[1] = 999;
190 
191  DetId detId_inner = innerHit->hit()->geographicalId();
192  DetId detId_outer = outerHit->hit()->geographicalId();
193 
194  stationCoded[0] = stationCode(innerHit);
195  stationCoded[1] = stationCode(outerHit);
196 
197  std::ostringstream os0 ;
198  std::ostringstream os1 ;
199  os0 << abs(stationCoded[0]);
200  os1 << abs(stationCoded[1]);
201 
202  //std::cout<<" st1 = "<<stationCoded[0]<<" st2 = "<<stationCoded[1]<<std::endl;
203  //std::cout<<" detId_inner = "<<detId_inner.rawId()<<" detId_outer = "<< detId_outer.rawId()<<std::endl;
204  std::string combination = "0";
205  std::string init_combination = combination;
206  bool singleSegment = false;
207  //if(detId_first == detId_second){
208  if( stationCoded[0] == stationCoded[1]){
209  // single segment - DT or CSC
210  singleSegment = true;
211  //eta = innerPoint.eta();
212  GlobalVector gv = innerHit->globalDirection();
213  double gvPerpPos = gv.x()*gv.x() + gv.y()*gv.y();
214  if (gvPerpPos < 1e-32) gvPerpPos = 1e-32; gvPerpPos=sqrt(gvPerpPos);
215  // Psi is angle between the segment origin and segment direction
216  // Use dot product between two vectors to get Psi in global x-y plane
217  double cosDpsi = (gv.x()*innerPoint.x() + gv.y()*innerPoint.y());
218  if (cosDpsi!=0){
219  cosDpsi /= sqrt(innerPoint.x()*innerPoint.x() + innerPoint.y()*innerPoint.y());
220  cosDpsi /= gvPerpPos;
221  cosDpsi = cosDpsi > 1 ? 1 : cosDpsi;
222  cosDpsi = cosDpsi < -1 ? -1 : cosDpsi;
223  }
224 
225  double axb = ( innerPoint.x()*gv.y() ) - ( innerPoint.y()*gv.x() ) ;
226  sign = (axb < 0.) ? 1 : -1;
227 
228  double dpsi = fabs(acos(cosDpsi)) ;
229  if ( dpsi > TMath::Pi()/2.) {
230  dpsi = TMath::Pi() - dpsi;
231  }
232  if (fabs(dpsi) < 0.00005) {
233  dpsi = 0.00005;
234  }
235  dPhi = dpsi;
236 
237  if(innerHit->isDT())
238  {
239  DTChamberId dtCh(detId_inner);
240  std::ostringstream os;
241  os << dtCh.station() << abs(dtCh.wheel());
242  combination = "SMB_"+os.str();
243  }
244  else if(innerHit->isCSC())
245  {
246  CSCDetId cscId(detId_inner);
247  std::ostringstream os;
248  int ring = cscId.ring();
249  if(ring == 4) ring = 1;
250  os << cscId.station() << ring;
251  combination = "SME_"+os.str();
252  }
253  else
254  {
255  throw cms::Exception("MuonSeedPtExtractor") << "Bad hit DetId";
256  }
257  }
258  else{
259  if(stationCoded[0]<0){
260  if(stationCoded[1]<0){
261  // DT-DT
262  combination = "DT_" + os0.str() + os1.str();
263  }
264  else{
265  // DT-CSC
266  eta = fabs(etaInner);
267  if(-1==stationCoded[0]){
268  switch (stationCoded[1]){
269  case 1:
270  combination = "OL_1213";
271  break;
272  case 2:
273  combination = "OL_1222";
274  break;
275  case 3:
276  combination = "OL_1232";
277  break;
278  default:
279  // can not be
280  break;
281  }
282  }
283  else if (-2==stationCoded[0]){
284  if(1==stationCoded[1]){
285  combination = "OL_2213";
286  }
287  else{
288  // can not be (not coded?)
289  combination = "OL_2222";// in case
290  }
291  }
292  else{
293  // can not be
294  }
295  }
296  }
297  else{
298  if(stationCoded[1]<0){
299  // CSC-DT
300  // can not be
301  }
302  else{
303  // CSC-CSC
304  combination = "CSC_" + os0.str() + os1.str();
305  if("CSC_04" == combination){
306  combination = "CSC_14";
307  }
308  }
309  }
310  }
311 
312  std::vector<double> pTestimate(2);//
313  //std::cout<<" combination = "<<combination<<std::endl;
314  if(init_combination!=combination){
315  //std::cout<<" combination = "<<combination<<" eta = "<<eta<<" dPhi = "<<dPhi<<std::endl;
316  ParametersMap::const_iterator parametersItr = theParametersForCombo.find(combination);
317  if(parametersItr == theParametersForCombo.end()) {
318  // edm::LogWarning("RecoMuon|MuonSeedGenerator|MuonSeedPtExtractor") << "Cannot find parameters for combo " << combination;
319  edm::LogWarning("BadSegmentCombination") << "Cannot find parameters for combo " << combination;
320  pTestimate[0] = pTestimate[1] = 100;
321  // throw cms::Exception("MuonSeedPtEstimator") << "Cannot find parameters for combo " << combination;
322  } else {
323 
324  if(scaleDT_ && outerHit->isDT() )
325  {
326  pTestimate = getPt(parametersItr->second, eta, dPhi, combination, detId_outer);
327  }
328  else
329  {
330  pTestimate = getPt(parametersItr->second, eta, dPhi);
331  }
332 
333  if(singleSegment){
334  pTestimate[0] = fabs(pTestimate[0]);
335  pTestimate[1] = fabs(pTestimate[1]);
336  }
337  pTestimate[0] *= double(sign);
338  }
339  }
340  else{
341  // often a MB3 - ME1/3 seed
342  pTestimate[0] = pTestimate[1] = 100;
343  // hmm
344  }
345  // std::cout<<" pTestimate[0] = "<<pTestimate[0]<<" pTestimate[1] = "<<pTestimate[1]<<std::endl;
346  /*
347  MuonRecHitContainer_clusters[cluster][iHit+1]->isDT());
348  if(specialCase){
349  DTChamberId dtCh(detId);
350  DTChamberId dtCh_2(detId_2);
351  specialCase = (dtCh.station() == dtCh_2.station());
352  }
353  */
354  //return vPara;
355  return pTestimate;
356 }
357 
358 
360 {
361  DetId detId(hit->hit()->geographicalId());
362  int result = -999;
363  if(hit->isDT() ){
364  DTChamberId dtCh(detId);
365  //std::cout<<"first (DT) St/W/S = "<<dtCh.station()<<"/"<<dtCh.wheel()<<"/"<<dtCh.sector()<<"/"<<std::endl;
366  result = -1 * dtCh.station();
367  }
368  else if( hit->isCSC() ){
369  CSCDetId cscID(detId);
370  //std::cout<<"first (CSC) E/S/R/C = "<<cscID.endcap()<<"/"<<cscID.station()<<"/"<<cscID.ring()<<"/"<<cscID.chamber()<<std::endl;
371  result = cscID.station();
372  if(result == 1 && (1 == cscID.ring() || 4 == cscID.ring()) )
373  result = 0;
374  }
375  else if(hit->isRPC()){
376  }
377  return result;
378 }
379 
380 
381 std::vector<double> MuonSeedPtExtractor::getPt(const std::vector<double> & vPara, double eta, double dPhi ) const {
382  //std::cout<<" eta = "<<eta<<" dPhi = "<<dPhi<<" vPara[0] = "<<vPara[0]<<" vPara[1] = "<<vPara[1]<<" vPara[2] = "<<vPara[2]<<std::endl;
383  double h = fabs(eta);
384  double estPt = ( vPara[0] + vPara[1]*h + vPara[2]*h*h ) / dPhi;
385  double estSPt = ( vPara[3] + vPara[4]*h + vPara[5]*h*h ) / dPhi;
386  // std::cout<<"estPt = "<<estPt<<std::endl;
387  std::vector<double> paraPt ;
388  paraPt.push_back( estPt );
389  paraPt.push_back( estSPt ) ;
390  return paraPt ;
391 }
392 
393 // combined pt parameterization and pt scale
394 std::vector<double> MuonSeedPtExtractor::getPt(const std::vector<double> & vPara, double eta, double dPhi,
395  const std::string & combination, const DTChamberId & outerDetId ) const {
396  double h = fabs(eta);
397  double estPt = ( vPara[0] + vPara[1]*h + vPara[2]*h*h ) / dPhi;
398  // changed by S.C.
399  double estSPt = ( vPara[3] + vPara[4]*h + vPara[5]*h*h ) / dPhi;
400 
401  // scale the pt and spt , changed by S.C.
402  int wheel = 0;
403  if(combination[0] == 'D') {
404  wheel = abs(outerDetId.wheel());
405  }
406 
407  std::ostringstream os;
408  os << combination << "_" << wheel << "_scale";
409 
410  ScalesMap::const_iterator scalesItr = theScalesForCombo.find(os.str());
411  if ( scalesItr != theScalesForCombo.end()) {
412  double t1 = scalesItr->second[0];
413  double scaleFactor = 1./( 1. + t1/( estPt + 10. ) ) ;
414 
415  estSPt = estSPt * scaleFactor ;
416  estPt = estPt * scaleFactor ;
417  }
418 
419  // std::cout<<"estPt = "<<estPt<<std::endl;
420  std::vector<double> paraPt ;
421  paraPt.push_back( estPt );
422  paraPt.push_back( estSPt ) ;
423  return paraPt ;
424 }
425 
virtual ~MuonSeedPtExtractor()
Destructor.
const double Pi
T getParameter(std::string const &) const
T perp() const
Definition: PV3DBase.h:72
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
void fillParametersForCombo(const std::string &name, const edm::ParameterSet &pset)
double sign(double x)
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
void init(const edm::ParameterSet &par)
void fillScalesForCombo(const std::string &name, const edm::ParameterSet &pset)
tuple result
Definition: mps_fire.py:83
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
int stationCode(MuonTransientTrackingRecHit::ConstMuonRecHitPointer hit) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int ring() const
Definition: CSCDetId.h:75
ParametersMap theParametersForCombo
Definition: DetId.h:18
T eta() const
Definition: PV3DBase.h:76
MuonSeedPtExtractor(const edm::ParameterSet &)
Constructor with Parameter set and MuonServiceProxy.
int station() const
Definition: CSCDetId.h:86
std::shared_ptr< MuonTransientTrackingRecHit const > ConstMuonRecHitPointer
virtual std::vector< double > pT_extract(MuonTransientTrackingRecHit::ConstMuonRecHitPointer firstHit, MuonTransientTrackingRecHit::ConstMuonRecHitPointer secondHit) const
int station() const
Return the station number.
Definition: DTChamberId.h:51
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
T x() const
Definition: PV3DBase.h:62
std::vector< double > getPt(const std::vector< double > &vPara, double eta, double dPhi) const