00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020 #include "RecoVertex/V0Producer/interface/V0Fitter.h"
00021 #include "CommonTools/CandUtils/interface/AddFourMomenta.h"
00022
00023 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00024 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00025 #include "TrackingTools/PatternTools/interface/ClosestApproachInRPhi.h"
00026 #include "Geometry/CommonDetUnit/interface/GlobalTrackingGeometry.h"
00027 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00028 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00029 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
00030
00031 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00032
00033 #include <Math/Functions.h>
00034 #include <Math/SVector.h>
00035 #include <Math/SMatrix.h>
00036 #include <typeinfo>
00037
00038
00039
00040 const double piMass = 0.13957018;
00041 const double piMassSquared = piMass*piMass;
00042 const double protonMass = 0.938272013;
00043 const double protonMassSquared = protonMass*protonMass;
00044 const double kShortMass = 0.497614;
00045 const double lambdaMass = 1.115683;
00046
00047
00048 V0Fitter::V0Fitter(const edm::ParameterSet& theParameters,
00049 const edm::Event& iEvent, const edm::EventSetup& iSetup) {
00050 using std::string;
00051
00052
00053 recoAlg = theParameters.getParameter<edm::InputTag>("trackRecoAlgorithm");
00054
00055
00056
00057
00058
00059 useRefTrax = theParameters.getParameter<bool>(string("useSmoothing"));
00060
00061
00062 doKshorts = theParameters.getParameter<bool>(string("selectKshorts"));
00063
00064 doLambdas = theParameters.getParameter<bool>(string("selectLambdas"));
00065
00066
00067 chi2Cut = theParameters.getParameter<double>(string("vtxChi2Cut"));
00068 tkChi2Cut = theParameters.getParameter<double>(string("tkChi2Cut"));
00069 tkNhitsCut = theParameters.getParameter<int>(string("tkNhitsCut"));
00070 rVtxCut = theParameters.getParameter<double>(string("rVtxCut"));
00071 vtxSigCut = theParameters.getParameter<double>(string("vtxSignificance2DCut"));
00072 collinCut = theParameters.getParameter<double>(string("collinearityCut"));
00073 kShortMassCut = theParameters.getParameter<double>(string("kShortMassCut"));
00074 lambdaMassCut = theParameters.getParameter<double>(string("lambdaMassCut"));
00075 impactParameterSigCut = theParameters.getParameter<double>(string("impactParameterSigCut"));
00076 mPiPiCut = theParameters.getParameter<double>(string("mPiPiCut"));
00077 tkDCACut = theParameters.getParameter<double>(string("tkDCACut"));
00078 vtxFitter = theParameters.getParameter<edm::InputTag>("vertexFitter");
00079 innerHitPosCut = theParameters.getParameter<double>(string("innerHitPosCut"));
00080 std::vector<std::string> qual = theParameters.getParameter<std::vector<std::string> >("trackQualities");
00081 for (unsigned int ndx = 0; ndx < qual.size(); ndx++) {
00082 qualities.push_back(reco::TrackBase::qualityByName(qual[ndx]));
00083 }
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093 fitAll(iEvent, iSetup);
00094
00095
00096
00097
00098
00099 }
00100
00101 V0Fitter::~V0Fitter() {
00102 }
00103
00104
00105 void V0Fitter::fitAll(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
00106
00107 using std::vector;
00108 using std::cout;
00109 using std::endl;
00110 using namespace reco;
00111 using namespace edm;
00112
00113
00114
00115 std::vector<TrackRef> theTrackRefs;
00116 std::vector<TransientTrack> theTransTracks;
00117
00118
00119 Handle<reco::TrackCollection> theTrackHandle;
00120 Handle<reco::BeamSpot> theBeamSpotHandle;
00121 ESHandle<MagneticField> bFieldHandle;
00122 ESHandle<TrackerGeometry> trackerGeomHandle;
00123 ESHandle<GlobalTrackingGeometry> globTkGeomHandle;
00124
00125
00126
00127
00128 iEvent.getByLabel(recoAlg, theTrackHandle);
00129 iEvent.getByLabel(std::string("offlineBeamSpot"), theBeamSpotHandle);
00130 if( !theTrackHandle->size() ) return;
00131 iSetup.get<IdealMagneticFieldRecord>().get(bFieldHandle);
00132 iSetup.get<TrackerDigiGeometryRecord>().get(trackerGeomHandle);
00133 iSetup.get<GlobalTrackingGeometryRecord>().get(globTkGeomHandle);
00134
00135 trackerGeom = trackerGeomHandle.product();
00136 magField = bFieldHandle.product();
00137
00138
00139 for(unsigned int indx = 0; indx < theTrackHandle->size(); indx++) {
00140 TrackRef tmpRef( theTrackHandle, indx );
00141 bool quality_ok = true;
00142 if (qualities.size()!=0) {
00143 quality_ok = false;
00144 for (unsigned int ndx_ = 0; ndx_ < qualities.size(); ndx_++) {
00145 if (tmpRef->quality(qualities[ndx_])){
00146 quality_ok = true;
00147 break;
00148 }
00149 }
00150 }
00151 if( !quality_ok ) continue;
00152
00153
00154 if( tmpRef->normalizedChi2() < tkChi2Cut &&
00155 tmpRef->numberOfValidHits() >= tkNhitsCut ) {
00156 TransientTrack tmpTk( *tmpRef, &(*bFieldHandle), globTkGeomHandle );
00157 TrajectoryStateTransform theTransform;
00158 FreeTrajectoryState initialFTS = theTransform.initialFreeState(*tmpRef, magField);
00159 TSCBLBuilderNoMaterial blsBuilder;
00160 TrajectoryStateClosestToBeamLine tscb( blsBuilder(initialFTS, *theBeamSpotHandle) );
00161
00162 if( tscb.isValid() ) {
00163 if( tscb.transverseImpactParameter().significance() > impactParameterSigCut ) {
00164 theTrackRefs.push_back( tmpRef );
00165 theTransTracks.push_back( tmpTk );
00166 }
00167 }
00168 }
00169 }
00170
00171
00172
00173
00174
00175 for(unsigned int trdx1 = 0; trdx1 < theTrackRefs.size(); trdx1++) {
00176
00177 for(unsigned int trdx2 = trdx1 + 1; trdx2 < theTrackRefs.size(); trdx2++) {
00178
00179
00180 std::vector<TransientTrack> transTracks;
00181
00182 TrackRef positiveTrackRef;
00183 TrackRef negativeTrackRef;
00184 TransientTrack* posTransTkPtr = 0;
00185 TransientTrack* negTransTkPtr = 0;
00186
00187
00188
00189
00190 if(theTrackRefs[trdx1]->charge() < 0. &&
00191 theTrackRefs[trdx2]->charge() > 0.) {
00192 negativeTrackRef = theTrackRefs[trdx1];
00193 positiveTrackRef = theTrackRefs[trdx2];
00194 negTransTkPtr = &theTransTracks[trdx1];
00195 posTransTkPtr = &theTransTracks[trdx2];
00196 }
00197 else if(theTrackRefs[trdx1]->charge() > 0. &&
00198 theTrackRefs[trdx2]->charge() < 0.) {
00199 negativeTrackRef = theTrackRefs[trdx2];
00200 positiveTrackRef = theTrackRefs[trdx1];
00201 negTransTkPtr = &theTransTracks[trdx2];
00202 posTransTkPtr = &theTransTracks[trdx1];
00203 }
00204
00205
00206 else continue;
00207
00208
00209 transTracks.push_back(*posTransTkPtr);
00210 transTracks.push_back(*negTransTkPtr);
00211
00212
00213 FreeTrajectoryState posState = posTransTkPtr->impactPointTSCP().theState();
00214 FreeTrajectoryState negState = negTransTkPtr->impactPointTSCP().theState();
00215
00216 if( !posTransTkPtr->impactPointTSCP().isValid() || !negTransTkPtr->impactPointTSCP().isValid() ) continue;
00217
00218
00219 ClosestApproachInRPhi cApp;
00220 cApp.calculate(posState, negState);
00221 if( !cApp.status() ) continue;
00222 float dca = fabs( cApp.distance() );
00223 GlobalPoint cxPt = cApp.crossingPoint();
00224
00225 if (dca < 0. || dca > tkDCACut) continue;
00226 if (sqrt( cxPt.x()*cxPt.x() + cxPt.y()*cxPt.y() ) > 120.
00227 || std::abs(cxPt.z()) > 300.) continue;
00228
00229
00230 TrajectoryStateClosestToPoint posTSCP =
00231 posTransTkPtr->trajectoryStateClosestToPoint( cxPt );
00232 TrajectoryStateClosestToPoint negTSCP =
00233 negTransTkPtr->trajectoryStateClosestToPoint( cxPt );
00234
00235 if( !posTSCP.isValid() || !negTSCP.isValid() ) continue;
00236
00237
00238
00239
00240
00241
00242
00243 double totalE = sqrt( posTSCP.momentum().mag2() + piMassSquared ) +
00244 sqrt( negTSCP.momentum().mag2() + piMassSquared );
00245 double totalESq = totalE*totalE;
00246 double totalPSq =
00247 ( posTSCP.momentum() + negTSCP.momentum() ).mag2();
00248 double mass = sqrt( totalESq - totalPSq);
00249
00250
00251
00252 if( mass > mPiPiCut ) continue;
00253
00254
00255 TransientVertex theRecoVertex;
00256 if(vtxFitter == std::string("KalmanVertexFitter")) {
00257 KalmanVertexFitter theKalmanFitter(useRefTrax == 0 ? false : true);
00258 theRecoVertex = theKalmanFitter.vertex(transTracks);
00259 }
00260 else if (vtxFitter == std::string("AdaptiveVertexFitter")) {
00261 useRefTrax = false;
00262 AdaptiveVertexFitter theAdaptiveFitter;
00263 theRecoVertex = theAdaptiveFitter.vertex(transTracks);
00264 }
00265
00266
00267
00268 if( !theRecoVertex.isValid() || theRecoVertex.totalChiSquared() < 0. ) {
00269 continue;
00270 }
00271
00272
00273 reco::Vertex theVtx = theRecoVertex;
00274
00275
00276 std::vector<TransientTrack> refittedTrax;
00277 if( theRecoVertex.hasRefittedTracks() ) {
00278 refittedTrax = theRecoVertex.refittedTracks();
00279 }
00280
00281
00282
00283
00284
00285 typedef ROOT::Math::SMatrix<double, 3, 3, ROOT::Math::MatRepSym<double, 3> > SMatrixSym3D;
00286 typedef ROOT::Math::SVector<double, 3> SVector3;
00287
00288 GlobalPoint vtxPos(theVtx.x(), theVtx.y(), theVtx.z());
00289
00290 GlobalPoint beamSpotPos(theBeamSpotHandle->position().x(),
00291 theBeamSpotHandle->position().y(),
00292 theBeamSpotHandle->position().z());
00293
00294 SMatrixSym3D totalCov = theBeamSpotHandle->covariance3D() + theVtx.covariance();
00295 SVector3 distanceVector(vtxPos.x() - beamSpotPos.x(),
00296 vtxPos.y() - beamSpotPos.y(),
00297 0.);
00298
00299
00300 double rVtxMag = ROOT::Math::Mag(distanceVector);
00301 double sigmaRvtxMag = sqrt(ROOT::Math::Similarity(totalCov, distanceVector)) / rVtxMag;
00302
00303
00304
00305
00306 if( innerHitPosCut > 0. && positiveTrackRef->innerOk() ) {
00307 reco::Vertex::Point posTkHitPos = positiveTrackRef->innerPosition();
00308 double posTkHitPosD2 =
00309 (posTkHitPos.x()-beamSpotPos.x())*(posTkHitPos.x()-beamSpotPos.x()) +
00310 (posTkHitPos.y()-beamSpotPos.y())*(posTkHitPos.y()-beamSpotPos.y());
00311 if( sqrt( posTkHitPosD2 ) < ( rVtxMag - sigmaRvtxMag*innerHitPosCut )
00312 ) {
00313 continue;
00314 }
00315 }
00316 if( innerHitPosCut > 0. && negativeTrackRef->innerOk() ) {
00317 reco::Vertex::Point negTkHitPos = negativeTrackRef->innerPosition();
00318 double negTkHitPosD2 =
00319 (negTkHitPos.x()-beamSpotPos.x())*(negTkHitPos.x()-beamSpotPos.x()) +
00320 (negTkHitPos.y()-beamSpotPos.y())*(negTkHitPos.y()-beamSpotPos.y());
00321 if( sqrt( negTkHitPosD2 ) < ( rVtxMag - sigmaRvtxMag*innerHitPosCut )
00322 ) {
00323 continue;
00324 }
00325 }
00326
00327 if( theVtx.normalizedChi2() > chi2Cut ||
00328 rVtxMag < rVtxCut ||
00329 rVtxMag / sigmaRvtxMag < vtxSigCut ) {
00330 continue;
00331 }
00332
00333
00334
00335 TrajectoryStateClosestToPoint* trajPlus;
00336 TrajectoryStateClosestToPoint* trajMins;
00337
00338 if( useRefTrax && refittedTrax.size() > 1 ) {
00339
00340 std::vector<TransientTrack>::iterator traxIter = refittedTrax.begin(),
00341 traxEnd = refittedTrax.end();
00342
00343
00344
00345 TransientTrack* thePositiveRefTrack = 0;
00346 TransientTrack* theNegativeRefTrack = 0;
00347
00348 for( ; traxIter != traxEnd; ++traxIter) {
00349 if( traxIter->track().charge() > 0. ) {
00350 thePositiveRefTrack = new TransientTrack(*traxIter);
00351 }
00352 else if (traxIter->track().charge() < 0.) {
00353 theNegativeRefTrack = new TransientTrack(*traxIter);
00354 }
00355 }
00356 trajPlus = new TrajectoryStateClosestToPoint(
00357 thePositiveRefTrack->trajectoryStateClosestToPoint(vtxPos));
00358 trajMins = new TrajectoryStateClosestToPoint(
00359 theNegativeRefTrack->trajectoryStateClosestToPoint(vtxPos));
00360 delete thePositiveRefTrack;
00361 delete theNegativeRefTrack;
00362 }
00363 else {
00364 trajPlus = new TrajectoryStateClosestToPoint(
00365 posTransTkPtr->trajectoryStateClosestToPoint(vtxPos));
00366 trajMins = new TrajectoryStateClosestToPoint(
00367 negTransTkPtr->trajectoryStateClosestToPoint(vtxPos));
00368
00369 }
00370
00371 if( !trajPlus->isValid() || !trajMins->isValid() ) continue;
00372
00373 posTransTkPtr = negTransTkPtr = 0;
00374
00375 GlobalVector positiveP(trajPlus->momentum());
00376 GlobalVector negativeP(trajMins->momentum());
00377 GlobalVector totalP(positiveP + negativeP);
00378
00379
00380 delete trajPlus;
00381 delete trajMins;
00382 trajPlus = trajMins = 0;
00383
00384
00385
00386 double piPlusE = sqrt( positiveP.mag2() + piMassSquared );
00387 double piMinusE = sqrt( negativeP.mag2() + piMassSquared );
00388 double protonE = sqrt( positiveP.mag2() + protonMassSquared );
00389 double antiProtonE = sqrt( negativeP.mag2() + protonMassSquared );
00390 double kShortETot = piPlusE + piMinusE;
00391 double lambdaEtot = protonE + piMinusE;
00392 double lambdaBarEtot = antiProtonE + piPlusE;
00393
00394 using namespace reco;
00395
00396
00397 const Particle::LorentzVector kShortP4(totalP.x(),
00398 totalP.y(), totalP.z(),
00399 kShortETot);
00400 const Particle::LorentzVector lambdaP4(totalP.x(),
00401 totalP.y(), totalP.z(),
00402 lambdaEtot);
00403 const Particle::LorentzVector lambdaBarP4(totalP.x(),
00404 totalP.y(), totalP.z(),
00405 lambdaBarEtot);
00406
00407 Particle::Point vtx(theVtx.x(), theVtx.y(), theVtx.z());
00408 const Vertex::CovarianceMatrix vtxCov(theVtx.covariance());
00409 double vtxChi2(theVtx.chi2());
00410 double vtxNdof(theVtx.ndof());
00411
00412
00413 VertexCompositeCandidate* theKshort = 0;
00414 VertexCompositeCandidate* theLambda = 0;
00415 VertexCompositeCandidate* theLambdaBar = 0;
00416
00417 if( doKshorts ) {
00418 theKshort = new VertexCompositeCandidate(0, kShortP4, vtx, vtxCov, vtxChi2, vtxNdof);
00419 }
00420 if( doLambdas ) {
00421 if( positiveP.mag() > negativeP.mag() ) {
00422 theLambda =
00423 new VertexCompositeCandidate(0, lambdaP4, vtx, vtxCov, vtxChi2, vtxNdof);
00424 }
00425 else {
00426 theLambdaBar =
00427 new VertexCompositeCandidate(0, lambdaBarP4, vtx, vtxCov, vtxChi2, vtxNdof);
00428 }
00429 }
00430
00431
00432 RecoChargedCandidate
00433 thePiPlusCand(1, Particle::LorentzVector(positiveP.x(),
00434 positiveP.y(), positiveP.z(),
00435 piPlusE), vtx);
00436 thePiPlusCand.setTrack(positiveTrackRef);
00437
00438 RecoChargedCandidate
00439 thePiMinusCand(-1, Particle::LorentzVector(negativeP.x(),
00440 negativeP.y(), negativeP.z(),
00441 piMinusE), vtx);
00442 thePiMinusCand.setTrack(negativeTrackRef);
00443
00444
00445 RecoChargedCandidate
00446 theProtonCand(1, Particle::LorentzVector(positiveP.x(),
00447 positiveP.y(), positiveP.z(),
00448 protonE), vtx);
00449 theProtonCand.setTrack(positiveTrackRef);
00450
00451 RecoChargedCandidate
00452 theAntiProtonCand(-1, Particle::LorentzVector(negativeP.x(),
00453 negativeP.y(), negativeP.z(),
00454 antiProtonE), vtx);
00455 theAntiProtonCand.setTrack(negativeTrackRef);
00456
00457
00458 AddFourMomenta addp4;
00459
00460
00461 if( doKshorts ) {
00462 theKshort->addDaughter(thePiPlusCand);
00463 theKshort->addDaughter(thePiMinusCand);
00464 theKshort->setPdgId(310);
00465 addp4.set( *theKshort );
00466 if( theKshort->mass() < kShortMass + kShortMassCut &&
00467 theKshort->mass() > kShortMass - kShortMassCut ) {
00468 theKshorts.push_back( *theKshort );
00469 }
00470 }
00471
00472 if( doLambdas && theLambda ) {
00473 theLambda->addDaughter(theProtonCand);
00474 theLambda->addDaughter(thePiMinusCand);
00475 theLambda->setPdgId(3122);
00476 addp4.set( *theLambda );
00477 if( theLambda->mass() < lambdaMass + lambdaMassCut &&
00478 theLambda->mass() > lambdaMass - lambdaMassCut ) {
00479 theLambdas.push_back( *theLambda );
00480 }
00481 }
00482 else if ( doLambdas && theLambdaBar ) {
00483 theLambdaBar->addDaughter(theAntiProtonCand);
00484 theLambdaBar->addDaughter(thePiPlusCand);
00485 theLambdaBar->setPdgId(-3122);
00486 addp4.set( *theLambdaBar );
00487 if( theLambdaBar->mass() < lambdaMass + lambdaMassCut &&
00488 theLambdaBar->mass() > lambdaMass - lambdaMassCut ) {
00489 theLambdas.push_back( *theLambdaBar );
00490 }
00491 }
00492
00493 if(theKshort) delete theKshort;
00494 if(theLambda) delete theLambda;
00495 if(theLambdaBar) delete theLambdaBar;
00496 theKshort = theLambda = theLambdaBar = 0;
00497
00498 }
00499 }
00500 }
00501
00502
00503 const reco::VertexCompositeCandidateCollection& V0Fitter::getKshorts() const {
00504 return theKshorts;
00505 }
00506
00507 const reco::VertexCompositeCandidateCollection& V0Fitter::getLambdas() const {
00508 return theLambdas;
00509 }
00510
00511
00512
00513 double V0Fitter::findV0MassError(const GlobalPoint &vtxPos, std::vector<reco::TransientTrack> dauTracks) {
00514 return -1.;
00515 }
00516
00517
00518
00519
00520
00521
00522
00523
00524
00525
00526
00527
00528
00529
00530
00531
00532
00533
00534
00535
00536
00537
00538
00539
00540
00541
00542
00543
00544
00545
00546
00547
00548
00549
00550
00551
00552
00553
00554
00555