CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PhotonConversionTrajectorySeedProducerFromSingleLegAlgo.cc
Go to the documentation of this file.
3 
4 //#define debugTSPFSLA
5 
6 inline double sqr(double a){return a*a;}
7 
10  :_conf(conf),seedCollection(0),seedCollectionOfSourceTracks(0),
11  hitsfactoryPSet(conf.getParameter<edm::ParameterSet>("OrderedHitsFactoryPSet")),
12  creatorPSet(conf.getParameter<edm::ParameterSet>("SeedCreatorPSet")),
13  regfactoryPSet(conf.getParameter<edm::ParameterSet>("RegionFactoryPSet")),
14  theClusterCheck(conf.getParameter<edm::ParameterSet>("ClusterCheckPSet")),
15  theSilentOnClusterCheck(conf.getParameter<edm::ParameterSet>("ClusterCheckPSet").getUntrackedParameter<bool>("silentClusterCheck",false)),
16  _vtxMinDoF(conf.getParameter<double>("vtxMinDoF")),
17  _maxDZSigmas(conf.getParameter<double>("maxDZSigmas")),
18  _maxNumSelVtx(conf.getParameter<uint32_t>("maxNumSelVtx")),
19  _applyTkVtxConstraint(conf.getParameter<bool>("applyTkVtxConstraint")),
20  _countSeedTracks(0),
21  _primaryVtxInputTag(conf.getParameter<edm::InputTag>("primaryVerticesTag")),
22  _beamSpotInputTag(conf.getParameter<edm::InputTag>("beamSpotInputTag"))
23 {
24  init();
25 }
26 
30  delete theHitsGenerator;
31  if(theSeedCreator!=NULL)
32  delete theSeedCreator;
34  delete theRegionProducer;
35 }
36 
38 init(){
42 }
43 
46 
47  myEsetup = &setup;
48  myEvent = &event;
49 
50  if(seedCollection!=0)
51  delete seedCollection;
52 
55 
58 
59  size_t clustsOrZero = theClusterCheck.tooManyClusters(event);
60  if (clustsOrZero){
62  edm::LogError("TooManyClusters") << "Found too many clusters (" << clustsOrZero << "), bailing out.\n";
63  return ;
64  }
65 
66 
67  edm::ESHandle<MagneticField> handleMagField;
68  setup.get<IdealMagneticFieldRecord>().get(handleMagField);
69  magField = handleMagField.product();
70  if (unlikely(magField->inTesla(GlobalPoint(0.,0.,0.)).z()<0.01)) return;
71 
73 
74 
75  event.getByLabel(_primaryVtxInputTag, vertexHandle);
76  if (!vertexHandle.isValid() || vertexHandle->empty()){
77  edm::LogError("PhotonConversionFinderFromTracks") << "Error! Can't get the product primary Vertex Collection "<< _primaryVtxInputTag << "\n";
78  return;
79  }
80 
81  event.getByLabel(_beamSpotInputTag,recoBeamSpotHandle);
82 
83 
84  regions = theRegionProducer->regions(event,setup);
85 
86  //Do the analysis
87  loopOnTracks();
88 
89 
90 #ifdef debugTSPFSLA
91  std::stringstream ss;
92  ss.str("");
93  ss << "\n++++++++++++++++++\n";
94  ss << "seed collection size " << seedCollection->size();
95  BOOST_FOREACH(TrajectorySeed tjS,*seedCollection){
96  po.print(ss, tjS);
97  }
98  edm::LogInfo("debugTrajSeedFromSingleLeg") << ss.str();
99  //-------------------------------------------------
100 #endif
101 
102  // clear memory
104  for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) delete (*ir);
105 
106 }
107 
108 
111 
112  //--- Get Tracks
114 
115  if(trackCollectionH.isValid()==0){
116  edm::LogError("MissingInput")<<" could not find track collecion:"<<_conf.getParameter<edm::InputTag>("TrackRefitter");
117  return;
118  }
119  size_t idx=0, sel=0;
121 
122  ss.str("");
123 
124  for( reco::TrackCollection::const_iterator tr = trackCollectionH->begin();
125  tr != trackCollectionH->end(); tr++, idx++) {
126 
127  // #ifdef debugTSPFSLA
128  // ss << "\nStuding track Nb " << idx;
129  // #endif
130 
131  if(rejectTrack(*tr)) continue;
132  std::vector<reco::Vertex> selectedPriVtxCompatibleWithTrack;
134  selectedPriVtxCompatibleWithTrack.push_back(*(vertexHandle->begin())); //Same approach as before
135  }else{
136  if(!selectPriVtxCompatibleWithTrack(*tr,selectedPriVtxCompatibleWithTrack)) continue;
137  }
138 
139  sel++;
140  loopOnPriVtx(*tr,selectedPriVtxCompatibleWithTrack);
141  }
142 #ifdef debugTSPFSLA
143  edm::LogInfo("debugTrajSeedFromSingleLeg") << ss.str();
144  edm::LogInfo("debugTrajSeedFromSingleLeg") << "Inspected " << sel << " tracks over " << idx << " tracks. \t # tracks providing at least one seed " << _countSeedTracks ;
145 #endif
146 }
147 
149 selectPriVtxCompatibleWithTrack(const reco::Track& tk, std::vector<reco::Vertex>& selectedPriVtxCompatibleWithTrack){
150 
151  std::vector< std::pair< double, short> > idx;
152  short count=-1;
153 
154  double cosPhi=tk.px()/tk.pt();
155  double sinPhi=tk.py()/tk.pt();
156  double sphi2=tk.covariance(2,2);
157  double stheta2=tk.covariance(1,1);
158 
159  BOOST_FOREACH(const reco::Vertex& vtx, *vertexHandle){
160  count++;
161  if(vtx.ndof()<= _vtxMinDoF) continue;
162 
163  double _dz= tk.dz(vtx.position());
164  double _dzError=tk.dzError();
165 
166  double cotTheta=tk.pz()/tk.pt();
167  double dx = vtx.position().x();
168  double dy = vtx.position().y();
169  double sx2=vtx.covariance(0,0);
170  double sy2=vtx.covariance(1,1);
171 
172  double sxy2= sqr(cosPhi*cotTheta)*sx2+
173  sqr(sinPhi*cotTheta)*sy2+
174  sqr(cotTheta*(-dx*sinPhi+dy*cosPhi))*sphi2+
175  sqr((1+cotTheta*cotTheta)*(dx*cosPhi+dy*sinPhi))*stheta2;
176 
177  _dzError=sqrt(_dzError*_dzError+vtx.covariance(2,2)+sxy2); //there is a missing component, related to the element (vtx.x*px+vtx.y*py)/pt * pz/pt. since the tk ref point is at the point of closest approach, this scalar product should be almost zero.
178 
179 #ifdef debugTSPFSLA
180  ss << " primary vtx " << vtx.position() << " \tk vz " << tk.vz() << " vx " << tk.vx() << " vy " << tk.vy() << " pz/pt " << tk.pz()/tk.pt() << " \t dz " << _dz << " \t " << _dzError << " sxy2 "<< sxy2<< " \t dz/dzErr " << _dz/_dzError<< std::endl;
181 #endif
182 
183  if(fabs(_dz)/_dzError > _maxDZSigmas) continue;
184 
185  idx.push_back(std::pair<double,short>(fabs(_dz),count));
186  }
187  if(idx.size()==0) {
188 #ifdef debugTSPFSLA
189  ss << "no vertex selected " << std::endl;
190 #endif
191  return false;
192 }
193 
194  std::stable_sort(idx.begin(),idx.end(),lt_);
195  for(size_t i=0;i<_maxNumSelVtx && i<idx.size();++i){
196  selectedPriVtxCompatibleWithTrack.push_back((*vertexHandle)[idx[i].second]);
197 #ifdef debugTSPFSLA
198  ss << "selected vtx dz " << idx[0].first << " position" << idx[0].second << std::endl;
199 #endif
200  }
201 
202  return true;
203 }
204 
206 loopOnPriVtx(const reco::Track& tk, const std::vector<reco::Vertex>& selectedPriVtxCompatibleWithTrack){
207 
208  bool foundAtLeastASeedCand=false;
209  BOOST_FOREACH(const reco::Vertex vtx, selectedPriVtxCompatibleWithTrack){
210 
211  math::XYZPoint primaryVertexPoint=math::XYZPoint(vtx.position());
212 
213  for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) {
214  const TrackingRegion & region = **ir;
215 
216 #ifdef debugTSPFSLA
217  ss << "[PrintRegion] " << region.print() << std::endl;
218 #endif
219 
220  //This if is just for the _countSeedTracks. otherwise
221  //inspectTrack(&tk,region, primaryVertexPoint);
222  //would be enough
223 
224  if(
225  inspectTrack(&tk,region, primaryVertexPoint)
226  and
227  !foundAtLeastASeedCand
228  ){
229  foundAtLeastASeedCand=true;
230  _countSeedTracks++;
231  }
232 
233  }
234  }
235 }
236 
238 rejectTrack(const reco::Track& track){
239 
242  beamSpot = math::XYZVector(recoBeamSpotHandle->position());
243 
244  _IdealHelixParameters.setData(&track,beamSpot);
246  //this case means a null results on the _IdealHelixParameters side
247  return true;
248  }
249 
250  float rMin=2.; //cm
251  if(_IdealHelixParameters.GetTangentPoint().rho()<rMin){
252  //this case means a track that has the tangent point nearby the primary vertex
253  // if the track is primary, this number tends to be the primary vertex itself
254  //Rejecting all the potential photon conversions having a "vertex" inside the beampipe
255  //We should not miss too much, seen that the conversions at the beam pipe are the better reconstructed
256  return true;
257  }
258  }
259 
260  //-------------------------------------------------------
261  /*
262  float maxPt2=64.; //Cut on pt^2 Indeed doesn't do almost nothing
263  if(track.momentum().Perp2() > maxPt2)
264  return true;
265  */
266  //-------------------------------------------------------
267  //Cut in the barrel eta region FIXME: to be extended to endcaps
268  /*
269  float maxEta=1.3;
270  if(fabs(track.eta()) > maxEta)
271  return true;
272  */
273  //-------------------------------------------------------
274  //Reject tracks that have a first valid hit in the pixel barrel/endcap layer/disk 1
275  //assume that the hits are aligned along momentum
276  /*
277  const reco::HitPattern& p=track.hitPattern();
278  for (int i=0; i<p.numberOfHits(); i++) {
279  uint32_t hit = p.getHitPattern(i);
280  // if the hit is valid and in pixel barrel, print out the layer
281  if (! p.validHitFilter(hit) ) continue;
282  if( (p.pixelBarrelHitFilter(hit) || p.pixelEndcapHitFilter(hit))
283  &&
284  p.getLayer(hit) == 1
285  )
286  return true;
287  else
288  break; //because the first valid hit is in a different layer
289  }
290  */
291  //-------------------------------------------------------
292 
293 
294  return false;
295 }
296 
297 
299 inspectTrack(const reco::Track* track, const TrackingRegion & region, math::XYZPoint& primaryVertexPoint){
300 
301  _IdealHelixParameters.setData(track,primaryVertexPoint);
302 
305  //this case means a null results on the _IdealHelixParameters side
306  return false;
307  }
308 
309  float rMin=3.; //cm
310  if(_IdealHelixParameters.GetTangentPoint().rho()<rMin){
311  //this case means a track that has the tangent point nearby the primary vertex
312  // if the track is primary, this number tends to be the primary vertex itself
313  //Rejecting all the potential photon conversions having a "vertex" inside the beampipe
314  //We should not miss too much, seen that the conversions at the beam pipe are the better reconstructed
315  return false;
316  }
317 
318  float ptmin = 0.5;
319  float originRBound = 3;
320  float originZBound = 3.;
321 
322  GlobalPoint originPos;
326  );
327  float cotTheta;
330  }else{
332  cotTheta=99999.f;
333  else
334  cotTheta=-99999.f;
335  }
336  GlobalVector originBounds(originRBound,originRBound,originZBound);
337 
338  GlobalPoint pvtxPoint(primaryVertexPoint.x(),
339  primaryVertexPoint.y(),
340  primaryVertexPoint.z()
341  );
342  ConversionRegion convRegion(originPos, pvtxPoint, cotTheta, track->thetaError(), -1*track->charge());
343 
344 #ifdef debugTSPFSLA
345  ss << "\nConversion Point " << originPos << " " << originPos.perp() << "\n";
346 #endif
347 
348  const OrderedSeedingHits & hitss = theHitsGenerator->run(convRegion, region, *myEvent, *myEsetup);
349 
350  unsigned int nHitss = hitss.size();
351 
352  if(nHitss==0)
353  return false;
354 
355 #ifdef debugTSPFSLA
356  ss << "\n nHitss " << nHitss << "\n";
357 #endif
358 
359  if (seedCollection->empty()) seedCollection->reserve(nHitss); // don't do multiple reserves in the case of multiple regions: it would make things even worse
360  // as it will cause N re-allocations instead of the normal log(N)/log(2)
361  for (unsigned int iHits = 0; iHits < nHitss; ++iHits) {
362 
363 #ifdef debugTSPFSLA
364  ss << "\n iHits " << iHits << "\n";
365 #endif
366  const SeedingHitSet & hits = hitss[iHits];
367  //if (!theComparitor || theComparitor->compatible( hits, es) ) {
368  //try{
369  theSeedCreator->trajectorySeed(*seedCollection,hits, originPos, originBounds, ptmin, *myEsetup,convRegion.cotTheta(),ss);
370  //}catch(cms::Exception& er){
371  // edm::LogError("SeedingConversion") << " Problem in the Single Leg Seed creator " <<er.what()<<std::endl;
372  //}catch(std::exception& er){
373  // edm::LogError("SeedingConversion") << " Problem in the Single Leg Seed creator " << er.what()<<std::endl;
374  //}
375  }
376  return true;
377 }
378 
379 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
T perp() const
Definition: PV3DBase.h:71
void setMagnField(const MagneticField *magnField)
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
#define abs(x)
Definition: mlp_lapack.h:159
#define NULL
Definition: scimark2.h:8
virtual unsigned int size() const =0
virtual const TrajectorySeed * trajectorySeed(TrajectorySeedCollection &seedCollection, const SeedingHitSet &hits, const GlobalPoint &vertex, const GlobalVector &vertexBounds, float ptmin, const edm::EventSetup &es, float cotTheta, std::stringstream &ss)
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:110
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:133
const Point & position() const
position
Definition: Vertex.h:93
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:182
math::XYZVector GetMomentumAtTangentPoint() const
U second(std::pair< T, U > const &p)
#define unlikely(x)
Definition: Likely.h:21
const OrderedHitPairs & run(const ConversionRegion &convRegion, const TrackingRegion &region, const edm::Event &ev, const edm::EventSetup &es)
math::XYZVector GetTangentPoint() const
std::vector< TrajectorySeed > TrajectorySeedCollection
bool isnan(float x)
Definition: math.h:13
T sqrt(T t)
Definition: SSEVec.h:46
double pt() const
track transverse momentum
Definition: TrackBase.h:131
T z() const
Definition: PV3DBase.h:63
size_t tooManyClusters(const edm::Event &e) const
double f[11][100]
void print(std::stringstream &ss, const SiStripCluster &clus)
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool isValid() const
Definition: HandleBase.h:76
bool inspectTrack(const reco::Track *track, const TrackingRegion &region, math::XYZPoint &primaryVertexPoint)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
double ndof() const
Definition: Vertex.h:89
tuple conf
Definition: dbtoconf.py:185
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:137
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:127
double dzError() const
error on dz
Definition: TrackBase.h:215
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:147
void loopOnPriVtx(const reco::Track &tk, const std::vector< reco::Vertex > &selectedPriVtxCompatibleWithTrack)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
double ptmin
Definition: HydjetWrapper.h:86
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:145
virtual std::vector< TrackingRegion * > regions(const edm::Event &ev, const edm::EventSetup &) const
double a
Definition: hdecay.h:121
Square< F >::type sqr(const F &f)
Definition: Square.h:13
bool lt_(std::pair< double, short > a, std::pair< double, short > b)
virtual std::string print() const =0
int charge() const
track electric charge
Definition: TrackBase.h:113
void setData(const reco::Track *track, math::XYZVector refPoint=math::XYZVector(0, 0, 0))
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
bool selectPriVtxCompatibleWithTrack(const reco::Track &tk, std::vector< reco::Vertex > &selectedPriVtxCompatibleWithTrack)
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:135
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:143
double thetaError() const
error on theta
Definition: TrackBase.h:201