CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SETSeedFinder.cc
Go to the documentation of this file.
1 
15 
16 #include "TMath.h"
17 
18 
19 using namespace edm;
20 using namespace std;
21 
22 const string metname = "Muon|RecoMuon|SETMuonSeedFinder";
23 
26 {
27  // Parameter set for the Builder
28  ParameterSet trajectoryBuilderParameters = parameterSet.getParameter<ParameterSet>("SETTrajBuilderParameters");
29  apply_prePruning = trajectoryBuilderParameters.getParameter<bool>("Apply_prePruning");
30  // load pT seed parameters
31  thePtExtractor = new MuonSeedPtExtractor(trajectoryBuilderParameters);
32 
33 }
34 
35 
37  std::vector<TrajectorySeed> & result)
38 {
39 }
40 
41 
42 
43 // there is an existing sorter somewhere in the CMSSW code (I think) - delete that
44 struct sorter{
47  return (hit_1->globalPosition().mag2()<hit_2->globalPosition().mag2());
48 
49  }
50 } sortSegRadius;// smaller first
51 
52 
53 
54 std::vector<SETSeedFinder::MuonRecHitContainer>
56 {
57  stable_sort(cluster.begin(), cluster.end(),sortSegRadius);
58  //---- group hits in detector layers (if in same layer); the idea is that
59  //---- some hits could not belong to a track simultaneously - these will be in a
60  //---- group; two hits from one and the same group will not go to the same track
61  std::vector< MuonRecHitContainer > MuonRecHitContainer_perLayer;
62  if(cluster.size()){
63  int iHit =0;
64  MuonRecHitContainer hitsInThisLayer;
65  hitsInThisLayer.push_back(cluster[iHit]);
66  DetId detId = cluster[iHit]->hit()->geographicalId();
67  const GeomDet* geomDet = theService->trackingGeometry()->idToDet( detId );
68  while(iHit<int(cluster.size())-1){
69  DetId detId_2 = cluster[iHit+1]->hit()->geographicalId();
70  const GlobalPoint gp_nextHit = cluster[iHit+1]->globalPosition();
71 
72  // this is the distance of the "second" hit to the "first" detector (containing the "first hit")
73  float distanceToDetector = fabs(geomDet->surface().localZ(gp_nextHit));
74 
75  //---- hits from DT and CSC could be very close in angle but incosistent with
76  //---- belonging to a common track (and these are different surfaces);
77  //---- also DT (and CSC now - 090822) hits from a station (in a pre-cluster) should be always in a group together;
78  //---- take this into account and put such hits in a group together
79 
80  bool specialCase = false;
81  if( detId.subdetId() == MuonSubdetId::DT &&
82  detId_2.subdetId() == MuonSubdetId::DT ){
83  DTChamberId dtCh(detId);
84  DTChamberId dtCh_2(detId_2);
85  specialCase = (dtCh.station() == dtCh_2.station());
86  }
87  else if(detId.subdetId() == MuonSubdetId::CSC &&
88  detId_2.subdetId() == MuonSubdetId::CSC){
89  CSCDetId cscCh(detId);
90  CSCDetId cscCh_2(detId_2);
91  specialCase = (cscCh.station() == cscCh_2.station() && cscCh.ring() == cscCh_2.ring());
92  }
93 
94  if(distanceToDetector<0.001 || true==specialCase){ // hardcoded value - remove!
95  hitsInThisLayer.push_back(cluster[iHit+1]);
96  }
97  else{
98  specialCase = false;
99  if(( (cluster[iHit]->isDT() &&
100  cluster[iHit+1]->isCSC()) ||
101  (cluster[iHit]->isCSC() &&
102  cluster[iHit+1]->isDT())) &&
103  //---- what is the minimal distance between a DT and a CSC hit belonging
104  //---- to a common track? (well, with "reasonable" errors; put 10 cm for now)
105  fabs(cluster[iHit+1]->globalPosition().mag() -
106  cluster[iHit]->globalPosition().mag())<10.){
107  hitsInThisLayer.push_back(cluster[iHit+1]);
108  // change to Stoyan - now we also update the detID here... give it a try. IBL 080905
109  detId = cluster[iHit+1]->hit()->geographicalId();
110  geomDet = theService->trackingGeometry()->idToDet( detId );
111  }
112  else if(!specialCase){
113  //---- put the group of hits in the vector (containing the groups of hits)
114  //---- and continue with next layer (group)
115  MuonRecHitContainer_perLayer.push_back(hitsInThisLayer);
116  hitsInThisLayer.clear();
117  hitsInThisLayer.push_back(cluster[iHit+1]);
118  detId = cluster[iHit+1]->hit()->geographicalId();
119  geomDet = theService->trackingGeometry()->idToDet( detId );
120  }
121  }
122  ++iHit;
123  }
124  MuonRecHitContainer_perLayer.push_back(hitsInThisLayer);
125  }
126  return MuonRecHitContainer_perLayer;
127 }
128 //
129 void SETSeedFinder::limitCombinatorics(std::vector< MuonRecHitContainer > & MuonRecHitContainer_perLayer){
130  const int maximumNumberOfCombinations = 1000000;
131  unsigned nLayers = MuonRecHitContainer_perLayer.size();
132  if(1==nLayers){
133  return ;
134  }
135  // maximal number of (segment) layers would be upto ~12; see next function
136  // below is just a quick fix for a rare "overflow"
137  if(MuonRecHitContainer_perLayer.size()>15){
138  MuonRecHitContainer_perLayer.resize(1);
139  return;
140  }
141 
142  std::vector <double> sizeOfLayer(nLayers);
143  //std::cout<<" nLayers = "<<nLayers<<std::endl;
144  double nAllCombinations = 1.;
145  for(unsigned int i = 0;i<nLayers;++i ){
146  //std::cout<<" i = "<<i<<" size = "<<MuonRecHitContainer_perLayer.at(i).size()<<std::endl;
147  sizeOfLayer.at(i) = MuonRecHitContainer_perLayer.at(i).size();
148  nAllCombinations*=MuonRecHitContainer_perLayer.at(i).size();
149  }
150  //std::cout<<"nAllCombinations = "<<nAllCombinations<<std::endl;
151  //---- Erase most busy detector layers until we get less than maximumNumberOfCombinations combinations
152  int iCycle = 0;
153  while(nAllCombinations > float(maximumNumberOfCombinations)){
154  ++iCycle;
155  std::vector <double>::iterator maxEl_it = max_element(sizeOfLayer.begin(),sizeOfLayer.end());
156  int maxEl = maxEl_it - sizeOfLayer.begin();
157  nAllCombinations/=MuonRecHitContainer_perLayer.at(maxEl).size();
158  //std::cout<<" iCycle = "<<iCycle<<" nAllCombinations = "<<nAllCombinations<<std::endl;
159  MuonRecHitContainer_perLayer.erase(MuonRecHitContainer_perLayer.begin()+maxEl);
160  sizeOfLayer.erase(sizeOfLayer.begin()+maxEl);
161  }
162  return;
163 }
164 //
165 std::vector<SETSeedFinder::MuonRecHitContainer>
166 SETSeedFinder::findAllValidSets(const std::vector<SETSeedFinder::MuonRecHitContainer> & MuonRecHitContainer_perLayer)
167 {
168  std::vector <MuonRecHitContainer> allValidSets;
169  // build all possible combinations (i.e valid sets; the algorithm name is after this feature -
170  // SET algorithm)
171  //
172  // ugly... use recursive function?!
173  // or implement Ingo's suggestion (a la ST)
174  unsigned nLayers = MuonRecHitContainer_perLayer.size();
175  if(1==nLayers){
176  return allValidSets;
177  }
178  MuonRecHitContainer validSet;
179  unsigned int iPos0 = 0;
180  std::vector <unsigned int> iLayer(12);// could there be more than 11 layers?
181  std::vector <unsigned int> size(12);
182  if(iPos0<nLayers){
183  size.at(iPos0) = MuonRecHitContainer_perLayer.at(iPos0).size();
184  for(iLayer[iPos0] = 0; iLayer[iPos0]<size[iPos0];++iLayer[iPos0]){
185  validSet.clear();
186  validSet.push_back(MuonRecHitContainer_perLayer[iPos0][iLayer[iPos0]]);
187  unsigned int iPos1 = 1;
188  if(iPos1<nLayers){
189  size.at(iPos1) = MuonRecHitContainer_perLayer.at(iPos1).size();
190  for(iLayer[iPos1] = 0; iLayer[iPos1]<size[iPos1];++iLayer[iPos1]){
191  validSet.resize(iPos1);
192  validSet.push_back(MuonRecHitContainer_perLayer[iPos1][iLayer[iPos1]]);
193  unsigned int iPos2 = 2;
194  if(iPos2<nLayers){
195  size.at(iPos2) = MuonRecHitContainer_perLayer.at(iPos2).size();
196  for(iLayer[iPos2] = 0; iLayer[iPos2]<size[iPos2];++iLayer[iPos2]){
197  validSet.resize(iPos2);
198  validSet.push_back(MuonRecHitContainer_perLayer[iPos2][iLayer[iPos2]]);
199  unsigned int iPos3 = 3;
200  if(iPos3<nLayers){
201  size.at(iPos3) = MuonRecHitContainer_perLayer.at(iPos3).size();
202  for(iLayer[iPos3] = 0; iLayer[iPos3]<size[iPos3];++iLayer[iPos3]){
203  validSet.resize(iPos3);
204  validSet.push_back(MuonRecHitContainer_perLayer[iPos3][iLayer[iPos3]]);
205  unsigned int iPos4 = 4;
206  if(iPos4<nLayers){
207  size.at(iPos4) = MuonRecHitContainer_perLayer.at(iPos4).size();
208  for(iLayer[iPos4] = 0; iLayer[iPos4]<size[iPos4];++iLayer[iPos4]){
209  validSet.resize(iPos4);
210  validSet.push_back(MuonRecHitContainer_perLayer[iPos4][iLayer[iPos4]]);
211  unsigned int iPos5 = 5;
212  if(iPos5<nLayers){
213  size.at(iPos5) = MuonRecHitContainer_perLayer.at(iPos5).size();
214  for(iLayer[iPos5] = 0; iLayer[iPos5]<size[iPos5];++iLayer[iPos5]){
215  validSet.resize(iPos5);
216  validSet.push_back(MuonRecHitContainer_perLayer[iPos5][iLayer[iPos5]]);
217  unsigned int iPos6 = 6;
218  if(iPos6<nLayers){
219  size.at(iPos6) = MuonRecHitContainer_perLayer.at(iPos6).size();
220  for(iLayer[iPos6] = 0; iLayer[iPos6]<size[iPos6];++iLayer[iPos6]){
221  validSet.resize(iPos6);
222  validSet.push_back(MuonRecHitContainer_perLayer[iPos6][iLayer[iPos6]]);
223  unsigned int iPos7 = 7;
224  if(iPos7<nLayers){
225  size.at(iPos7) = MuonRecHitContainer_perLayer.at(iPos7).size();
226  for(iLayer[iPos7] = 0; iLayer[iPos7]<size[iPos7];++iLayer[iPos7]){
227  validSet.resize(iPos7);
228  validSet.push_back(MuonRecHitContainer_perLayer[iPos7][iLayer[iPos7]]);
229  unsigned int iPos8 = 8;
230  if(iPos8<nLayers){
231  size.at(iPos8) = MuonRecHitContainer_perLayer.at(iPos8).size();
232  for(iLayer[iPos8] = 0; iLayer[iPos8]<size[iPos8];++iLayer[iPos8]){
233  validSet.resize(iPos8);
234  validSet.push_back(MuonRecHitContainer_perLayer[iPos8][iLayer[iPos8]]);
235  unsigned int iPos9 = 9;
236  if(iPos9<nLayers){
237  size.at(iPos9) = MuonRecHitContainer_perLayer.at(iPos9).size();
238  for(iLayer[iPos9] = 0; iLayer[iPos9]<size[iPos9];++iLayer[iPos9]){
239  validSet.resize(iPos9);
240  validSet.push_back(MuonRecHitContainer_perLayer[iPos9][iLayer[iPos9]]);
241  unsigned int iPos10 = 10;
242  if(iPos10<nLayers){
243  size.at(iPos10) = MuonRecHitContainer_perLayer.at(iPos10).size();
244  for(iLayer[iPos10] = 0; iLayer[iPos10]<size[iPos10];++iLayer[iPos10]){
245  validSet.resize(iPos10);
246  validSet.push_back(MuonRecHitContainer_perLayer[iPos10][iLayer[iPos10]]);
247  unsigned int iPos11 = 11;// more?
248  if(iPos11<nLayers){
249  size.at(iPos11) = MuonRecHitContainer_perLayer.at(iPos11).size();
250  for(iLayer[iPos11] = 0; iLayer[iPos11]<size[iPos11];++iLayer[iPos11]){
251  }
252  }
253  else{
254  allValidSets.push_back(validSet);
255 
256  }
257  }
258  }
259  else{
260  allValidSets.push_back(validSet);
261  }
262  }
263  }
264  else{
265  allValidSets.push_back(validSet);
266  }
267  }
268  }
269  else{
270  allValidSets.push_back(validSet);
271  }
272  }
273  }
274  else{
275  allValidSets.push_back(validSet);
276  }
277  }
278  }
279  else{
280  allValidSets.push_back(validSet);
281  }
282  }
283  }
284  else{
285  allValidSets.push_back(validSet);
286  }
287  }
288  }
289  else{
290  allValidSets.push_back(validSet);
291  }
292  }
293  }
294  else{
295  allValidSets.push_back(validSet);
296  }
297  }
298  }
299  else{
300  allValidSets.push_back(validSet);
301  }
302  }
303  }
304  else{
305  allValidSets.push_back(validSet);
306  }
307  }
308  }
309  else{
310  allValidSets.push_back(validSet);
311  }
312  return allValidSets;
313 }
314 
315 
316 std::pair <int, int> // or <bool, bool>
317 SETSeedFinder::checkAngleDeviation(double dPhi_1, double dPhi_2) const
318 {
319  // Two conditions:
320  // a) deviations should be only to one side (above some absolute value cut to avoid
321  // material effects; this should be refined)
322  // b) deviatiation in preceding steps should be bigger due to higher magnetic field
323  // (again - a minimal value cut should be in place; this also should account for
324  // the small (Z) distances in overlaping CSC chambers)
325 
326  double mult = dPhi_1 * dPhi_2;
327  int signVal = 1;
328  if(fabs(dPhi_1)<fabs(dPhi_2)){
329  signVal = -1;
330  }
331  int signMult = -1;
332  if(mult>0) signMult = 1;
333  std::pair <int, int> sign;
334  sign = make_pair (signVal, signMult);
335 
336  return sign;
337 }
338 
339 
340 void SETSeedFinder::validSetsPrePruning(std::vector<SETSeedFinder::MuonRecHitContainer> & allValidSets)
341 {
342  //---- this actually is a pre-pruning; it does not include any fit information;
343  //---- it is intended to remove only very "wild" segments from a set;
344  //---- no "good" segment is to be lost (otherwise - widen the parameters)
345 
346  for(unsigned int iSet = 0;iSet<allValidSets.size();++iSet){
347  pre_prune(allValidSets[iSet]);
348  }
349 }
350 
351 
353 {
354  unsigned nHits = validSet.size();
355  if(nHits>3){ // to decide we need at least 4 measurements
356  // any information could be used to make a decision for pruning
357  // maybe dPhi (delta Phi) is enough
358  std::vector <double> dPhi;
359  double dPhi_tmp;
360  bool wildCandidate;
361  int pruneHit_tmp;
362 
363  for(unsigned int iHit = 1;iHit<nHits;++iHit){
364  dPhi_tmp = validSet[iHit]->globalPosition().phi() -
365  validSet[iHit-1]->globalPosition().phi();
366  dPhi.push_back(dPhi_tmp);
367  }
368  std::vector <int> pruneHit;
369  //---- loop over all the hits in a set
370 
371  for(unsigned int iHit = 0;iHit<nHits;++iHit){
372  double dPHI_MIN = 0.02;//?? hardcoded - remove it
373  if(iHit){
374  // if we have to remove the very first hit (iHit == 0) then
375  // we'll probably be in trouble already
376  wildCandidate = false;
377  // actually 2D is bad only if not r-phi... Should I refine it?
378  // a hit is a candidate for pruning only if dPhi > dPHI_MIN;
379  // pruning decision is based on combination of hits characteristics
380  if(4==validSet[iHit-1]->dimension() && 4 == validSet[iHit]->dimension() &&
381  fabs(validSet[iHit]->globalPosition().phi() -
382  validSet[iHit-1]->globalPosition().phi())>dPHI_MIN ){
383  wildCandidate = true;
384  }
385  pruneHit_tmp = -1;
386  if(wildCandidate){
387  // OK - this couple doesn't look good (and is from 4D segments); proceed...
388  if(1==iHit){// the first and the last hits are special case
389  if(4==validSet[iHit+1]->dimension() && 4 == validSet[iHit+2]->dimension()){//4D?
390  // is the picture better if we remove the second hit?
391  dPhi_tmp = validSet[iHit+1]->globalPosition().phi() -
392  validSet[iHit-1]->globalPosition().phi();
393  // is the deviation what we expect (sign, not magnitude)?
394  std::pair <int, int> sign = checkAngleDeviation(dPhi_tmp, dPhi[2]);
395  if( 1==sign.first && 1==sign.second){
396  pruneHit_tmp = iHit;// mark the hit 1 for removing
397  }
398  }
399  }
400  else if(iHit>1 && iHit<validSet.size()-1){
401  if(4 == validSet[0]->dimension() && // we rely on the first (most important) couple
402  4 == validSet[1]->dimension() &&
403  pruneHit.back()!=int(iHit-1) && pruneHit.back()!=1){// check if hits are already marked
404  // decide which of the two hits should be removed (if any; preferably the outer one i.e.
405  // iHit rather than iHit-1); here - check what we get by removing iHit
406  dPhi_tmp = validSet[iHit+1]->globalPosition().phi() -
407  validSet[iHit-1]->globalPosition().phi();
408  // first couple is most important anyway so again compare to it
409  std::pair <int, int> sign = checkAngleDeviation(dPhi[0],dPhi_tmp);
410  if( 1==sign.first && 1==sign.second){
411  pruneHit_tmp = iHit; // mark the hit iHit for removing
412  }
413  else{ // iHit is not to be removed; proceed...
414  // what if we remove (iHit - 1) instead of iHit?
415  dPhi_tmp = validSet[iHit+1]->globalPosition().phi() -
416  validSet[iHit]->globalPosition().phi();
417  std::pair <int, int> sign = checkAngleDeviation(dPhi[0],dPhi_tmp);
418  if( 1==sign.first && 1==sign.second){
419  pruneHit_tmp = iHit-1;// mark the hit (iHit -1) for removing
420  }
421  }
422  }
423  }
424  else{
425  // the last hit: if picture is not good - remove it
426  if(pruneHit.size()>1 && pruneHit[pruneHit.size()-1]<0 && pruneHit[pruneHit.size()-2]<0){
427  std::pair <int, int> sign = checkAngleDeviation(dPhi[dPhi.size()-2], dPhi[dPhi.size()-1]);
428  if( -1==sign.first && -1==sign.second){// here logic is a bit twisted
429  pruneHit_tmp = iHit; // mark the last hit for removing
430  }
431  }
432  }
433  }
434  pruneHit.push_back(pruneHit_tmp);
435  }
436  }
437  // }
438  // actual pruning
439  for(unsigned int iHit = 1;iHit<nHits;++iHit){
440  int count = 0;
441  if(pruneHit[iHit-1]>0){
442  validSet.erase(validSet.begin()+pruneHit[iHit-1]-count);
443  ++count;
444  }
445  }
446  }
447 }
448 
449 
450 std::vector <SeedCandidate> SETSeedFinder::
451 fillSeedCandidates(std::vector <MuonRecHitContainer> & allValidSets){
452  //---- we have the valid sets constructed; transform the information in an
453  //---- apropriate form; meanwhile - estimate the momentum for a given set
454 
455  // RPCs should not be used (no parametrization)
456  std::vector <SeedCandidate> seedCandidates_inCluster;
457  // calculate and fill the inputs needed
458  // loop over all valid sets
459  for(unsigned int iSet = 0;iSet<allValidSets.size();++iSet){
460  //
461  //std::cout<<" This is SET number : "<<iSet<<std::endl;
462  //for(unsigned int iHit = 0;iHit<allValidSets[iSet].size();++iHit){
463  //std::cout<<" measurements in the SET: iHit = "<<iHit<<" pos = "<<allValidSets[iSet][iHit]->globalPosition()<<
464  //" dim = "<<allValidSets[iSet][iHit]->dimension()<<std::endl;
465  //}
466 
467 
468  CLHEP::Hep3Vector momEstimate;
469  int chargeEstimate;
470  estimateMomentum(allValidSets[iSet], momEstimate, chargeEstimate);
471  MuonRecHitContainer MuonRecHitContainer_theSet_prep;
472  // currently hardcoded - will be in proper loop of course:
473 
474  SeedCandidate seedCandidates_inCluster_prep;
475  seedCandidates_inCluster_prep.theSet = allValidSets[iSet];
476  seedCandidates_inCluster_prep.momentum = momEstimate;
477  seedCandidates_inCluster_prep.charge = chargeEstimate;
478  seedCandidates_inCluster.push_back(seedCandidates_inCluster_prep);
479  // END estimateMomentum
480  }
481  return seedCandidates_inCluster;
482 }
483 
485  CLHEP::Hep3Vector & momEstimate, int & charge) const
486 {
487  int firstMeasurement = -1;
488  int lastMeasurement = -1;
489 
490  // don't use 2D measurements for momentum estimation
491 
492  //if( 4==allValidSets[iSet].front()->dimension() &&
493  //(allValidSets[iSet].front()->isCSC() || allValidSets[iSet].front()->isDT())){
494  //firstMeasurement = 0;
495  //}
496  //else{
497  // which is the "first" hit (4D)?
498  for(unsigned int iMeas = 0;iMeas<validSet.size();++iMeas){
499  if(4==validSet[iMeas]->dimension() &&
500  (validSet[iMeas]->isCSC() || validSet[iMeas]->isDT())){
501  firstMeasurement = iMeas;
502  break;
503  }
504  }
505  //}
506 
507  std::vector<double> momentum_estimate;
508  double pT = 0.;
511  // which is the second hit?
512  for(int loop = 0; loop<2; ++loop){// it is actually not used; to be removed
513  // this is the last measurement
514  if(!loop){// this is what is used currently
515  // 23.04.09 : it becomes a problem with introduction of ME42 chambers -
516  // the initial pT parametrization is incorrect for them
517  for(int iMeas = validSet.size()-1;iMeas>-1;--iMeas){
518  if(4==validSet[iMeas]->dimension() &&
519  (validSet[iMeas]->isCSC() || validSet[iMeas]->isDT()) &&
520  // below is a fix saying "don't use ME4 chambers for initial pT estimation";
521  // not using ME41 should not be a big loss too (and is more "symmetric" solution)
522  fabs(validSet[iMeas]->globalPosition().z())<1000.){
523  lastMeasurement = iMeas;
524  break;
525  }
526  }
527  }
528  else{
529  // this is the second measurement
530  for(unsigned int iMeas = 1;iMeas<validSet.size();++iMeas){
531  if(4==validSet[iMeas]->dimension() &&
532  (validSet[iMeas]->isCSC() || validSet[iMeas]->isDT())){
533  lastMeasurement = iMeas;
534  break;
535  }
536  }
537  }
538  // only 2D measurements (it should have been already abandoned)
539  if(-1==lastMeasurement && -1==firstMeasurement){
540  firstMeasurement = 0;
541  lastMeasurement = validSet.size()-1;
542  }
543  // because of the ME42 above lastMeasurement could be -1
544  else if(-1==lastMeasurement){
545  lastMeasurement = firstMeasurement;
546  }
547  else if(-1==firstMeasurement){
548  firstMeasurement = lastMeasurement;
549  }
550 
551  firstHit = validSet[firstMeasurement];
552  secondHit = validSet[lastMeasurement];
553  if(firstHit->isRPC() && secondHit->isRPC()){ // remove all RPCs from here?
554  momentum_estimate.push_back(300.);
555  momentum_estimate.push_back(300.);
556  }
557  else{
558  if(firstHit->isRPC()){
559  firstHit = secondHit;
560  }
561  else if(secondHit->isRPC()){
562  secondHit = firstHit;
563  }
564  //---- estimate pT given two hits
565  //std::cout<<" hits for initial pT estimate: first -> dim = "<<firstHit->dimension()<<" pos = "<<firstHit->globalPosition()<<
566  //" , second -> "<<" dim = "<<secondHit->dimension()<<" pos = "<<secondHit->globalPosition()<<std::endl;
567  //---- pT throws exception if hits are MB4
568  // (no coding for them - 2D hits in the outer station)
569  if(2==firstHit->dimension() && 2==secondHit->dimension()){
570  momentum_estimate.push_back(999999999.);
571  momentum_estimate.push_back(999999999.);
572  }
573  else{
574  momentum_estimate = thePtExtractor->pT_extract(firstHit, secondHit);
575  }
576  }
577  pT = fabs(momentum_estimate[0]);
578  if(1 || pT>40.){ //it is skipped; we have to look at least into number of hits in the chamber actually...
579  // and then decide which segment to use
580  // use the last measurement, otherwise use the second; this is to be investigated
581  break;
582  }
583  }
584 
585  const float pT_min = 1.99;// many hardcoded - remove them!
586  if(pT>3000.){
587  pT=3000.;
588  }
589  else if(pT<pT_min ){
590  if(pT>0){
591  pT=pT_min ;
592  }
593  else if(pT>(-1)*pT_min ){
594  pT=(-1)*pT_min ;
595  }
596  else if (pT<-3000.){
597  pT= -3000;
598  }
599  }
600  //std::cout<<" THE pT from the parametrization: "<<momentum_estimate[0]<<std::endl;
601  // estimate the charge of the track candidate from the delta phi of two segments:
602  //int charge = dPhi > 0 ? 1 : -1; // what we want is: dphi < 0 => charge = -1
603  charge = momentum_estimate[0]> 0 ? 1 : -1;
604 
605  // we have the pT - get the 3D momentum estimate as well
606 
607  // this is already final info:
608  double xHit = validSet[firstMeasurement]->globalPosition().x();
609  double yHit = validSet[firstMeasurement]->globalPosition().y();
610  double rHit = TMath::Sqrt(pow(xHit,2) + pow(yHit,2));
611 
612  double thetaInner = validSet[firstMeasurement]->globalPosition().theta();
613  // if some of the segments is missing r-phi measurement then we should
614  // use only the 4D phi estimate (also use 4D eta estimate only)
615  // the direction is not so important (it will be corrected)
616 
617  double rTrack = (pT /(0.3*3.8))*100.; //times 100 for conversion to cm!
618 
619  double par = -1.*(2./charge)*(TMath::ASin(rHit/(2*rTrack)));
620  double sinPar = TMath::Sin( par );
621  double cosPar = TMath::Cos( par );
622 
623  // calculate phi at coordinate origin (0,0,0).
624  double sinPhiH = 1./(2.*charge*rTrack)*(xHit + ((sinPar)/(cosPar-1.))*yHit);
625  double cosPhiH = -1./(2.*charge*rTrack)*(((sinPar)/(1.-cosPar))*xHit + yHit);
626 
627  // finally set the return vector
628 
629  // try out the reco info:
630  momEstimate = CLHEP::Hep3Vector(pT*cosPhiH, pT*sinPhiH, pT/TMath::Tan(thetaInner)); // should used into to theta directly here (rather than tan(atan2(...)))
631  //Hep3Vector momEstimate(6.97961, 5.89732, -50.0855);
632  const float minMomenum = 5.; //hardcoded - remove it! same in SETFilter
633  if (momEstimate.mag()<minMomenum){
634  int sign = (pT<0.) ? -1: 1;
635  pT = sign * (fabs(pT)+1);
636  CLHEP::Hep3Vector momEstimate2(pT*cosPhiH, pT*sinPhiH, pT/TMath::Tan(thetaInner));
637  momEstimate = momEstimate2;
638  if (momEstimate.mag()<minMomenum){
639  pT = sign * (fabs(pT)+1);
640  CLHEP::Hep3Vector momEstimate3(pT*cosPhiH, pT*sinPhiH, pT/TMath::Tan(thetaInner));
641  momEstimate = momEstimate3;
642  if (momEstimate.mag()<minMomenum){
643  pT = sign * (fabs(pT)+1);
644  CLHEP::Hep3Vector momEstimate4(pT*cosPhiH, pT*sinPhiH, pT/TMath::Tan(thetaInner));
645  momEstimate = momEstimate4;
646  }
647  }
648  }
649 }
650 
651 
652 
655 {
656  edm::OwnVector<TrackingRecHit> recHitsContainer;
657  for(unsigned int iHit = 0;iHit < hits.size();++iHit){
658  recHitsContainer.push_back(hits.at(iHit)->hit()->clone());
659  }
662  dir = alongMomentum;// why forward (for rechits) later?
663  }
664 
665  PTrajectoryStateOnDet const & seedTSOS =
666  trajectoryStateTransform::persistentState( firstTSOS, hits.at(0)->geographicalId().rawId());
667  TrajectorySeed seed(seedTSOS,recHitsContainer,dir);
668 
669  //MuonPatternRecoDumper debug;
670  //std::cout<<" firstTSOS = "<<debug.dumpTSOS(firstTSOS)<<std::endl;
671  //std::cout<<" iTraj = ???"<<" hits = "<<range.second-range.first<<std::endl;
672  //std::cout<<" nhits = "<<hits.size()<<std::endl;
673  //for(unsigned int iRH=0;iRH<hits.size();++iRH){
674  //std::cout<<" RH = "<<iRH+1<<" globPos = "<<hits.at(iRH)->globalPosition()<<std::endl;
675  //}
676  return seed;
677 }
678 
679 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool apply_prePruning
Definition: SETSeedFinder.h:61
const std::string metname
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
CLHEP::Hep3Vector momentum
Definition: SETFilter.h:36
bool operator()(TransientTrackingRecHit::ConstRecHitPointer hit_1, TransientTrackingRecHit::ConstRecHitPointer hit_2)
Definition: SETFilter.cc:31
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
std::pair< int, int > checkAngleDeviation(double dPhi_1, double dPhi_2) const
float localZ(const GlobalPoint &gp) const
Fast access to distance from plane for a point.
Definition: Plane.h:52
void pre_prune(MuonRecHitContainer &validSet) const
bool useSegmentsInTrajectory
Definition: SETSeedFinder.h:62
PropagationDirection
TrajectorySeed makeSeed(const TrajectoryStateOnSurface &tsos, const TransientTrackingRecHit::ConstRecHitContainer &hits) const
double charge(const std::vector< uint8_t > &Ampls)
std::vector< MuonRecHitContainer > sortByLayer(MuonRecHitContainer &cluster) const
void estimateMomentum(const MuonRecHitContainer &validSet, CLHEP::Hep3Vector &momentum, int &charge) const
void push_back(D *&d)
Definition: OwnVector.h:273
static const int CSC
Definition: MuonSubdetId.h:15
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
tuple result
Definition: query.py:137
MuonTransientTrackingRecHit::MuonRecHitContainer MuonRecHitContainer
Definition: SETSeedFinder.h:15
struct sorter sortSegRadius
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
virtual void seeds(const MuonRecHitContainer &cluster, std::vector< TrajectorySeed > &result)
MuonSeedPtExtractor * thePtExtractor
std::vector< ConstRecHitPointer > ConstRecHitContainer
void limitCombinatorics(std::vector< MuonRecHitContainer > &MuonRecHitContainer_perLayer)
int ring() const
Definition: CSCDetId.h:77
Definition: DetId.h:20
std::vector< MuonRecHitContainer > findAllValidSets(const std::vector< MuonRecHitContainer > &MuonRecHitContainer_perLayer)
MuonServiceProxy * theService
Definition: SETSeedFinder.h:59
const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
void validSetsPrePruning(std::vector< MuonRecHitContainer > &allValidSets)
MuonTransientTrackingRecHit::MuonRecHitContainer theSet
Definition: SETFilter.h:35
int station() const
Definition: CSCDetId.h:88
static const int DT
Definition: MuonSubdetId.h:14
dbl *** dir
Definition: mlp_gen.cc:35
std::vector< SeedCandidate > fillSeedCandidates(std::vector< MuonRecHitContainer > &allValidSets)
virtual std::vector< double > pT_extract(MuonTransientTrackingRecHit::ConstMuonRecHitPointer firstHit, MuonTransientTrackingRecHit::ConstMuonRecHitPointer secondHit) const
int station() const
Return the station number.
Definition: DTChamberId.h:53
uint32_t dimension(pat::CandKinResolution::Parametrization parametrization)
Returns the number of free parameters in a parametrization (3 or 4)
SETSeedFinder(const edm::ParameterSet &pset)
tuple size
Write out results.
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
Definition: DDAxes.h:10