Go to the documentation of this file.00001 #include "RecoParticleFlow/PFTracking/interface/PFDisplacedVertexCandidateFinder.h"
00002
00003 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00004 #include "DataFormats/TrackReco/interface/Track.h"
00005
00006 #include "MagneticField/Engine/interface/MagneticField.h"
00007
00008 using namespace std;
00009 using namespace reco;
00010
00011
00012
00013
00014 PFDisplacedVertexCandidateFinder::PFDisplacedVertexCandidateFinder() :
00015 vertexCandidates_( new PFDisplacedVertexCandidateCollection ),
00016 dcaCut_(1000),
00017 primaryVertexCut2_(0.0),
00018 dcaPInnerHitCut2_(1000.0),
00019 vertexCandidatesSize_(50),
00020 debug_(false) {
00021
00022 TwoTrackMinimumDistance theMinimum(TwoTrackMinimumDistance::SlowMode);
00023 theMinimum_ = theMinimum;
00024 }
00025
00026
00027 PFDisplacedVertexCandidateFinder::~PFDisplacedVertexCandidateFinder() {
00028
00029 #ifdef PFLOW_DEBUG
00030 if(debug_)
00031 cout<<"~PFDisplacedVertexCandidateFinder - number of remaining elements: "
00032 <<eventTracks_.size()<<endl;
00033 #endif
00034
00035 }
00036
00037
00038
00039 void
00040 PFDisplacedVertexCandidateFinder::setInput(const edm::Handle<TrackCollection>& trackh,
00041 const MagneticField* magField) {
00042
00043 magField_ = magField;
00044 trackMask_.clear();
00045 trackMask_.resize(trackh->size());
00046 eventTrackTrajectories_.clear();
00047 eventTrackTrajectories_.resize(trackh->size());
00048
00049 for(unsigned i=0;i<trackMask_.size(); i++) trackMask_[i] = true;
00050
00051 eventTracks_.clear();
00052 if(trackh.isValid()) {
00053 for(unsigned i=0;i<trackh->size(); i++) {
00054
00055 TrackRef tref( trackh, i);
00056 TrackBaseRef tbref(tref);
00057
00058 if( !isSelected( tbref ) ) {
00059 trackMask_[i] = false;
00060 continue;
00061 }
00062
00063 const Track* trk = tref.get();
00064 eventTracks_.push_back( tbref );
00065 eventTrackTrajectories_[i] = getGlobalTrajectoryParameters(trk);
00066 }
00067 }
00068 }
00069
00070
00071
00072
00073 void
00074 PFDisplacedVertexCandidateFinder::findDisplacedVertexCandidates() {
00075
00076 if (debug_) cout << "========= Start Finding Displaced Vertex Candidates =========" << endl;
00077
00078
00079 if(vertexCandidates_.get() )vertexCandidates_->clear();
00080 else
00081 vertexCandidates_.reset( new PFDisplacedVertexCandidateCollection );
00082
00083 vertexCandidates_->reserve(vertexCandidatesSize_);
00084 for(IE ie = eventTracks_.begin();
00085 ie != eventTracks_.end();) {
00086
00087
00088
00089
00090 PFDisplacedVertexCandidate tempVertexCandidate;
00091
00092 ie = associate( eventTracks_.end() , ie, tempVertexCandidate);
00093
00094
00095 if(tempVertexCandidate.isValid()) {
00096 packLinks( tempVertexCandidate );
00097 vertexCandidates_->push_back( tempVertexCandidate );
00098 }
00099 }
00100
00101 if(debug_) cout << "========= End Finding Displaced Vertex Candidates =========" << endl;
00102
00103 }
00104
00105
00106
00107
00108 PFDisplacedVertexCandidateFinder::IE
00109 PFDisplacedVertexCandidateFinder::associate(IE last, IE next,
00110 PFDisplacedVertexCandidate& tempVertexCandidate) {
00111
00112
00113 #ifdef PFLOW_DEBUG
00114 if(debug_ ) cout<<"== Start the association procedure =="<<endl;
00115 #endif
00116
00117 if( last!= eventTracks_.end() ) {
00118 double dist = -1;
00119 GlobalPoint P(0,0,0);
00120 PFDisplacedVertexCandidate::VertexLinkTest linktest;
00121 link( (*last), (*next), dist, P, linktest);
00122
00123
00124 if(dist<-0.5) {
00125 #ifdef PFLOW_DEBUG
00126 if(debug_ ) cout<<"link failed"<<endl;
00127 #endif
00128 return ++next;
00129 }
00130 else {
00131
00132 tempVertexCandidate.addElement( (*next));
00133 trackMask_[(*next).key()] = false;
00134 #ifdef PFLOW_DEBUG
00135 if(debug_ )
00136 cout<<"link parameters "
00137 << " *next = " << (*next).key()
00138 << " *last = " << (*last).key()
00139 << " dist = " << dist
00140 << " P.x = " << P.x()
00141 << " P.y = " << P.y()
00142 << " P.z = " << P.z()
00143 << endl;
00144 #endif
00145 }
00146 }
00147 else {
00148
00149 #ifdef PFLOW_DEBUG
00150 if(debug_ ) cout<<"adding to block element "
00151 << (*next).key()
00152 <<endl;
00153 #endif
00154 tempVertexCandidate.addElement( (*next));
00155 trackMask_[(*next).key()] = false;
00156
00157
00158 }
00159
00160
00161 #ifdef PFLOW_DEBUG
00162 if(debug_ ) {
00163 for(unsigned i=0; i<trackMask_.size(); i++)
00164 cout << " Mask[" << i << "] = " << trackMask_[i];
00165 cout << "" << endl;
00166 }
00167 #endif
00168
00169 for(IE ie = eventTracks_.begin();
00170 ie != eventTracks_.end();) {
00171
00172 if( ie == last || ie == next ) { ++ie; continue;}
00173
00174
00175 if( !trackMask_[(*ie).key()] ) { ++ie; continue;}
00176
00177 #ifdef PFLOW_DEBUG
00178 if(debug_ ) cout<<"calling associate "
00179 << (*next).key()
00180 <<" & "
00181 << (*ie).key()
00182 <<endl;
00183 #endif
00184 ie = associate(next, ie, tempVertexCandidate);
00185
00186 }
00187
00188 #ifdef PFLOW_DEBUG
00189 if(debug_ ) {
00190 cout<<"**** removing element "<<endl;
00191 }
00192 #endif
00193
00194 IE iteratorToNextFreeElement = eventTracks_.erase( next );
00195
00196 #ifdef PFLOW_DEBUG
00197 if(debug_ ) cout<< "== End the association procedure ==" <<endl;
00198 #endif
00199
00200 return iteratorToNextFreeElement;
00201 }
00202
00203
00204
00205 void
00206 PFDisplacedVertexCandidateFinder::link(const TrackBaseRef& el1,
00207 const TrackBaseRef& el2,
00208 double& dist, GlobalPoint& P,
00209 PFDisplacedVertexCandidate::VertexLinkTest& vertexLinkTest
00210 ){
00211
00212 if ( fabs(el1->eta()-el2->eta()) > 1) {dist = -1; return;}
00213 if ( el1->pt()>2 && el2->pt()>2 && fabs(el1->phi()-el2->phi()) > 1) {dist = -1; return;}
00214
00215 const GlobalTrajectoryParameters& gt1 = eventTrackTrajectories_[el1.key()];
00216 const GlobalTrajectoryParameters& gt2 = eventTrackTrajectories_[el2.key()];
00217
00218
00219
00220 theMinimum_.calculate(gt1,gt2);
00221
00222
00223 dist = theMinimum_.distance();
00224 P = theMinimum_.crossingPoint();
00225
00226 vertexLinkTest = PFDisplacedVertexCandidate::LINKTEST_DCA;
00227
00228
00229
00230 if (dist > dcaCut_) {dist = -1; return;}
00231
00232
00233 double rho2 = P.x()*P.x()+P.y()*P.y();
00234
00235 if ( rho2 < primaryVertexCut2_) {dist = -1; return;}
00236
00237
00238
00239 double tob_rho_limit2 = 10000;
00240 double tec_z_limit = 270;
00241
00242 if ( rho2 > tob_rho_limit2) {dist = -1; return;}
00243 if ( fabs(P.z()) > tec_z_limit) {dist = -1; return;}
00244
00245 return;
00246
00247
00248
00249
00250
00251
00252
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271
00272
00273 }
00274
00275
00276
00277
00278 void
00279 PFDisplacedVertexCandidateFinder::packLinks( PFDisplacedVertexCandidate& vertexCandidate) {
00280
00281
00282 const vector < TrackBaseRef >& els = vertexCandidate.elements();
00283
00284
00285 for( unsigned i1=0; i1<els.size(); i1++ ) {
00286 for( unsigned i2=i1+1; i2<els.size(); i2++ ) {
00287
00288
00289 if( i1==i2 ) continue;
00290
00291 double dist = -1;
00292 GlobalPoint P(0,0,0);
00293 PFDisplacedVertexCandidate::VertexLinkTest linktest;
00294
00295 link( els[i1], els[i2], dist, P, linktest);
00296
00297
00298 #ifdef PFLOW_DEBUG
00299 if( debug_ )
00300 cout << "Setting link between elements " << i1 << " key " << els[i1].key()
00301 << " and " << i2 << " key " << els[i2].key()
00302 << " of dist =" << dist << " computed from link test "
00303 << linktest << endl;
00304 #endif
00305
00306 if(dist >-0.5) vertexCandidate.setLink( i1, i2, dist, P, linktest );
00307 }
00308 }
00309
00310 }
00311
00312
00313
00314
00315
00316
00317
00318
00319 GlobalTrajectoryParameters
00320 PFDisplacedVertexCandidateFinder::getGlobalTrajectoryParameters
00321 (const Track* track) const
00322 {
00323
00324 const GlobalPoint position(track->vx(),
00325 track->vy(),
00326 track->vz());
00327
00328 const GlobalVector momentum(track->momentum().x(),
00329 track->momentum().y(),
00330 track->momentum().z());
00331
00332 GlobalTrajectoryParameters gtp(position,momentum,
00333 track->charge(),magField_);
00334
00335 return gtp;
00336 }
00337
00338
00339
00340 bool
00341 PFDisplacedVertexCandidateFinder::goodPtResolution( const TrackBaseRef& trackref) const {
00342
00343 double nChi2 = trackref->normalizedChi2();
00344 double pt = trackref->pt();
00345 double dpt = trackref->ptError();
00346
00347 double pt_error = dpt/pt*100;
00348
00349 if (debug_) cout << " PFDisplacedVertexFinder: PFrecTrack->Track Pt= "
00350 << pt << " dPt/Pt = " << pt_error << "% nChi2 = " << nChi2 << endl;
00351 if (nChi2 > nChi2_max_ || pt < pt_min_){
00352
00353 if (debug_) cout << " PFBlockAlgo: skip badly measured or low pt track"
00354 << " nChi2_cut = " << 5
00355 << " pt_cut = " << 0.2 << endl;
00356 return false;
00357 }
00358
00359 return true;
00360 }
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372 ostream& operator<<(std::ostream& out, const PFDisplacedVertexCandidateFinder& a) {
00373 if(! out) return out;
00374
00375 out<<"====== Particle Flow Block Algorithm ======= ";
00376 out<<endl;
00377 out<<"number of unassociated elements : "<<a.eventTracks_.size()<<endl;
00378 out<<endl;
00379
00380 for(PFDisplacedVertexCandidateFinder::IEC ie = a.eventTracks_.begin();
00381 ie != a.eventTracks_.end(); ie++) {
00382
00383 double pt = (*ie).get()->pt();
00384
00385 math::XYZPoint Pi = (*ie).get()->innerPosition();
00386 math::XYZPoint Po = (*ie).get()->outerPosition();
00387
00388 double innermost_radius = sqrt(Pi.x()*Pi.x() + Pi.y()*Pi.y() + Pi.z()*Pi.z());
00389 double outermost_radius = sqrt(Po.x()*Po.x() + Po.y()*Po.y() + Po.z()*Po.z());
00390 double innermost_rho = sqrt(Pi.x()*Pi.x() + Pi.y()*Pi.y());
00391 double outermost_rho = sqrt(Po.x()*Po.x() + Po.y()*Po.y());
00392
00393 out<<"ie = " << (*ie).key()
00394 <<" pt = " << pt
00395 <<" innermost hit radius = " << innermost_radius << " rho = " << innermost_rho
00396 <<" outermost hit radius = " << outermost_radius << " rho = " << outermost_rho
00397 <<endl;
00398 }
00399
00400
00401 const std::auto_ptr< reco::PFDisplacedVertexCandidateCollection >& vertexCandidates
00402 = a.vertexCandidates();
00403
00404 if(!vertexCandidates.get() ) {
00405 out<<"vertexCandidates already transfered"<<endl;
00406 }
00407 else {
00408 out<<"number of vertexCandidates : "<<vertexCandidates->size()<<endl;
00409 out<<endl;
00410
00411
00412 for(PFDisplacedVertexCandidateFinder::IBC ib=vertexCandidates->begin();
00413 ib != vertexCandidates->end(); ib++)
00414 ib->Dump();
00415
00416
00417
00418 }
00419
00420 return out;
00421 }