00001
00002 #include "CommonTools/RecoUtils/interface/PF_PU_AssoMapAlgos.h"
00003
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005
00006 #include "DataFormats/Math/interface/deltaR.h"
00007
00008 #include "DataFormats/GeometryCommonDetAlgo/interface/Measurement1D.h"
00009 #include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
00010 #include "RecoVertex/VertexPrimitives/interface/ConvertToFromReco.h"
00011 #include "RecoVertex/VertexPrimitives/interface/VertexState.h"
00012
00013 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
00014
00015 #include "TrackingTools/IPTools/interface/IPTools.h"
00016 #include "RecoEgamma/EgammaTools/interface/ConversionTools.h"
00017
00018 using namespace edm;
00019 using namespace std;
00020 using namespace reco;
00021
00022
00023
00024
00025
00026 PF_PU_AssoMapAlgos::PF_PU_AssoMapAlgos(const edm::ParameterSet& iConfig)
00027 : maxNumWarnings_(3),
00028 numWarnings_(0)
00029 {
00030
00031 input_MaxNumAssociations_ = iConfig.getParameter<int>("MaxNumberOfAssociations");
00032
00033 input_VertexCollection_= iConfig.getParameter<InputTag>("VertexCollection");
00034
00035 input_BeamSpot_= iConfig.getParameter<InputTag>("BeamSpot");
00036
00037 input_doReassociation_= iConfig.getParameter<bool>("doReassociation");
00038 cleanedColls_ = iConfig.getParameter<bool>("GetCleanedCollections");
00039
00040 ConversionsCollection_= iConfig.getParameter<InputTag>("ConversionsCollection");
00041
00042 KshortCollection_= iConfig.getParameter<InputTag>("V0KshortCollection");
00043 LambdaCollection_= iConfig.getParameter<InputTag>("V0LambdaCollection");
00044
00045 NIVertexCollection_= iConfig.getParameter<InputTag>("NIVertexCollection");
00046
00047 input_FinalAssociation_= iConfig.getUntrackedParameter<int>("FinalAssociation", 0);
00048
00049 ignoremissingpfcollection_ = iConfig.getParameter<bool>("ignoreMissingCollection");
00050
00051 input_nTrack_ = iConfig.getParameter<double>("nTrackWeight");
00052
00053 }
00054
00055
00056
00057
00058
00059 void
00060 PF_PU_AssoMapAlgos::GetInputCollections(edm::Event& iEvent, const edm::EventSetup& iSetup)
00061 {
00062
00063
00064 iEvent.getByLabel(input_BeamSpot_, beamspotH);
00065
00066
00067 iEvent.getByLabel(ConversionsCollection_, convCollH);
00068 cleanedConvCollP = PF_PU_AssoMapAlgos::GetCleanedConversions(convCollH,beamspotH,cleanedColls_);
00069
00070
00071 iEvent.getByLabel(KshortCollection_, vertCompCandCollKshortH);
00072 cleanedKshortCollP = PF_PU_AssoMapAlgos::GetCleanedKshort(vertCompCandCollKshortH,beamspotH,cleanedColls_);
00073
00074
00075 iEvent.getByLabel(LambdaCollection_, vertCompCandCollLambdaH);
00076 cleanedLambdaCollP = PF_PU_AssoMapAlgos::GetCleanedLambda(vertCompCandCollLambdaH,beamspotH,cleanedColls_);
00077
00078
00079
00080 missingColls = false;
00081 if(!iEvent.getByLabel(NIVertexCollection_,displVertexCollH)){
00082 if (ignoremissingpfcollection_){
00083
00084 missingColls = true;
00085
00086 if ( numWarnings_ < maxNumWarnings_ ) {
00087 LogWarning("PF_PU_AssoMapAlgos::GetInputCollections")
00088 << "No Extra objects available in input file --> skipping reconstruction of displaced vertices !!" << endl;
00089 ++numWarnings_;
00090 }
00091
00092 }
00093 } else {
00094
00095 cleanedNICollP = PF_PU_AssoMapAlgos::GetCleanedNI(displVertexCollH,beamspotH,true);
00096
00097 }
00098
00099
00100 iEvent.getByLabel(input_VertexCollection_, vtxcollH);
00101
00102 iSetup.get<IdealMagneticFieldRecord>().get(bFieldH);
00103
00104 }
00105
00106
00107
00108
00109 std::auto_ptr<TrackToVertexAssMap>
00110 PF_PU_AssoMapAlgos::CreateTrackToVertexMap(edm::Handle<reco::TrackCollection> trkcollH, const edm::EventSetup& iSetup)
00111 {
00112
00113 auto_ptr<TrackToVertexAssMap> track2vertex(new TrackToVertexAssMap());
00114
00115 int num_vertices = vtxcollH->size();
00116 if ( num_vertices < input_MaxNumAssociations_) input_MaxNumAssociations_ = num_vertices;
00117
00118
00119 for ( size_t idxTrack = 0; idxTrack < trkcollH->size(); ++idxTrack ) {
00120
00121 TrackRef trackref = TrackRef(trkcollH, idxTrack);
00122
00123 TransientTrack transtrk(trackref, &(*bFieldH) );
00124 transtrk.setBeamSpot(*beamspotH);
00125 transtrk.setES(iSetup);
00126
00127 vector<VertexRef>* vtxColl_help = CreateVertexVector(vtxcollH);
00128
00129 for ( int assoc_ite = 0; assoc_ite < input_MaxNumAssociations_; ++assoc_ite ) {
00130
00131 VertexStepPair assocVtx = FindAssociation(trackref, vtxColl_help, bFieldH, iSetup, beamspotH, assoc_ite);
00132 int step = assocVtx.second;
00133 double distance = ( IPTools::absoluteImpactParameter3D( transtrk, *(assocVtx.first) ) ).second.value();
00134
00135 int quality = DefineQuality(assoc_ite, step, distance);
00136
00137
00138
00139
00140
00141
00142
00143 track2vertex->insert( assocVtx.first, make_pair(trackref, quality) );
00144
00145 PF_PU_AssoMapAlgos::EraseVertex(vtxColl_help, assocVtx.first);
00146
00147 }
00148
00149 delete vtxColl_help;
00150
00151 }
00152
00153 return track2vertex;
00154
00155 }
00156
00157
00158
00159
00160
00161 std::auto_ptr<VertexToTrackAssMap>
00162 PF_PU_AssoMapAlgos::CreateVertexToTrackMap(edm::Handle<reco::TrackCollection> trkcollH, const edm::EventSetup& iSetup)
00163 {
00164
00165 auto_ptr<VertexToTrackAssMap> vertex2track(new VertexToTrackAssMap());
00166
00167 int num_vertices = vtxcollH->size();
00168 if ( num_vertices < input_MaxNumAssociations_) input_MaxNumAssociations_ = num_vertices;
00169
00170
00171 for ( size_t idxTrack = 0; idxTrack < trkcollH->size(); ++idxTrack ) {
00172
00173 TrackRef trackref = TrackRef(trkcollH, idxTrack);
00174
00175 TransientTrack transtrk(trackref, &(*bFieldH) );
00176 transtrk.setBeamSpot(*beamspotH);
00177 transtrk.setES(iSetup);
00178
00179 vector<VertexRef>* vtxColl_help = CreateVertexVector(vtxcollH);
00180
00181 for ( int assoc_ite = 0; assoc_ite < input_MaxNumAssociations_; ++assoc_ite ) {
00182
00183 VertexStepPair assocVtx = FindAssociation(trackref, vtxColl_help, bFieldH, iSetup, beamspotH, assoc_ite);
00184 int step = assocVtx.second;
00185 double distance = ( IPTools::absoluteImpactParameter3D( transtrk, *(assocVtx.first) ) ).second.value();
00186
00187 int quality = DefineQuality(assoc_ite, step, distance);
00188
00189
00190
00191
00192
00193
00194 vertex2track->insert( trackref, make_pair(assocVtx.first, quality) );
00195
00196 PF_PU_AssoMapAlgos::EraseVertex(vtxColl_help, assocVtx.first);
00197
00198 }
00199
00200 delete vtxColl_help;
00201
00202 }
00203
00204 return vertex2track;
00205
00206 }
00207
00208
00209
00210
00211
00212 auto_ptr<TrackToVertexAssMap>
00213 PF_PU_AssoMapAlgos::SortAssociationMap(TrackToVertexAssMap* trackvertexassInput)
00214 {
00215
00216 auto_ptr<TrackToVertexAssMap> trackvertexassOutput(new TrackToVertexAssMap() );
00217
00218
00219 VertexPtsumVector vertexptsumvector;
00220
00221
00222 for(TrackToVertexAssMap::const_iterator assomap_ite=trackvertexassInput->begin(); assomap_ite!=trackvertexassInput->end(); assomap_ite++){
00223
00224 const VertexRef assomap_vertexref = assomap_ite->key;
00225 const TrackQualityPairVector trckcoll = assomap_ite->val;
00226
00227 float ptsum = 0;
00228
00229 TrackRef trackref;
00230
00231
00232 for(unsigned int trckcoll_ite=0; trckcoll_ite<trckcoll.size(); trckcoll_ite++){
00233
00234 trackref = trckcoll[trckcoll_ite].first;
00235 int quality = trckcoll[trckcoll_ite].second;
00236
00237 if ( quality<=2 ) continue;
00238
00239 double man_pT = trackref->pt() - trackref->ptError();
00240 if(man_pT>0.) ptsum+=man_pT*man_pT;
00241
00242 }
00243
00244 vertexptsumvector.push_back(make_pair(assomap_vertexref,ptsum));
00245
00246 }
00247
00248 while (vertexptsumvector.size()!=0){
00249
00250 VertexRef vertexref_highestpT;
00251 float highestpT = 0.;
00252 int highestpT_index = 0;
00253
00254 for(unsigned int vtxptsumvec_ite=0; vtxptsumvec_ite<vertexptsumvector.size(); vtxptsumvec_ite++){
00255
00256 if(vertexptsumvector[vtxptsumvec_ite].second>highestpT){
00257
00258 vertexref_highestpT = vertexptsumvector[vtxptsumvec_ite].first;
00259 highestpT = vertexptsumvector[vtxptsumvec_ite].second;
00260 highestpT_index = vtxptsumvec_ite;
00261
00262 }
00263
00264 }
00265
00266
00267 for(TrackToVertexAssMap::const_iterator assomap_ite=trackvertexassInput->begin(); assomap_ite!=trackvertexassInput->end(); assomap_ite++){
00268
00269 const VertexRef assomap_vertexref = assomap_ite->key;
00270 const TrackQualityPairVector trckcoll = assomap_ite->val;
00271
00272
00273
00274 if(assomap_vertexref==vertexref_highestpT)
00275 for(unsigned int trckcoll_ite=0; trckcoll_ite<trckcoll.size(); trckcoll_ite++)
00276 trackvertexassOutput->insert(assomap_vertexref,trckcoll[trckcoll_ite]);
00277
00278 }
00279
00280 vertexptsumvector.erase(vertexptsumvector.begin()+highestpT_index);
00281
00282 }
00283
00284 return trackvertexassOutput;
00285
00286 }
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298 std::vector<reco::VertexRef>*
00299 PF_PU_AssoMapAlgos::CreateVertexVector(edm::Handle<reco::VertexCollection> vtxcollH)
00300 {
00301
00302 vector<VertexRef>* output = new vector<VertexRef>();
00303
00304 for(unsigned int index_vtx=0; index_vtx<vtxcollH->size(); ++index_vtx){
00305
00306 VertexRef vertexref(vtxcollH,index_vtx);
00307
00308 output->push_back(vertexref);
00309
00310 }
00311
00312 return output;
00313
00314 }
00315
00316
00317
00318
00319
00320 void
00321 PF_PU_AssoMapAlgos::EraseVertex(std::vector<reco::VertexRef>* vtxcollV, reco::VertexRef toErase)
00322 {
00323
00324 for(unsigned int index_vtx=0; index_vtx<vtxcollV->size(); ++index_vtx){
00325
00326 VertexRef vertexref = vtxcollV->at(index_vtx);
00327
00328 if ( vertexref == toErase ){
00329 vtxcollV->erase(vtxcollV->begin()+index_vtx);
00330 break;
00331 }
00332
00333 }
00334
00335 }
00336
00337
00338
00339
00340
00341
00342 VertexRef
00343 PF_PU_AssoMapAlgos::FindClosestZ(const reco::TrackRef trkref, std::vector<reco::VertexRef>* vtxcollV, double tWeight)
00344 {
00345
00346 double ztrack = trkref->vertex().z();
00347
00348 VertexRef foundVertexRef = vtxcollV->at(0);
00349
00350 double dzmin = 1e5;
00351
00352
00353 for(unsigned int index_vtx=0; index_vtx<vtxcollV->size(); ++index_vtx){
00354
00355 VertexRef vertexref = vtxcollV->at(index_vtx);
00356
00357 double nTracks = sqrt(vertexref->tracksSize());
00358
00359 double z_distance = fabs(ztrack - vertexref->z());
00360
00361 double weightedDistance = z_distance-tWeight*nTracks;
00362
00363 if(weightedDistance<dzmin) {
00364 dzmin = weightedDistance;
00365 foundVertexRef = vertexref;
00366 }
00367
00368 }
00369
00370 return foundVertexRef;
00371 }
00372
00373
00374
00375
00376
00377
00378 VertexRef
00379 PF_PU_AssoMapAlgos::FindClosest3D(TransientTrack transtrk, std::vector<reco::VertexRef>* vtxcollV, double tWeight)
00380 {
00381
00382 VertexRef foundVertexRef = vtxcollV->at(0);
00383
00384 double d3min = 1e5;
00385
00386
00387 for(unsigned int index_vtx=0; index_vtx<vtxcollV->size(); ++index_vtx){
00388
00389 VertexRef vertexref = vtxcollV->at(index_vtx);
00390
00391 double nTracks = sqrt(vertexref->tracksSize());
00392
00393 double distance = 1e5;
00394 pair<bool,Measurement1D> IpPair = IPTools::absoluteImpactParameter3D(transtrk, *vertexref);
00395
00396 if(IpPair.first)
00397 distance = IpPair.second.value();
00398
00399 double weightedDistance = distance-tWeight*nTracks;
00400
00401 if(weightedDistance<d3min) {
00402 d3min = weightedDistance;
00403 foundVertexRef = vertexref;
00404 }
00405
00406 }
00407
00408 return foundVertexRef;
00409 }
00410
00411
00412
00413
00414
00415 double
00416 PF_PU_AssoMapAlgos::dR(const math::XYZPoint& vtx_pos, const math::XYZVector& vtx_mom, edm::Handle<reco::BeamSpot> bsH)
00417 {
00418
00419 double bs_x = bsH->x0();
00420 double bs_y = bsH->y0();
00421 double bs_z = bsH->z0();
00422
00423 double connVec_x = vtx_pos.x() - bs_x;
00424 double connVec_y = vtx_pos.y() - bs_y;
00425 double connVec_z = vtx_pos.z() - bs_z;
00426
00427 double connVec_r = sqrt(connVec_x*connVec_x + connVec_y*connVec_y + connVec_z*connVec_z);
00428 double connVec_theta = acos(connVec_z*1./connVec_r);
00429
00430 double connVec_eta = -1.*log(tan(connVec_theta*1./2.));
00431 double connVec_phi = atan2(connVec_y,connVec_x);
00432
00433 return deltaR(vtx_mom.eta(),vtx_mom.phi(),connVec_eta,connVec_phi);
00434
00435 }
00436
00437
00438
00439
00440
00441 auto_ptr<ConversionCollection>
00442 PF_PU_AssoMapAlgos::GetCleanedConversions(edm::Handle<reco::ConversionCollection> convCollH, Handle<BeamSpot> bsH, bool cleanedColl)
00443 {
00444 auto_ptr<ConversionCollection> cleanedConvColl(new ConversionCollection() );
00445
00446 for (unsigned int convcoll_idx=0; convcoll_idx<convCollH->size(); convcoll_idx++){
00447
00448 ConversionRef convref(convCollH,convcoll_idx);
00449
00450 if(!cleanedColl){
00451 cleanedConvColl->push_back(*convref);
00452 continue;
00453 }
00454
00455 if( (convref->nTracks()==2) &&
00456 (fabs(convref->pairInvariantMass())<=0.1) ){
00457
00458 cleanedConvColl->push_back(*convref);
00459
00460 }
00461
00462 }
00463
00464 return cleanedConvColl;
00465
00466 }
00467
00468
00469
00470
00471
00472
00473 bool
00474 PF_PU_AssoMapAlgos::ComesFromConversion(const TrackRef trackref, const ConversionCollection& cleanedConvColl, Conversion* gamma)
00475 {
00476
00477 for(unsigned int convcoll_ite=0; convcoll_ite<cleanedConvColl.size(); convcoll_ite++){
00478
00479 if(ConversionTools::matchesConversion(trackref,cleanedConvColl.at(convcoll_ite))){
00480
00481 *gamma = cleanedConvColl.at(convcoll_ite);
00482 return true;
00483
00484 }
00485
00486 }
00487
00488 return false;
00489 }
00490
00491
00492
00493
00494
00495
00496 VertexRef
00497 PF_PU_AssoMapAlgos::FindConversionVertex(const reco::TrackRef trackref, const reco::Conversion& gamma, ESHandle<MagneticField> bfH, const EventSetup& iSetup, edm::Handle<reco::BeamSpot> bsH, std::vector<reco::VertexRef>* vtxcollV, double tWeight)
00498 {
00499
00500 math::XYZPoint conv_pos = gamma.conversionVertex().position();
00501
00502 math::XYZVector conv_mom(gamma.refittedPair4Momentum().x(),
00503 gamma.refittedPair4Momentum().y(),
00504 gamma.refittedPair4Momentum().z());
00505
00506 Track photon(trackref->chi2(), trackref->ndof(), conv_pos, conv_mom, 0, trackref->covariance());
00507
00508 TransientTrack transpho(photon, &(*bfH) );
00509 transpho.setBeamSpot(*bsH);
00510 transpho.setES(iSetup);
00511
00512 return FindClosest3D(transpho, vtxcollV, tWeight);
00513
00514 }
00515
00516
00517
00518
00519
00520 auto_ptr<VertexCompositeCandidateCollection>
00521 PF_PU_AssoMapAlgos::GetCleanedKshort(Handle<VertexCompositeCandidateCollection> KshortsH, Handle<BeamSpot> bsH, bool cleanedColl)
00522 {
00523
00524 auto_ptr<VertexCompositeCandidateCollection> cleanedKaonColl(new VertexCompositeCandidateCollection() );
00525
00526 for (unsigned int kscoll_idx=0; kscoll_idx<KshortsH->size(); kscoll_idx++){
00527
00528 VertexCompositeCandidateRef ksref(KshortsH,kscoll_idx);
00529
00530 if(!cleanedColl){
00531 cleanedKaonColl->push_back(*ksref);
00532 continue;
00533 }
00534
00535 VertexDistance3D distanceComputer;
00536
00537 GlobalPoint dec_pos = RecoVertex::convertPos(ksref->vertex());
00538
00539 GlobalError decayVertexError = GlobalError(ksref->vertexCovariance(0,0), ksref->vertexCovariance(0,1), ksref->vertexCovariance(1,1), ksref->vertexCovariance(0,2), ksref->vertexCovariance(1,2), ksref->vertexCovariance(2,2));
00540
00541 math::XYZVector dec_mom(ksref->momentum().x(),
00542 ksref->momentum().y(),
00543 ksref->momentum().z());
00544
00545 GlobalPoint bsPosition = RecoVertex::convertPos(bsH->position());
00546 GlobalError bsError = RecoVertex::convertError(bsH->covariance3D());
00547
00548 double kaon_significance = (distanceComputer.distance(VertexState(bsPosition,bsError), VertexState(dec_pos, decayVertexError))).significance();
00549
00550 if ((ksref->vertex().rho()>=3.) &&
00551 (ksref->vertexNormalizedChi2()<=3.) &&
00552 (fabs(ksref->mass() - kMass)<=0.01) &&
00553 (kaon_significance>15.) &&
00554 (PF_PU_AssoMapAlgos::dR(ksref->vertex(),dec_mom,bsH)<=0.3) ){
00555
00556 cleanedKaonColl->push_back(*ksref);
00557
00558 }
00559
00560 }
00561
00562 return cleanedKaonColl;
00563 }
00564
00565
00566
00567
00568
00569 auto_ptr<VertexCompositeCandidateCollection>
00570 PF_PU_AssoMapAlgos::GetCleanedLambda(Handle<VertexCompositeCandidateCollection> LambdasH, Handle<BeamSpot> bsH, bool cleanedColl)
00571 {
00572
00573 auto_ptr<VertexCompositeCandidateCollection> cleanedLambdaColl(new VertexCompositeCandidateCollection() );
00574
00575 for (unsigned int lambdacoll_idx=0; lambdacoll_idx<LambdasH->size(); lambdacoll_idx++){
00576
00577 VertexCompositeCandidateRef lambdaref(LambdasH,lambdacoll_idx);
00578
00579 if(!cleanedColl){
00580 cleanedLambdaColl->push_back(*lambdaref);
00581 continue;
00582 }
00583
00584 VertexDistance3D distanceComputer;
00585
00586 GlobalPoint dec_pos = RecoVertex::convertPos(lambdaref->vertex());
00587
00588 GlobalError decayVertexError = GlobalError(lambdaref->vertexCovariance(0,0), lambdaref->vertexCovariance(0,1), lambdaref->vertexCovariance(1,1), lambdaref->vertexCovariance(0,2), lambdaref->vertexCovariance(1,2), lambdaref->vertexCovariance(2,2));
00589
00590 math::XYZVector dec_mom(lambdaref->momentum().x(),
00591 lambdaref->momentum().y(),
00592 lambdaref->momentum().z());
00593
00594 GlobalPoint bsPosition = RecoVertex::convertPos(bsH->position());
00595 GlobalError bsError = RecoVertex::convertError(bsH->covariance3D());
00596
00597 double lambda_significance = (distanceComputer.distance(VertexState(bsPosition,bsError), VertexState(dec_pos, decayVertexError))).significance();
00598
00599 if ((lambdaref->vertex().rho()>=3.) &&
00600 (lambdaref->vertexNormalizedChi2()<=3.) &&
00601 (fabs(lambdaref->mass() - lamMass)<=0.005) &&
00602 (lambda_significance>15.) &&
00603 (PF_PU_AssoMapAlgos::dR(lambdaref->vertex(),dec_mom,bsH)<=0.3) ){
00604
00605 cleanedLambdaColl->push_back(*lambdaref);
00606
00607 }
00608
00609 }
00610
00611 return cleanedLambdaColl;
00612 }
00613
00614
00615
00616
00617
00618 bool
00619 PF_PU_AssoMapAlgos::ComesFromV0Decay(const TrackRef trackref, const VertexCompositeCandidateCollection& cleanedKshort, const VertexCompositeCandidateCollection& cleanedLambda, VertexCompositeCandidate* V0)
00620 {
00621
00622
00623 for(VertexCompositeCandidateCollection::const_iterator iKS=cleanedKshort.begin(); iKS!=cleanedKshort.end(); iKS++){
00624
00625 const RecoChargedCandidate *dauCand1 = dynamic_cast<const RecoChargedCandidate*>(iKS->daughter(0));
00626 TrackRef dauTk1 = dauCand1->track();
00627 const RecoChargedCandidate *dauCand2 = dynamic_cast<const RecoChargedCandidate*>(iKS->daughter(1));
00628 TrackRef dauTk2 = dauCand2->track();
00629
00630 if((trackref==dauTk1) || (trackref==dauTk2)){
00631
00632 *V0 = *iKS;
00633 return true;
00634
00635 }
00636
00637 }
00638
00639
00640 for(VertexCompositeCandidateCollection::const_iterator iLambda=cleanedLambda.begin(); iLambda!=cleanedLambda.end(); iLambda++){
00641
00642 const RecoChargedCandidate *dauCand1 = dynamic_cast<const RecoChargedCandidate*>(iLambda->daughter(0));
00643 TrackRef dauTk1 = dauCand1->track();
00644 const RecoChargedCandidate *dauCand2 = dynamic_cast<const RecoChargedCandidate*>(iLambda->daughter(1));
00645 TrackRef dauTk2 = dauCand2->track();
00646
00647 if((trackref==dauTk1) || (trackref==dauTk2)){
00648
00649 *V0 = *iLambda;
00650 return true;
00651
00652 }
00653
00654 }
00655
00656 return false;
00657 }
00658
00659
00660
00661
00662
00663
00664 VertexRef
00665 PF_PU_AssoMapAlgos::FindV0Vertex(const TrackRef trackref, const VertexCompositeCandidate& V0_vtx, ESHandle<MagneticField> bFieldH, const EventSetup& iSetup, Handle<BeamSpot> bsH, std::vector<reco::VertexRef>* vtxcollV, double tWeight)
00666 {
00667
00668 math::XYZPoint dec_pos = V0_vtx.vertex();
00669
00670 math::XYZVector dec_mom(V0_vtx.momentum().x(),
00671 V0_vtx.momentum().y(),
00672 V0_vtx.momentum().z());
00673
00674 Track V0(trackref->chi2(), trackref->ndof(), dec_pos, dec_mom, 0, trackref->covariance());
00675
00676 TransientTrack transV0(V0, &(*bFieldH) );
00677 transV0.setBeamSpot(*bsH);
00678 transV0.setES(iSetup);
00679
00680 return FindClosest3D(transV0, vtxcollV, tWeight);
00681
00682 }
00683
00684
00685
00686
00687
00688
00689 auto_ptr<PFDisplacedVertexCollection>
00690 PF_PU_AssoMapAlgos::GetCleanedNI(Handle<PFDisplacedVertexCollection> NuclIntH, Handle<BeamSpot> bsH, bool cleanedColl)
00691 {
00692
00693 auto_ptr<PFDisplacedVertexCollection> cleanedNIColl(new PFDisplacedVertexCollection() );
00694
00695 for (PFDisplacedVertexCollection::const_iterator niref=NuclIntH->begin(); niref!=NuclIntH->end(); niref++){
00696
00697
00698 if( (niref->isFake()) || !(niref->isNucl()) ) continue;
00699
00700 if(!cleanedColl){
00701 cleanedNIColl->push_back(*niref);
00702 continue;
00703 }
00704
00705 VertexDistance3D distanceComputer;
00706
00707 GlobalPoint ni_pos = RecoVertex::convertPos(niref->position());
00708 GlobalError interactionVertexError = RecoVertex::convertError(niref->error());
00709
00710 math::XYZVector ni_mom(niref->primaryMomentum().x(),
00711 niref->primaryMomentum().y(),
00712 niref->primaryMomentum().z());
00713
00714 GlobalPoint bsPosition = RecoVertex::convertPos(bsH->position());
00715 GlobalError bsError = RecoVertex::convertError(bsH->covariance3D());
00716
00717 double nuclint_significance = (distanceComputer.distance(VertexState(bsPosition,bsError), VertexState(ni_pos, interactionVertexError))).significance();
00718
00719 if ((niref->position().rho()>=3.) &&
00720 (nuclint_significance>15.) &&
00721 (PF_PU_AssoMapAlgos::dR(niref->position(),ni_mom,bsH)<=0.3) ){
00722
00723 cleanedNIColl->push_back(*niref);
00724
00725 }
00726
00727 }
00728
00729 return cleanedNIColl;
00730 }
00731
00732
00733
00734
00735
00736 bool
00737 PF_PU_AssoMapAlgos::ComesFromNI(const TrackRef trackref, const PFDisplacedVertexCollection& cleanedNI, PFDisplacedVertex* displVtx)
00738 {
00739
00740
00741 for(PFDisplacedVertexCollection::const_iterator iDisplV=cleanedNI.begin(); iDisplV!=cleanedNI.end(); iDisplV++){
00742
00743 if(iDisplV->trackWeight(trackref)>1.e-5){
00744
00745 *displVtx = *iDisplV;
00746 return true;
00747
00748 }
00749
00750 }
00751
00752 return false;
00753 }
00754
00755
00756
00757
00758
00759
00760 VertexRef
00761 PF_PU_AssoMapAlgos::FindNIVertex(const TrackRef trackref, const PFDisplacedVertex& displVtx, ESHandle<MagneticField> bFieldH, const EventSetup& iSetup, Handle<BeamSpot> bsH, std::vector<reco::VertexRef>* vtxcollV, double tWeight)
00762 {
00763
00764 TrackCollection refittedTracks = displVtx.refittedTracks();
00765
00766 if((displVtx.isTherePrimaryTracks()) || (displVtx.isThereMergedTracks())){
00767
00768 for(TrackCollection::const_iterator trkcoll_ite=refittedTracks.begin(); trkcoll_ite!=refittedTracks.end(); trkcoll_ite++){
00769
00770 const TrackBaseRef retrackbaseref = displVtx.originalTrack(*trkcoll_ite);
00771
00772 if(displVtx.isIncomingTrack(retrackbaseref)){
00773
00774 VertexRef VOAssociation = TrackWeightAssociation(retrackbaseref, vtxcollV);
00775
00776 if( VOAssociation->trackWeight(retrackbaseref) >= 1.e-5 ){
00777 return VOAssociation;
00778 }
00779
00780 TransientTrack transIncom(*retrackbaseref, &(*bFieldH) );
00781 transIncom.setBeamSpot(*bsH);
00782 transIncom.setES(iSetup);
00783
00784 return FindClosest3D(transIncom, vtxcollV, tWeight);
00785
00786 }
00787
00788 }
00789
00790 }
00791
00792 math::XYZPoint ni_pos = displVtx.position();
00793
00794 math::XYZVector ni_mom(displVtx.primaryMomentum().x(),
00795 displVtx.primaryMomentum().y(),
00796 displVtx.primaryMomentum().z());
00797
00798 Track incom(trackref->chi2(), trackref->ndof(), ni_pos, ni_mom, 0, trackref->covariance());
00799
00800 TransientTrack transIncom(incom, &(*bFieldH) );
00801 transIncom.setBeamSpot(*bsH);
00802 transIncom.setES(iSetup);
00803
00804 return FindClosest3D(transIncom, vtxcollV, tWeight);
00805
00806 }
00807
00808
00809
00810
00811
00812 VertexRef
00813 PF_PU_AssoMapAlgos::TrackWeightAssociation(const TrackBaseRef& trackbaseRef, std::vector<reco::VertexRef>* vtxcollV)
00814 {
00815
00816 VertexRef bestvertexref = vtxcollV->at(0);
00817 float bestweight = 0.;
00818
00819
00820 for(unsigned int index_vtx=0; index_vtx<vtxcollV->size(); ++index_vtx){
00821
00822 VertexRef vertexref = vtxcollV->at(index_vtx);
00823
00824
00825 float weight = vertexref->trackWeight(trackbaseRef);
00826 if(weight>bestweight){
00827 bestweight = weight;
00828 bestvertexref = vertexref;
00829 }
00830
00831 }
00832
00833 return bestvertexref;
00834
00835 }
00836
00837
00838
00839
00840
00841 VertexStepPair
00842 PF_PU_AssoMapAlgos::FindAssociation(const reco::TrackRef& trackref, std::vector<reco::VertexRef>* vtxColl, edm::ESHandle<MagneticField> bfH, const edm::EventSetup& iSetup, edm::Handle<reco::BeamSpot> bsH, int assocNum)
00843 {
00844
00845 const TrackBaseRef& trackbaseRef = TrackBaseRef(trackref);
00846
00847 VertexRef foundVertex;
00848
00849
00850
00851 if ( assocNum>0 ) goto finalStep;
00852
00853
00854
00855 foundVertex = TrackWeightAssociation(trackbaseRef, vtxColl);
00856
00857 if ( foundVertex->trackWeight(trackbaseRef) >= 1.e-5 ){
00858 return make_pair( foundVertex, 0. );
00859 }
00860
00861
00862
00863
00864
00865
00866 if ( input_doReassociation_ ) {
00867
00868
00869
00870 Conversion gamma;
00871 if ( ComesFromConversion(trackref, *cleanedConvCollP, &gamma) ){
00872 foundVertex = FindConversionVertex(trackref, gamma, bfH, iSetup, bsH, vtxColl, input_nTrack_);
00873 return make_pair( foundVertex, 1. );
00874 }
00875
00876
00877
00878 VertexCompositeCandidate V0;
00879 if ( ComesFromV0Decay(trackref, *cleanedKshortCollP, *cleanedLambdaCollP, &V0) ) {
00880 foundVertex = FindV0Vertex(trackref, V0, bfH, iSetup, bsH, vtxColl, input_nTrack_);
00881 return make_pair( foundVertex, 1. );
00882 }
00883
00884 if ( !missingColls ) {
00885
00886
00887
00888 PFDisplacedVertex displVtx;
00889 if ( ComesFromNI(trackref, *cleanedNICollP, &displVtx) ){
00890 foundVertex = FindNIVertex(trackref, displVtx, bfH, iSetup, bsH, vtxColl, input_nTrack_);
00891 return make_pair( foundVertex, 1. );
00892 }
00893
00894 }
00895
00896 }
00897
00898
00899
00900
00901
00902
00903
00904 finalStep:
00905
00906 switch (input_FinalAssociation_) {
00907
00908 case 1:{
00909
00910
00911 foundVertex = FindClosestZ(trackref,vtxColl,input_nTrack_);
00912 break;
00913
00914
00915 }
00916
00917 case 2:{
00918
00919
00920 TransientTrack transtrk(trackref, &(*bfH) );
00921 transtrk.setBeamSpot(*bsH);
00922 transtrk.setES(iSetup);
00923
00924 foundVertex = FindClosest3D(transtrk,vtxColl,input_nTrack_);
00925 break;
00926
00927 }
00928
00929 default:{
00930
00931
00932 foundVertex = vtxColl->at(0);
00933 break;
00934
00935 }
00936
00937 }
00938
00939 return make_pair( foundVertex, 2. );
00940
00941 }
00942
00943
00944
00945
00946
00947 int
00948 PF_PU_AssoMapAlgos::DefineQuality(int assoc_ite, int step, double distance)
00949 {
00950
00951 int quality = 0;
00952
00953 switch (step) {
00954
00955 case 0:{
00956
00957
00958 if ( distance <= tw_90 ) {
00959 quality = 5;
00960 } else {
00961 if ( distance <= tw_70 ) {
00962 quality = 4;
00963 } else {
00964 if ( distance <= tw_50 ) {
00965 quality = 3;
00966 } else {
00967 quality = 2;
00968 }
00969 }
00970 }
00971 break;
00972
00973 }
00974
00975 case 1:{
00976
00977
00978 if ( distance <= sec_70 ) {
00979 quality = 4;
00980 } else {
00981 if ( distance <= sec_50 ) {
00982 quality = 3;
00983 } else {
00984 quality = 2;
00985 }
00986 }
00987 break;
00988
00989 }
00990
00991 case 2:{
00992
00993
00994 if ( assoc_ite == 1 ) {
00995 quality = 1;
00996 } else {
00997 if ( assoc_ite >= 2 ) {
00998 quality = 0;
00999 } else {
01000 if ( distance <= fin_70 ) {
01001 quality = 4;
01002 } else {
01003 if ( distance <= fin_50 ) {
01004 quality = 3;
01005 } else {
01006 quality = 2;
01007 }
01008 }
01009 }
01010 }
01011 break;
01012
01013 }
01014
01015 default:{
01016
01017 quality = -1;
01018 break;
01019 }
01020
01021 }
01022
01023 return quality;
01024
01025 }