CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrajectorySeedProducer.cc
Go to the documentation of this file.
1 #include <memory>
2 
7 
13 
16 
18 
21 
22 
26 
29 
32 
33 //Propagator withMaterial
35 //analyticalpropagator
36 //#include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
37 
38 //
39 
40 //for debug only
41 //#define FAMOS_DEBUG
42 
44 {
45 
46  // The input tag for the beam spot
47  theBeamSpot = conf.getParameter<edm::InputTag>("beamSpot");
48 
49  // The name of the TrajectorySeed Collections
50  seedingAlgo = conf.getParameter<std::vector<std::string> >("seedingAlgo");
51  for ( unsigned i=0; i<seedingAlgo.size(); ++i )
52  produces<TrajectorySeedCollection>(seedingAlgo[i]);
53 
54  // The smallest true pT for a track candidate
55  pTMin = conf.getParameter<std::vector<double> >("pTMin");
56  if ( pTMin.size() != seedingAlgo.size() )
57  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
58  << " WARNING : pTMin does not have the proper size "
59  << std::endl;
60 
61  for ( unsigned i=0; i<pTMin.size(); ++i )
62  pTMin[i] *= pTMin[i]; // Cut is done of perp2() - CPU saver
63 
64  // The smallest number of Rec Hits for a track candidate
65  minRecHits = conf.getParameter<std::vector<unsigned int> >("minRecHits");
66  if ( minRecHits.size() != seedingAlgo.size() )
67  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
68  << " WARNING : minRecHits does not have the proper size "
69  << std::endl;
70  // Set the overall number hits to be checked
71  absMinRecHits = 0;
72  for ( unsigned ialgo=0; ialgo<minRecHits.size(); ++ialgo )
73  if ( minRecHits[ialgo] > absMinRecHits ) absMinRecHits = minRecHits[ialgo];
74 
75  // The smallest true impact parameters (d0 and z0) for a track candidate
76  maxD0 = conf.getParameter<std::vector<double> >("maxD0");
77  if ( maxD0.size() != seedingAlgo.size() )
78  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
79  << " WARNING : maxD0 does not have the proper size "
80  << std::endl;
81 
82  maxZ0 = conf.getParameter<std::vector<double> >("maxZ0");
83  if ( maxZ0.size() != seedingAlgo.size() )
84  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
85  << " WARNING : maxZ0 does not have the proper size "
86  << std::endl;
87 
88  // The name of the hit producer
89  hitProducer = conf.getParameter<edm::InputTag>("HitProducer");
90 
91  // The cuts for seed cleaning
92  seedCleaning = conf.getParameter<bool>("seedCleaning");
93 
94  // Number of hits needed for a seed
95  numberOfHits = conf.getParameter<std::vector<unsigned int> >("numberOfHits");
96  if ( numberOfHits.size() != seedingAlgo.size() )
97  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
98  << " WARNING : numberOfHits does not have the proper size "
99  << std::endl;
100 
101  // Seeding based on muons
102  selectMuons = conf.getParameter<bool>("selectMuons");
103 
104  // Layers
105  newSyntax = conf.getParameter<bool>("newSyntax");
106  if (newSyntax) {
107  std::vector<std::string> layerList = conf.getParameter<std::vector<std::string> >("layerList");
108  //for (unsigned i=0; i<layerList.size();i++) std::cout << "------- Layers = " << layerList[i] << std::endl;
109 
110  for(std::vector<std::string>::const_iterator it=layerList.begin(); it < layerList.end(); ++it) {
111  std::vector<LayerSpec> tempResult;
112  std::string line = *it;
114  while (pos != std::string::npos) {
115  pos=line.find("+");
116  std::string layer = line.substr(0, pos);
117  //
118  LayerSpec layerSpec;
119  layerSpec.name = line.substr(0, pos);;
120  // Possible names: BPix(1-3) || FPix(1-2)(pos, neg) || TIB(1-4) || TID(1-3)(pos, neg) || TOB(1-6) || TEC(1-9)(pos, neg)
121  // BPix(1-3)
122  if (layerSpec.name.substr(0,4)=="BPix" ) {
123  layerSpec.subDet=PXB;
124  layerSpec.side=BARREL;
125  layerSpec.idLayer = std::atoi(layerSpec.name.substr(4,1).c_str());
126  if (layerSpec.idLayer>3 || layerSpec.idLayer==0) {
127  throw cms::Exception("FastSimulation/Tracking/python")
128  << "Bad data naming in IterativeInitialStep_cff.py iterativeInitialSeeds.layerList" << std::endl
129  << "Layers: " << layerSpec.name << ", number needs to be in range 1-3" << std::endl;
130  }
131  }
132  // FPix(1-2)(pos, neg)
133  else if (layerSpec.name.substr(0,4)=="FPix" ) {
134  layerSpec.subDet=PXD;
135  if(layerSpec.name.substr(layerSpec.name.size()-3)=="pos")
136  layerSpec.side = POS_ENDCAP;
137  else //no validation if it's not neg
138  layerSpec.side = NEG_ENDCAP;
139  layerSpec.idLayer = std::atoi(layerSpec.name.substr(4,1).c_str());
140  if (layerSpec.idLayer>2 || layerSpec.idLayer==0) {
141  throw cms::Exception("FastSimulation/Tracking/python")
142  << "Bad data naming in IterativeInitialStep_cff.py iterativeInitialSeeds.layerList" << std::endl
143  << "Layers: " << layerSpec.name << ", number needs to be in range 1-2" << std::endl;
144  }
145  }
146  // TIB(1-4)
147  else if (layerSpec.name.substr(0,3)=="TIB" || layerSpec.name.substr(0,4)=="MTIB") {
148  layerSpec.subDet=TIB;
149  layerSpec.side=BARREL;
150  layerSpec.idLayer = std::atoi(layerSpec.name.substr(layerSpec.name.find('B')+1,1).c_str());
151  if (layerSpec.idLayer>4 || layerSpec.idLayer==0) {
152  throw cms::Exception("FastSimulation/Tracking/python")
153  << "Bad data naming in IterativeInitialStep_cff.py iterativeInitialSeeds.layerList" << std::endl
154  << "Layers: " << layerSpec.name << ", number needs to be in range 1-4" << std::endl;
155  }
156  }
157  // TID(1-3)(pos, neg)
158  else if (layerSpec.name.substr(0,3)=="TID" || layerSpec.name.substr(0,4)=="MTID") {
159  layerSpec.subDet=TID;
160  if(layerSpec.name.substr(layerSpec.name.size()-3)=="pos")
161  layerSpec.side = POS_ENDCAP;
162  else
163  layerSpec.side = NEG_ENDCAP;
164  layerSpec.idLayer = std::atoi(layerSpec.name.substr(layerSpec.name.find('D')+1,1).c_str());
165  if (layerSpec.idLayer>3 || layerSpec.idLayer==0) {
166  throw cms::Exception("FastSimulation/Tracking/python")
167  << "Bad data naming in IterativeInitialStep_cff.py iterativeInitialSeeds.layerList" << std::endl
168  << "Layers: " << layerSpec.name << ", number needs to be in range 1-3" << std::endl;
169  }
170  }
171  // TOB(1-6)
172  else if (layerSpec.name.substr(0,3)=="TOB" ) {
173  layerSpec.subDet=TOB;
174  layerSpec.side=BARREL;
175  layerSpec.idLayer = std::atoi(layerSpec.name.substr(3,1).c_str());
176  if (layerSpec.idLayer>6 || layerSpec.idLayer==0) {
177  throw cms::Exception("FastSimulation/Tracking/python")
178  << "Bad data naming in IterativeInitialStep_cff.py iterativeInitialSeeds.layerList" << std::endl
179  << "Layers: " << layerSpec.name << ", number needs to be in range 1-6" << std::endl;
180  }
181  }
182  // TEC(1-9)(pos, neg)
183  else if (layerSpec.name.substr(0,3)=="TEC" || layerSpec.name.substr(0,4)=="MTEC") {
184  layerSpec.subDet=TEC;
185  if(layerSpec.name.substr(layerSpec.name.size()-3)=="pos")
186  layerSpec.side = POS_ENDCAP;
187  else
188  layerSpec.side = NEG_ENDCAP;
189  layerSpec.idLayer = std::atoi(layerSpec.name.substr(layerSpec.name.find('C')+1,1).c_str());
190  if (layerSpec.idLayer>9 || layerSpec.idLayer==0) {
191  throw cms::Exception("FastSimulation/Tracking/python")
192  << "Bad data naming in IterativeInitialStep_cff.py iterativeInitialSeeds.layerList" << std::endl
193  << "Layers: " << layerSpec.name << ", number needs to be in range 1-9" << std::endl;
194  }
195  }
196  else {
197  throw cms::Exception("FastSimulation/Tracking/python")
198  << "Bad data naming in IterativeInitialStep_cff.py iterativeInitialSeeds.layerList" << std::endl
199  << "Layer: " << layerSpec.name << ", shouldn't exist" << std::endl
200  << "Case sensitive names: BPix FPix TIB MTIB TID MTID TOB TEC MTEC" << std::endl;
201  }
202  //
203  tempResult.push_back(layerSpec);
204  line=line.substr(pos+1,std::string::npos);
205  }
206  theLayersInSets.push_back(tempResult);
207  }
208 
209  //prints theLayersInSets
210  /* for (std::vector<std::vector<LayerSpec> >::const_iterator it = theLayersInSets.begin(); it != theLayersInSets.end(); ++it) {
211  for (std::vector<LayerSpec>::const_iterator is = it->begin(); is != it->end(); ++is) {
212  std::cout << is->name << " | " << is->subDet << " | " << is->idLayer << " | " << is->side << std::endl;
213  }
214  std::cout << "---" << std::endl;
215  } */
216  } else {
217  // TO BE DELETED (AG)
219  conf.getParameter<std::vector<unsigned int> >("firstHitSubDetectorNumber");
220  if ( firstHitSubDetectorNumber.size() != seedingAlgo.size() )
221  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
222  << " WARNING : firstHitSubDetectorNumber does not have the proper size "
223  << std::endl;
224 
225  std::vector<unsigned int> firstSubDets =
226  conf.getParameter<std::vector<unsigned int> >("firstHitSubDetectors");
227  unsigned isub1 = 0;
228  unsigned check1 = 0;
229  firstHitSubDetectors.resize(seedingAlgo.size());
230  for ( unsigned ialgo=0; ialgo<firstHitSubDetectorNumber.size(); ++ialgo ) {
231  check1 += firstHitSubDetectorNumber[ialgo];
232  for ( unsigned idet=0; idet<firstHitSubDetectorNumber[ialgo]; ++idet ) {
233  firstHitSubDetectors[ialgo].push_back(firstSubDets[isub1++]);
234  }
235  }
236  if ( firstSubDets.size() != check1 )
237  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
238  << " WARNING : firstHitSubDetectors does not have the proper size (should be " << check1 << ")"
239  << std::endl;
240 
241 
243  conf.getParameter<std::vector<unsigned int> >("secondHitSubDetectorNumber");
244  if ( secondHitSubDetectorNumber.size() != seedingAlgo.size() )
245  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
246  << " WARNING : secondHitSubDetectorNumber does not have the proper size "
247  << std::endl;
248 
249  std::vector<unsigned int> secondSubDets =
250  conf.getParameter<std::vector<unsigned int> >("secondHitSubDetectors");
251  unsigned isub2 = 0;
252  unsigned check2 = 0;
253  secondHitSubDetectors.resize(seedingAlgo.size());
254  for ( unsigned ialgo=0; ialgo<secondHitSubDetectorNumber.size(); ++ialgo ) {
255  check2 += secondHitSubDetectorNumber[ialgo];
256  for ( unsigned idet=0; idet<secondHitSubDetectorNumber[ialgo]; ++idet ) {
257  secondHitSubDetectors[ialgo].push_back(secondSubDets[isub2++]);
258  }
259  }
260  if ( secondSubDets.size() != check2 )
261  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
262  << " WARNING : secondHitSubDetectors does not have the proper size (should be " << check2 << ")"
263  << std::endl;
264 
266  conf.getParameter<std::vector<unsigned int> >("thirdHitSubDetectorNumber");
267  if ( thirdHitSubDetectorNumber.size() != seedingAlgo.size() )
268  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
269  << " WARNING : thirdHitSubDetectorNumber does not have the proper size "
270  << std::endl;
271 
272  std::vector<unsigned int> thirdSubDets =
273  conf.getParameter<std::vector<unsigned int> >("thirdHitSubDetectors");
274  unsigned isub3 = 0;
275  unsigned check3 = 0;
276  thirdHitSubDetectors.resize(seedingAlgo.size());
277  for ( unsigned ialgo=0; ialgo<thirdHitSubDetectorNumber.size(); ++ialgo ) {
278  check3 += thirdHitSubDetectorNumber[ialgo];
279  for ( unsigned idet=0; idet<thirdHitSubDetectorNumber[ialgo]; ++idet ) {
280  thirdHitSubDetectors[ialgo].push_back(thirdSubDets[isub3++]);
281  }
282  }
283  if ( thirdSubDets.size() != check3 )
284  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
285  << " WARNING : thirdHitSubDetectors does not have the proper size (should be " << check3 << ")"
286  << std::endl;
287  }
288 
289  originRadius = conf.getParameter<std::vector<double> >("originRadius");
290  if ( originRadius.size() != seedingAlgo.size() )
291  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
292  << " WARNING : originRadius does not have the proper size "
293  << std::endl;
294 
295  originHalfLength = conf.getParameter<std::vector<double> >("originHalfLength");
296  if ( originHalfLength.size() != seedingAlgo.size() )
297  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
298  << " WARNING : originHalfLength does not have the proper size "
299  << std::endl;
300 
301  originpTMin = conf.getParameter<std::vector<double> >("originpTMin");
302  if ( originpTMin.size() != seedingAlgo.size() )
303  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
304  << " WARNING : originpTMin does not have the proper size "
305  << std::endl;
306 
307  primaryVertices = conf.getParameter<std::vector<edm::InputTag> >("primaryVertices");
308  if ( primaryVertices.size() != seedingAlgo.size() )
309  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
310  << " WARNING : primaryVertices does not have the proper size "
311  << std::endl;
312 
313  zVertexConstraint = conf.getParameter<std::vector<double> >("zVertexConstraint");
314  if ( zVertexConstraint.size() != seedingAlgo.size() )
315  throw cms::Exception("FastSimulation/TrajectorySeedProducer ")
316  << " WARNING : zVertexConstraint does not have the proper size "
317  << std::endl;
318  //removed - }
319 
320  // consumes
321  beamSpotToken = consumes<reco::BeamSpot>(theBeamSpot);
322  edm::InputTag _label("famosSimHits");
323  simTrackToken = consumes<edm::SimTrackContainer>(_label);
324  simVertexToken = consumes<edm::SimVertexContainer>(_label);
325  recHitToken = consumes<SiTrackerGSMatchedRecHit2DCollection>(hitProducer);
326  for ( unsigned ialgo=0; ialgo<seedingAlgo.size(); ++ialgo ) {
327  _label = edm::InputTag(primaryVertices[ialgo]);
328  recoVertexToken.push_back(consumes<reco::VertexCollection>(_label));
329  }
330 }
331 
332 // Virtual destructor needed.
334 
335  if(thePropagator) delete thePropagator;
336 
337  // do nothing
338 #ifdef FAMOS_DEBUG
339  std::cout << "TrajectorySeedProducer destructed" << std::endl;
340 #endif
341 
342 }
343 
344 void
346 
347  //services
348  // es.get<TrackerRecoGeometryRecord>().get(theGeomSearchTracker);
349 
353 
354 
355  es.get<IdealMagneticFieldRecord>().get(magField);
356  es.get<TrackerDigiGeometryRecord>().get(geometry);
357  es.get<MagneticFieldMapRecord>().get(magFieldMap);
358 
359  theMagField = &(*magField);
360  theGeometry = &(*geometry);
361  theFieldMap = &(*magFieldMap);
362 
364 
365  const GlobalPoint g(0.,0.,0.);
366 
367 }
368 
369  // Functions that gets called by framework every event
370 void
372 
373 
374  // if( seedingAlgo[0] == "FourthPixelLessPairs") std::cout << "Seed producer in 4th iteration " << std::endl;
375 
376 #ifdef FAMOS_DEBUG
377  std::cout << "################################################################" << std::endl;
378  std::cout << " TrajectorySeedProducer produce init " << std::endl;
379 #endif
380 
381  //Retrieve tracker topology from geometry
383  es.get<IdealGeometryRecord>().get(tTopoHand);
384  const TrackerTopology *tTopo=tTopoHand.product();
385 
386 
387  unsigned nSimTracks = 0;
388  unsigned nTracksWithHits = 0;
389  unsigned nTracksWithPT = 0;
390  unsigned nTracksWithD0Z0 = 0;
391  // unsigned nTrackCandidates = 0;
392  PTrajectoryStateOnDet initialState;
393 
394  // Output
395  std::vector<TrajectorySeedCollection*>
396  output(seedingAlgo.size(),static_cast<TrajectorySeedCollection*>(0));
397  for ( unsigned ialgo=0; ialgo<seedingAlgo.size(); ++ialgo ) {
398  // std::auto_ptr<TrajectorySeedCollection> p(new TrajectorySeedCollection );
399  output[ialgo] = new TrajectorySeedCollection;
400  }
401 
402  // Beam spot
403  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
404  e.getByToken(beamSpotToken,recoBeamSpotHandle);
405  math::XYZPoint BSPosition_ = recoBeamSpotHandle->position();
406 
407  //not used anymore. take the value from the py
408 
409  //double sigmaZ=recoBeamSpotHandle->sigmaZ();
410  //double sigmaZ0Error=recoBeamSpotHandle->sigmaZ0Error();
411  //double sigmaz0=sqrt(sigmaZ*sigmaZ+sigmaZ0Error*sigmaZ0Error);
412  x0 = BSPosition_.X();
413  y0 = BSPosition_.Y();
414  z0 = BSPosition_.Z();
415 
416  // SimTracks and SimVertices
418  e.getByToken(simTrackToken,theSimTracks);
419 
421  e.getByToken(simVertexToken,theSimVtx);
422 
423 #ifdef FAMOS_DEBUG
424  std::cout << " Step A: SimTracks found " << theSimTracks->size() << std::endl;
425 #endif
426 
427  // edm::Handle<SiTrackerGSRecHit2DCollection> theGSRecHits;
429  e.getByToken(recHitToken, theGSRecHits);
430 
431  // No tracking attempted if no hits (but put an empty collection in the event)!
432 #ifdef FAMOS_DEBUG
433  std::cout << " Step B: Full GS RecHits found " << theGSRecHits->size() << std::endl;
434 #endif
435  if(theGSRecHits->size() == 0) {
436  for ( unsigned ialgo=0; ialgo<seedingAlgo.size(); ++ialgo ) {
437  std::auto_ptr<TrajectorySeedCollection> p(output[ialgo]);
438  e.put(p,seedingAlgo[ialgo]);
439  }
440  return;
441  }
442 
443  // Primary vertices
444  vertices = std::vector<const reco::VertexCollection*>
445  (seedingAlgo.size(),static_cast<const reco::VertexCollection*>(0));
446  for ( unsigned ialgo=0; ialgo<seedingAlgo.size(); ++ialgo ) {
447  //PAT Attempt!!!!
448  //originHalfLength[ialgo] = 3.*sigmaz0; // Overrides the configuration
450  bool isVertexCollection = e.getByToken(recoVertexToken[ialgo],aHandle);
451  if (!isVertexCollection ) continue;
452  vertices[ialgo] = &(*aHandle);
453  }
454 
455 #ifdef FAMOS_DEBUG
456  std::cout << " Step C: Loop over the RecHits, track by track " << std::endl;
457 #endif
458 
459  // The vector of simTrack Id's carrying GSRecHits
460  const std::vector<unsigned> theSimTrackIds = theGSRecHits->ids();
461 
462  // loop over SimTrack Id's
463  for ( unsigned tkId=0; tkId != theSimTrackIds.size(); ++tkId ) {
464 
465 #ifdef FAMOS_DEBUG
466  std::cout << "Track number " << tkId << std::endl;
467 #endif
468 
469  ++nSimTracks;
470  unsigned simTrackId = theSimTrackIds[tkId];
471  const SimTrack& theSimTrack = (*theSimTracks)[simTrackId];
472 #ifdef FAMOS_DEBUG
473  std::cout << "Pt = " << std::sqrt(theSimTrack.momentum().Perp2())
474  << " eta " << theSimTrack.momentum().Eta()
475  << " pdg ID " << theSimTrack.type()
476  << std::endl;
477 #endif
478 
479  // Select only muons, if requested
480  if (selectMuons && abs(theSimTrack.type()) != 13) continue;
481 
482  // Check that the sim track comes from the main vertex (loose cut)
483  int vertexIndex = theSimTrack.vertIndex();
484  const SimVertex& theSimVertex = (*theSimVtx)[vertexIndex];
485 #ifdef FAMOS_DEBUG
486  std::cout << " o SimTrack " << theSimTrack << std::endl;
487  std::cout << " o SimVertex " << theSimVertex << std::endl;
488 #endif
489 
490  BaseParticlePropagator theParticle =
492  RawParticle(XYZTLorentzVector(theSimTrack.momentum().px(),
493  theSimTrack.momentum().py(),
494  theSimTrack.momentum().pz(),
495  theSimTrack.momentum().e()),
496  XYZTLorentzVector(theSimVertex.position().x(),
497  theSimVertex.position().y(),
498  theSimVertex.position().z(),
499  theSimVertex.position().t())),
500  0.,0.,4.);
501  theParticle.setCharge((*theSimTracks)[simTrackId].charge());
502 
503  SiTrackerGSMatchedRecHit2DCollection::range theRecHitRange = theGSRecHits->get(simTrackId);
504  SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorBegin = theRecHitRange.first;
505  SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorEnd = theRecHitRange.second;
510 
511  // Check the number of layers crossed
512  unsigned numberOfRecHits = 0;
513  TrackerRecHit previousHit, currentHit;
514  for ( iterRecHit = theRecHitRangeIteratorBegin;
515  iterRecHit != theRecHitRangeIteratorEnd;
516  ++iterRecHit) {
517  previousHit = currentHit;
518  currentHit = TrackerRecHit(&(*iterRecHit),theGeometry,tTopo);
519  if ( currentHit.isOnTheSameLayer(previousHit) ) continue;
520  ++numberOfRecHits;
521  if ( numberOfRecHits == absMinRecHits ) break;
522  }
523 
524  // Loop on the successive seedings
525  for ( unsigned int ialgo = 0; ialgo < seedingAlgo.size(); ++ialgo ) {
526 #ifdef FAMOS_DEBUG
527  std::cout << "Algo " << seedingAlgo[ialgo] << std::endl;
528 #endif
529 
530  // Request a minimum number of RecHits for the track to give a seed.
531 #ifdef FAMOS_DEBUG
532  std::cout << "The number of RecHits = " << numberOfRecHits << std::endl;
533 #endif
534  if ( numberOfRecHits < minRecHits[ialgo] ) continue;
535  ++nTracksWithHits;
536 
537  // Request a minimum pT for the sim track
538  if ( theSimTrack.momentum().Perp2() < pTMin[ialgo] ) continue;
539  ++nTracksWithPT;
540 
541  // Cut on sim track impact parameters
542  if ( theParticle.xyImpactParameter(x0,y0) > maxD0[ialgo] ) continue;
543  if ( fabs( theParticle.zImpactParameter(x0,y0) - z0 ) > maxZ0[ialgo] ) continue;
544  ++nTracksWithD0Z0;
545 
546  std::vector<TrackerRecHit >
547  theSeedHits(numberOfHits[ialgo],
548  static_cast<TrackerRecHit >(TrackerRecHit()));
549  TrackerRecHit& theSeedHits0 = theSeedHits[0];
550  TrackerRecHit& theSeedHits1 = theSeedHits[1];
551  TrackerRecHit& theSeedHits2 = theSeedHits[2];
552 
553  bool compatible = false;
554 
555  for ( iterRecHit1 = theRecHitRangeIteratorBegin; iterRecHit1 != theRecHitRangeIteratorEnd; ++iterRecHit1) {
556  //std::cout << (*iterRecHit1).localPosition().phi() << " | J - iterRecHit1"<< std::endl; //
557  theSeedHits[0] = TrackerRecHit(&(*iterRecHit1),theGeometry,tTopo);
558 #ifdef FAMOS_DEBUG
559  std::cout << "The first hit position = " << theSeedHits0.globalPosition() << std::endl;
560  std::cout << "The first hit subDetId = " << theSeedHits0.subDetId() << std::endl;
561  std::cout << "The first hit layer = " << theSeedHits0.layerNumber() << std::endl;
562 #endif
563  // Check if inside the requested detectors
564  bool isInside = true;
565  if (!selectMuons) {
566  //(newSyntax) ? std::cout << "TRUE " : std::cout << "FALSE "; //J-
567  if (newSyntax)
568  isInside = false; // AG placeholder true
569  else
570  isInside = theSeedHits0.subDetId() < firstHitSubDetectors[ialgo][0];
571  // bool isInside = theSeedHits0.subDetId() < firstHitSubDetectors[ialgo][0];
572  if ( isInside ) continue;
573  }
574 
575  // Check if on requested detectors
576  // bool isOndet = theSeedHits0.isOnRequestedDet(firstHitSubDetectors[ialgo]);
577  bool isOndet = true;
578  if (!selectMuons) {
579  if (newSyntax)
580  isOndet = theSeedHits[0].isOnRequestedDet(theLayersInSets);
581  else
582  isOndet = theSeedHits0.isOnRequestedDet(firstHitSubDetectors[ialgo], seedingAlgo[ialgo]);
583  //std::cout << firstHitSubDetectors[ialgo][0] << " | " << seedingAlgo[ialgo] << " " << std::endl; //seedingAlgo[iAlgo]: PixelTriplet, LowPtPixelTriplets, PixelPair, DetachedPixelTriplets, MixedTriplets, PixelLessPairs, TobTecLayerPairs......
584  // bool isOndet = theSeedHits0.isOnRequestedDet(firstHitSubDetectors[ialgo], seedingAlgo[ialgo]);
585  // if ( !isOndet ) break;
586  if ( !isOndet ) continue;
587  }
588 
589 #ifdef FAMOS_DEBUG
590  std::cout << "Apparently the first hit is on the requested detector! " << std::endl;
591 #endif
592  for ( iterRecHit2 = iterRecHit1+1; iterRecHit2 != theRecHitRangeIteratorEnd; ++iterRecHit2) {
593  theSeedHits[1] = TrackerRecHit(&(*iterRecHit2),theGeometry,tTopo);
594 #ifdef FAMOS_DEBUG
595  std::cout << "The second hit position = " << theSeedHits1.globalPosition() << std::endl;
596  std::cout << "The second hit subDetId = " << theSeedHits1.subDetId() << std::endl;
597  std::cout << "The second hit layer = " << theSeedHits1.layerNumber() << std::endl;
598 #endif
599 
600  if (!selectMuons) {
601  // Check if inside the requested detectors
602  if (newSyntax)
603  isInside = false; // AG placeholder true
604  else
605  isInside = theSeedHits1.subDetId() < secondHitSubDetectors[ialgo][0];
606  if ( isInside ) continue;
607 
608  // Check if on requested detectors
609  if (newSyntax)
610  isOndet = theSeedHits[0].isOnRequestedDet(theLayersInSets, theSeedHits[1]);
611  else
612  isOndet = theSeedHits1.isOnRequestedDet(secondHitSubDetectors[ialgo], seedingAlgo[ialgo]);
613  if ( !isOndet ) break;
614  }
615  // Check if on the same layer as previous hit
616  if ( theSeedHits1.isOnTheSameLayer(theSeedHits0) ) continue;
617 
618 #ifdef FAMOS_DEBUG
619  std::cout << "Apparently the second hit is on the requested detector! " << std::endl;
620 #endif
621  GlobalPoint gpos1 = theSeedHits0.globalPosition();
622  GlobalPoint gpos2 = theSeedHits1.globalPosition();
623  bool forward = theSeedHits0.isForward();
624  double error = std::sqrt(theSeedHits0.largerError()+theSeedHits1.largerError());
625  // compatible = compatibleWithVertex(gpos1,gpos2,ialgo);
626  //added out of desperation
627  if(seedingAlgo[0] == "PixelLess" || seedingAlgo[0] == "TobTecLayerPairs"){
628  compatible = true;
629  //std::cout << "Algo " << seedingAlgo[0] << "Det/layer = " << theSeedHits0.subDetId() << "/" << theSeedHits0.layerNumber() << std::endl;
630  } else {
631  compatible = compatibleWithBeamAxis(gpos1,gpos2,error,forward,ialgo);
632  }
633 
634 #ifdef FAMOS_DEBUG
635  std::cout << "Algo" << seedingAlgo[0] << "\t Are the two hits compatible with the PV? " << compatible << std::endl;
636 #endif
637  // Reject non suited pairs
638  if ( !compatible ) continue;
639 
640 #ifdef FAMOS_DEBUG
641  std::cout << "Pair kept! " << std::endl;
642 #endif
643 
644  // Leave here if only two hits are required.
645  if ( numberOfHits[ialgo] == 2 ) break;
646 
647  compatible = false;
648  // Check if there is a third satisfying hit otherwise
649  for ( iterRecHit3 = iterRecHit2+1; iterRecHit3 != theRecHitRangeIteratorEnd; ++iterRecHit3) {
650 
651  theSeedHits[2] = TrackerRecHit(&(*iterRecHit3),theGeometry,tTopo);
652 #ifdef FAMOS_DEBUG
653  std::cout << "The third hit position = " << theSeedHits2.globalPosition() << std::endl;
654  std::cout << "The third hit subDetId = " << theSeedHits2.subDetId() << std::endl;
655  std::cout << "The third hit layer = " << theSeedHits2.layerNumber() << std::endl;
656 #endif
657 
658  if (!selectMuons) {
659  // Check if inside the requested detectors
660  if (newSyntax)
661  isInside = false; // AG placeholder
662  else
663  isInside = theSeedHits2.subDetId() < thirdHitSubDetectors[ialgo][0];
664  if ( isInside ) continue;
665 
666  // Check if on requested detectors
667  if (newSyntax)
668  isOndet = theSeedHits[0].isOnRequestedDet(theLayersInSets, theSeedHits[1], theSeedHits[2]);
669  else
670  isOndet = theSeedHits2.isOnRequestedDet(thirdHitSubDetectors[ialgo], seedingAlgo[ialgo]);
671  // if ( !isOndet ) break;
672  if ( !isOndet ) continue;
673  }
674 
675  // Check if on the same layer as previous hit
676  compatible = !(theSeedHits2.isOnTheSameLayer(theSeedHits1));
677 
678 #ifdef FAMOS_DEBUG
679  if ( compatible )
680  std::cout << "Apparently the third hit is on the requested detector! " << std::endl;
681 #endif
682 
683  if ( compatible ) break;
684 
685  }
686 
687  if ( compatible ) break;
688 
689  }
690 
691  if ( compatible ) break;
692 
693  }
694 
695  // There is no compatible seed for this track with this seeding algorithm
696  // Go to next algo
697  if ( !compatible ) continue;
698 #ifdef FAMOS_DEBUG
699  if ( compatible )
700  std::cout << "@@@ There is at least a compatible seed" << std::endl;
701  else
702  std::cout << "@@@ There is no compatible seed" << std::endl;
703 #endif
704 
705 
706 #ifdef FAMOS_DEBUG
707  std::cout << "Preparing to create the TrajectorySeed" << std::endl;
708 #endif
709  // The seed is validated -> include in the collection
710  // 1) Create the vector of RecHits
712  for ( unsigned ih=0; ih<theSeedHits.size(); ++ih ) {
713  TrackingRecHit* aTrackingRecHit = theSeedHits[ih].hit()->clone();
714  recHits.push_back(aTrackingRecHit);
715  }
716 #ifdef FAMOS_DEBUG
717  std::cout << "with " << recHits.size() << " hits." << std::endl;
718 #endif
719 
720  // 2) Create the initial state
721  // a) origin vertex
722  GlobalPoint position((*theSimVtx)[vertexIndex].position().x(),
723  (*theSimVtx)[vertexIndex].position().y(),
724  (*theSimVtx)[vertexIndex].position().z());
725 
726  // b) initial momentum
727  GlobalVector momentum( (*theSimTracks)[simTrackId].momentum().x() ,
728  (*theSimTracks)[simTrackId].momentum().y() ,
729  (*theSimTracks)[simTrackId].momentum().z() );
730  // c) electric charge
731  float charge = (*theSimTracks)[simTrackId].charge();
732  // -> inital parameters
733  GlobalTrajectoryParameters initialParams(position,momentum,(int)charge,theMagField);
734  // -> large initial errors
735  AlgebraicSymMatrix55 errorMatrix= AlgebraicMatrixID();
736  // errorMatrix = errorMatrix * 10;
737 
738  //this line help the fit succeed in the case of pixelless tracks (4th and 5th iteration)
739  //for the future: probably the best thing is to use the mini-kalmanFilter
740  if(theSeedHits0.subDetId() !=1 || theSeedHits0.subDetId() !=2) errorMatrix = errorMatrix * 0.0000001;
741 
742 
743 
744 #ifdef FAMOS_DEBUG
745  std::cout << "TrajectorySeedProducer: SimTrack parameters " << std::endl;
746  std::cout << "\t\t pT = " << (*theSimTracks)[simTrackId].momentum().Pt() << std::endl;
747  std::cout << "\t\t eta = " << (*theSimTracks)[simTrackId].momentum().Eta() << std::endl;
748  std::cout << "\t\t phi = " << (*theSimTracks)[simTrackId].momentum().Phi() << std::endl;
749  std::cout << "TrajectorySeedProducer: AlgebraicSymMatrix " << errorMatrix << std::endl;
750 #endif
751  CurvilinearTrajectoryError initialError(errorMatrix);
752  // -> initial state
753  FreeTrajectoryState initialFTS(initialParams, initialError);
754 #ifdef FAMOS_DEBUG
755  std::cout << "TrajectorySeedProducer: FTS momentum " << initialFTS.momentum() << std::endl;
756 #endif
757  // const GeomDetUnit* initialLayer = theGeometry->idToDetUnit( recHits.front().geographicalId() );
758  const GeomDet* initialLayer = theGeometry->idToDet( recHits.front().geographicalId() );
759 
760  //this is wrong because the FTS is defined at vertex, and it need to be properly propagated.
761  // const TrajectoryStateOnSurface initialTSOS(initialFTS, initialLayer->surface());
762 
763  const TrajectoryStateOnSurface initialTSOS = thePropagator->propagate(initialFTS,initialLayer->surface()) ;
764  if (!initialTSOS.isValid()) continue;
765 
766 #ifdef FAMOS_DEBUG
767  std::cout << "TrajectorySeedProducer: TSOS global momentum " << initialTSOS.globalMomentum() << std::endl;
768  std::cout << "\t\t\tpT = " << initialTSOS.globalMomentum().perp() << std::endl;
769  std::cout << "\t\t\teta = " << initialTSOS.globalMomentum().eta() << std::endl;
770  std::cout << "\t\t\tphi = " << initialTSOS.globalMomentum().phi() << std::endl;
771  std::cout << "TrajectorySeedProducer: TSOS local momentum " << initialTSOS.localMomentum() << std::endl;
772  std::cout << "TrajectorySeedProducer: TSOS local error " << initialTSOS.localError().positionError() << std::endl;
773  std::cout << "TrajectorySeedProducer: TSOS local error matrix " << initialTSOS.localError().matrix() << std::endl;
774  std::cout << "TrajectorySeedProducer: TSOS surface side " << initialTSOS.surfaceSide() << std::endl;
775 #endif
776  stateOnDet(initialTSOS,
777  recHits.front().geographicalId().rawId(),
778  initialState);
779  // Create a new Trajectory Seed
780  output[ialgo]->push_back(TrajectorySeed(initialState, recHits, alongMomentum));
781 #ifdef FAMOS_DEBUG
782  std::cout << "Trajectory seed created ! " << std::endl;
783 #endif
784  break;
785  // End of the loop over seeding algorithms
786  }
787  // End on the loop over simtracks
788  }
789 
790  for ( unsigned ialgo=0; ialgo<seedingAlgo.size(); ++ialgo ) {
791  std::auto_ptr<TrajectorySeedCollection> p(output[ialgo]);
792  e.put(p,seedingAlgo[ialgo]);
793  }
794 
795 }
796 
797 
798 // This is a copy of a method in
799 // TrackingTools/TrajectoryState/src/TrajectoryStateTransform.cc
800 // but it does not return a pointer (thus avoiding a memory leak)
801 // In addition, it's also CPU more efficient, because
802 // ts.localError().matrix() is not copied
803 void
805  unsigned int detid,
806  PTrajectoryStateOnDet& pts) const
807 {
808 
809  const AlgebraicSymMatrix55& m = ts.localError().matrix();
810 
811  int dim = 5;
812 
813  float localErrors[15];
814  int k = 0;
815  for (int i=0; i<dim; ++i) {
816  for (int j=0; j<=i; ++j) {
817  localErrors[k++] = m(i,j);
818  }
819  }
820  int surfaceSide = static_cast<int>(ts.surfaceSide());
821 
823  localErrors, detid,
824  surfaceSide);
825 }
826 
827 bool
829  GlobalPoint& gpos2,
830  double error,
831  bool forward,
832  unsigned algo) const {
833 
834  if ( !seedCleaning ) return true;
835 
836  // The hits 1 and 2 positions, in HepLorentzVector's
837  XYZTLorentzVector thePos1(gpos1.x(),gpos1.y(),gpos1.z(),0.);
838  XYZTLorentzVector thePos2(gpos2.x(),gpos2.y(),gpos2.z(),0.);
839 #ifdef FAMOS_DEBUG
840  std::cout << "ThePos1 = " << thePos1 << std::endl;
841  std::cout << "ThePos2 = " << thePos2 << std::endl;
842 #endif
843 
844 
845  // Create new particles that pass through the second hit with pT = ptMin
846  // and charge = +/-1
847 
848  // The momentum direction is by default joining the two hits
849  XYZTLorentzVector theMom2 = (thePos2-thePos1);
850 
851  // The corresponding RawParticle, with an (irrelevant) electric charge
852  // (The charge is determined in the next step)
853  ParticlePropagator myPart(theMom2,thePos2,1.,theFieldMap);
854 
858  bool intersect = myPart.propagateToBeamCylinder(thePos1,originRadius[algo]*1.);
859  if ( !intersect ) return false;
860 
861 #ifdef FAMOS_DEBUG
862  std::cout << "MyPart R = " << myPart.R() << "\t Z = " << myPart.Z()
863  << "\t pT = " << myPart.Pt() << std::endl;
864 #endif
865 
866  // Check if the constraints are satisfied
867  // 1. pT at cylinder with radius originRadius
868  if ( myPart.Pt() < originpTMin[algo] ) return false;
869 
870  // 2. Z compatible with beam spot size
871  if ( fabs(myPart.Z()-z0) > originHalfLength[algo] ) return false;
872 
873  // 3. Z compatible with one of the primary vertices (always the case if no primary vertex)
874  const reco::VertexCollection* theVertices = vertices[algo];
875  if (!theVertices) return true;
876  unsigned nVertices = theVertices->size();
877  if ( !nVertices || zVertexConstraint[algo] < 0. ) return true;
878  // Radii of the two hits with respect to the beam spot position
879  double R1 = std::sqrt ( (thePos1.X()-x0)*(thePos1.X()-x0)
880  + (thePos1.Y()-y0)*(thePos1.Y()-y0) );
881  double R2 = std::sqrt ( (thePos2.X()-x0)*(thePos2.X()-x0)
882  + (thePos2.Y()-y0)*(thePos2.Y()-y0) );
883  // Loop on primary vertices
884  for ( unsigned iv=0; iv<nVertices; ++iv ) {
885  // Z position of the primary vertex
886  double zV = (*theVertices)[iv].z();
887  // Constraints on the inner hit
888  double checkRZ1 = forward ?
889  (thePos1.Z()-zV+zVertexConstraint[algo]) / (thePos2.Z()-zV+zVertexConstraint[algo]) * R2 :
890  -zVertexConstraint[algo] + R1/R2*(thePos2.Z()-zV+zVertexConstraint[algo]);
891  double checkRZ2 = forward ?
892  (thePos1.Z()-zV-zVertexConstraint[algo])/(thePos2.Z()-zV-zVertexConstraint[algo]) * R2 :
893  +zVertexConstraint[algo] + R1/R2*(thePos2.Z()-zV-zVertexConstraint[algo]);
894  double checkRZmin = std::min(checkRZ1,checkRZ2)-3.*error;
895  double checkRZmax = std::max(checkRZ1,checkRZ2)+3.*error;
896  // Check if the innerhit is within bounds
897  bool compat = forward ?
898  checkRZmin < R1 && R1 < checkRZmax :
899  checkRZmin < thePos1.Z()-zV && thePos1.Z()-zV < checkRZmax;
900  // If it is, just return ok
901  if ( compat ) return compat;
902  }
903  // Otherwise, return not ok
904  return false;
905 
906 }
907 
908 
bool compatibleWithBeamAxis(GlobalPoint &gpos1, GlobalPoint &gpos2, double error, bool forward, unsigned algo) const
void setCharge(float q)
set the MEASURED charge
Definition: RawParticle.cc:139
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::vector< unsigned > minRecHits
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
double zImpactParameter(double x0=0, double y0=0.) const
Longitudinal impact parameter.
std::vector< unsigned int > firstHitSubDetectorNumber
T perp() const
Definition: PV3DBase.h:72
virtual void beginRun(edm::Run const &run, const edm::EventSetup &es) override
std::vector< std::vector< LayerSpec > > theLayersInSets
TrajectorySeedProducer(const edm::ParameterSet &conf)
const LocalTrajectoryParameters & localParameters() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
ROOT::Math::SMatrixIdentity AlgebraicMatrixID
std::vector< double > pTMin
size_type size() const
Definition: OwnVector.h:248
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
std::vector< edm::EDGetTokenT< reco::VertexCollection > > recoVertexToken
bool propagateToBeamCylinder(const XYZTLorentzVector &v, double radius=0.)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::vector< double > originRadius
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
const MagneticFieldMap * theFieldMap
void stateOnDet(const TrajectoryStateOnSurface &ts, unsigned int detid, PTrajectoryStateOnDet &pts) const
A mere copy (without memory leak) of an existing tracking method.
std::vector< double > originpTMin
uint16_t size_type
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
double charge(const std::vector< uint8_t > &Ampls)
GlobalPoint globalPosition() const
The global position.
Definition: TrackerRecHit.h:99
edm::EDGetTokenT< SiTrackerGSMatchedRecHit2DCollection > recHitToken
double largerError()
LocalError positionError() const
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
float float float z
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
bool isOnTheSameLayer(const TrackerRecHit &other) const
Check if two hits are on the same layer of the same subdetector.
std::vector< unsigned int > numberOfHits
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
LocalVector localMomentum() const
PropagatorWithMaterial * thePropagator
double R() const
vertex radius
Definition: RawParticle.h:278
static const double pts[33]
Definition: Constants.h:30
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
void push_back(D *&d)
Definition: OwnVector.h:274
std::vector< std::vector< unsigned int > > secondHitSubDetectors
virtual void produce(edm::Event &e, const edm::EventSetup &es) override
std::vector< TrajectorySeed > TrajectorySeedCollection
const T & max(const T &a, const T &b)
std::vector< std::vector< unsigned int > > thirdHitSubDetectors
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
T sqrt(T t)
Definition: SSEVec.h:48
std::vector< double > zVertexConstraint
SurfaceSide surfaceSide() const
Position relative to material, defined relative to momentum vector.
T z() const
Definition: PV3DBase.h:64
double Z() const
z of vertex
Definition: RawParticle.h:276
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
const AlgebraicSymMatrix55 & matrix() const
const math::XYZTLorentzVectorD & position() const
Definition: CoreSimVertex.h:26
const LocalTrajectoryError & localError() const
const TrackerGeometry * theGeometry
edm::EDGetTokenT< edm::SimVertexContainer > simVertexToken
virtual const GeomDet * idToDet(DetId) const
GlobalVector momentum() const
tuple conf
Definition: dbtoconf.py:185
std::vector< unsigned int > thirdHitSubDetectorNumber
int k[5][pyjets_maxn]
virtual TrackingRecHit * clone() const =0
int vertIndex() const
index of the vertex in the Event container (-1 if no vertex)
Definition: SimTrack.h:29
bool isForward() const
Is it a forward hit ?
Definition: TrackerRecHit.h:93
edm::EDGetTokenT< edm::SimTrackContainer > simTrackToken
virtual TrackingRecHit const * hit() const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const T & get() const
Definition: EventSetup.h:55
std::vector< unsigned int > secondHitSubDetectorNumber
T const * product() const
Definition: ESHandle.h:62
unsigned int layerNumber() const
The Layer Number.
Definition: TrackerRecHit.h:84
T eta() const
Definition: PV3DBase.h:76
ESHandle< TrackerGeometry > geometry
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:40
bool isOnRequestedDet(const std::vector< unsigned int > &whichDet, const std::string &seedingAlgo) const
Check if the hit is on one of the requested detector.
GlobalVector globalMomentum() const
std::vector< const reco::VertexCollection * > vertices
std::vector< double > originHalfLength
const math::XYZTLorentzVectorD & momentum() const
particle info...
Definition: CoreSimTrack.h:36
std::vector< double > maxZ0
unsigned int subDetId() const
The subdet Id.
Definition: TrackerRecHit.h:81
std::vector< std::string > seedingAlgo
tuple cout
Definition: gather_cfg.py:121
std::vector< edm::InputTag > primaryVertices
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken
std::vector< std::vector< unsigned int > > firstHitSubDetectors
std::vector< double > maxD0
DetId geographicalId() const
Definition: DDAxes.h:10
T x() const
Definition: PV3DBase.h:62
reference front()
Definition: OwnVector.h:349
Definition: Run.h:41
double xyImpactParameter(double x0=0., double y0=0.) const
Transverse impact parameter.
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15
const MagneticField * theMagField