CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoEgamma/EgammaPhotonProducers/src/ConversionProducer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    ConversionProducer
00004 // Class:      ConversionProducer
00005 //
00013 //
00014 // Original Authors:  Hongliang Liu
00015 //         Created:  Thu Mar 13 17:40:48 CDT 2008
00016 // $Id: ConversionProducer.cc,v 1.14 2011/12/23 08:13:35 innocent Exp $
00017 //
00018 //
00019 
00020 
00021 // system include files
00022 #include <memory>
00023 #include <map>
00024 
00025 
00026 #include "RecoEgamma/EgammaPhotonProducers/interface/ConversionProducer.h"
00027 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00028 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00029 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00030 #include "DataFormats/TrajectoryState/interface/PTrajectoryStateOnDet.h"
00031 
00032 #include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00033 #include "TrackingTools/MaterialEffects/interface/PropagatorWithMaterial.h"
00034 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00035 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00036 
00037 #include "TrackingTools/TransientTrack/interface/TrackTransientTrack.h"
00038 #include "DataFormats/EgammaTrackReco/interface/ConversionTrack.h"
00039 
00040 
00041 #include "DataFormats/VertexReco/interface/Vertex.h"
00042 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00043 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00044 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
00045 
00046 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00047 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00048 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00049 
00050 #include "TrackingTools/PatternTools/interface/ClosestApproachInRPhi.h"
00051 #include "TrackingTools/PatternTools/interface/TwoTrackMinimumDistance.h"
00052 
00053 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00054 #include "DataFormats/Math/interface/deltaPhi.h"
00055 #include "RecoEgamma/EgammaPhotonAlgos/interface/ConversionVertexFinder.h"
00056 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00057 #include "RecoEgamma/EgammaPhotonAlgos/interface/TangentApproachInRPhi.h"
00058 #include "RecoEgamma/EgammaPhotonAlgos/interface/ConversionHitChecker.h"
00059 
00060 
00061 //Kinematic constraint vertex fitter
00062 #include "RecoVertex/KinematicFitPrimitives/interface/ParticleMass.h"
00063 #include "RecoVertex/KinematicFitPrimitives/interface/MultiTrackKinematicConstraint.h"
00064 #include <RecoVertex/KinematicFitPrimitives/interface/KinematicParticleFactoryFromTransientTrack.h>
00065 #include "RecoVertex/KinematicFit/interface/KinematicConstrainedVertexFitter.h"
00066 #include "RecoVertex/KinematicFit/interface/TwoTrackMassKinematicConstraint.h"
00067 #include "RecoVertex/KinematicFit/interface/KinematicParticleVertexFitter.h"
00068 #include "RecoVertex/KinematicFit/interface/KinematicParticleFitter.h"
00069 #include "RecoVertex/KinematicFit/interface/MassKinematicConstraint.h"
00070 #include "RecoVertex/KinematicFit/interface/ColinearityKinematicConstraint.h"
00071 
00072 
00073 
00074 ConversionProducer::ConversionProducer(const edm::ParameterSet& iConfig):
00075   theVertexFinder_(0)
00076 
00077 {
00078   algoName_ = iConfig.getParameter<std::string>( "AlgorithmName" );
00079 
00080   src_ = iConfig.getParameter<edm::InputTag>("src");
00081 
00082   maxNumOfTrackInPU_ = iConfig.getParameter<int>("maxNumOfTrackInPU");
00083   maxTrackRho_ = iConfig.getParameter<double>("maxTrackRho");
00084   maxTrackZ_ = iConfig.getParameter<double>("maxTrackZ");
00085 
00086   allowTrackBC_ = iConfig.getParameter<bool>("AllowTrackBC");
00087   allowD0_ = iConfig.getParameter<bool>("AllowD0");
00088   allowDeltaPhi_ = iConfig.getParameter<bool>("AllowDeltaPhi");
00089   allowDeltaCot_ = iConfig.getParameter<bool>("AllowDeltaCot");
00090   allowMinApproach_ = iConfig.getParameter<bool>("AllowMinApproach");
00091   allowOppCharge_ = iConfig.getParameter<bool>("AllowOppCharge");
00092 
00093   allowVertex_ = iConfig.getParameter<bool>("AllowVertex");
00094 
00095   bypassPreselGsf_ = iConfig.getParameter<bool>("bypassPreselGsf");
00096   bypassPreselEcal_ = iConfig.getParameter<bool>("bypassPreselEcal");
00097   bypassPreselEcalEcal_ = iConfig.getParameter<bool>("bypassPreselEcalEcal");
00098   
00099   deltaEta_ = iConfig.getParameter<double>("deltaEta");
00100   
00101   halfWayEta_ = iConfig.getParameter<double>("HalfwayEta");//open angle to search track matches with BC
00102 
00103   d0Cut_ = iConfig.getParameter<double>("d0");
00104     
00105   usePvtx_ = iConfig.getParameter<bool>("UsePvtx");//if use primary vertices
00106 
00107   vertexProducer_   = iConfig.getParameter<std::string>("primaryVertexProducer");
00108   
00109 
00110   //Track-cluster matching eta and phi cuts
00111   dEtaTkBC_ = iConfig.getParameter<double>("dEtaTrackBC");//TODO research on cut endcap/barrel
00112   dPhiTkBC_ = iConfig.getParameter<double>("dPhiTrackBC");
00113   
00114   bcBarrelCollection_     = iConfig.getParameter<edm::InputTag>("bcBarrelCollection");
00115   bcEndcapCollection_     = iConfig.getParameter<edm::InputTag>("bcEndcapCollection");
00116   
00117   scBarrelProducer_       = iConfig.getParameter<edm::InputTag>("scBarrelProducer");
00118   scEndcapProducer_       = iConfig.getParameter<edm::InputTag>("scEndcapProducer");
00119   
00120   energyBC_               = iConfig.getParameter<double>("EnergyBC");//BC energy threshold
00121   energyTotalBC_          = iConfig.getParameter<double>("EnergyTotalBC");//BC energy threshold
00122   minSCEt_                = iConfig.getParameter<double>("minSCEt");//super cluster energy threshold
00123   dEtacutForSCmatching_     = iConfig.getParameter<double>("dEtacutForSCmatching");// dEta between conversion momentum direction and SC position
00124   dPhicutForSCmatching_     = iConfig.getParameter<double>("dPhicutForSCmatching");// dPhi between conversion momentum direction and SC position
00125 
00126   
00127    
00128 
00129   //Track cuts on left right track: at least one leg reaches ECAL
00130   //Left track: must exist, must reach Ecal and match BC, so loose cut on Chi2 and tight on hits
00131   //Right track: not necessary to exist (if allowSingleLeg_), not necessary to reach ECAL or match BC, so tight cut on Chi2 and loose on hits
00132   maxChi2Left_ =  iConfig.getParameter<double>("MaxChi2Left");
00133   maxChi2Right_ =  iConfig.getParameter<double>("MaxChi2Right");
00134   minHitsLeft_ = iConfig.getParameter<int>("MinHitsLeft");
00135   minHitsRight_ = iConfig.getParameter<int>("MinHitsRight");
00136 
00137   //Track Open angle cut on delta cot(theta) and delta phi
00138   deltaCotTheta_ = iConfig.getParameter<double>("DeltaCotTheta");
00139   deltaPhi_ = iConfig.getParameter<double>("DeltaPhi");
00140   minApproachLow_ = iConfig.getParameter<double>("MinApproachLow");
00141   minApproachHigh_ = iConfig.getParameter<double>("MinApproachHigh");
00142   
00143 
00144   // if allow single track collection, by default False
00145   allowSingleLeg_ = iConfig.getParameter<bool>("AllowSingleLeg");
00146   rightBC_ = iConfig.getParameter<bool>("AllowRightBC");
00147 
00148   //track inner position dz cut, need RECO
00149   dzCut_ = iConfig.getParameter<double>("dz");
00150   //track analytical cross cut
00151   r_cut = iConfig.getParameter<double>("rCut");
00152   vtxChi2_ = iConfig.getParameter<double>("vtxChi2");
00153 
00154 
00155   theVertexFinder_ = new ConversionVertexFinder ( iConfig );
00156 
00157   thettbuilder_ = 0;
00158 
00159   //output
00160   ConvertedPhotonCollection_     = iConfig.getParameter<std::string>("convertedPhotonCollection");
00161 
00162   produces< reco::ConversionCollection >(ConvertedPhotonCollection_);
00163 
00164 }
00165 
00166 
00167 ConversionProducer::~ConversionProducer()
00168 {
00169 
00170   // do anything here that needs to be done at desctruction time
00171   // (e.g. close files, deallocate resources etc.)
00172   delete theVertexFinder_;
00173 }
00174 
00175 
00176 // ------------ method called to produce the data  ------------
00177 void
00178 ConversionProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00179 {
00180   using namespace edm;
00181 
00182   reco::ConversionCollection outputConvPhotonCollection;
00183   std::auto_ptr<reco::ConversionCollection> outputConvPhotonCollection_p(new reco::ConversionCollection);
00184 
00185   //std::cout << " ConversionProducer::produce " << std::endl;
00186   //Read multiple track input collections
00187 
00188   edm::Handle<edm::View<reco::ConversionTrack> > trackCollectionHandle;
00189   iEvent.getByLabel(src_,trackCollectionHandle);    
00190 
00191   //build map of ConversionTracks ordered in eta
00192   std::multimap<float, edm::Ptr<reco::ConversionTrack> > convTrackMap;
00193   for (edm::PtrVector<reco::ConversionTrack>::const_iterator tk_ref = trackCollectionHandle->ptrVector().begin(); tk_ref != trackCollectionHandle->ptrVector().end(); ++tk_ref ){
00194     convTrackMap.insert(std::make_pair((*tk_ref)->track()->eta(),*tk_ref));
00195   }
00196 
00197   edm::Handle<reco::VertexCollection> vertexHandle;
00198   reco::VertexCollection vertexCollection;
00199   if (usePvtx_){
00200     iEvent.getByLabel(vertexProducer_, vertexHandle);
00201     if (!vertexHandle.isValid()) {
00202       edm::LogError("ConversionProducer") << "Error! Can't get the product primary Vertex Collection "<< "\n";
00203       usePvtx_ = false;
00204     }
00205     if (usePvtx_)
00206       vertexCollection = *(vertexHandle.product());
00207   }
00208     
00209   edm::ESHandle<TransientTrackBuilder> hTransientTrackBuilder;
00210   iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",hTransientTrackBuilder);
00211   thettbuilder_ = hTransientTrackBuilder.product();
00212     
00213   reco::Vertex the_pvtx;
00214   //because the priamry vertex is sorted by quality, the first one is the best
00215   if (!vertexCollection.empty())
00216     the_pvtx = *(vertexCollection.begin());
00217     
00218   if (trackCollectionHandle->size()> maxNumOfTrackInPU_){
00219     iEvent.put( outputConvPhotonCollection_p, ConvertedPhotonCollection_);
00220     return;
00221   }
00222     
00223     
00224   // build Super and Basic cluster geometry map to search in eta bounds for clusters
00225   std::multimap<double, reco::CaloClusterPtr> basicClusterPtrs;
00226   std::multimap<double, reco::CaloClusterPtr> superClusterPtrs;
00227 
00228   
00229   buildSuperAndBasicClusterGeoMap(iEvent,basicClusterPtrs,superClusterPtrs);
00230       
00231   buildCollection( iEvent, iSetup, convTrackMap,  superClusterPtrs, basicClusterPtrs, the_pvtx, outputConvPhotonCollection);//allow empty basicClusterPtrs
00232     
00233   outputConvPhotonCollection_p->assign(outputConvPhotonCollection.begin(), outputConvPhotonCollection.end());
00234   iEvent.put( outputConvPhotonCollection_p, ConvertedPhotonCollection_);
00235     
00236 }
00237 
00238 
00239 void ConversionProducer::buildSuperAndBasicClusterGeoMap(const edm::Event& iEvent,  
00240                                                          std::multimap<double, reco::CaloClusterPtr>& basicClusterPtrs,
00241                                                          std::multimap<double, reco::CaloClusterPtr>& superClusterPtrs){
00242 
00243   // Get the Super Cluster collection in the Barrel
00244   edm::Handle<edm::View<reco::CaloCluster> > scBarrelHandle;
00245   iEvent.getByLabel(scBarrelProducer_,scBarrelHandle);
00246   if (!scBarrelHandle.isValid()) {
00247     edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the product "<<scBarrelProducer_;
00248   }
00249     
00250   // Get the Super Cluster collection in the Endcap
00251   edm::Handle<edm::View<reco::CaloCluster> > scEndcapHandle;
00252   iEvent.getByLabel(scEndcapProducer_,scEndcapHandle);
00253   if (!scEndcapHandle.isValid()) {
00254     edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the product "<<scEndcapProducer_;
00255   }
00256     
00257     
00258   edm::Handle<edm::View<reco::CaloCluster> > bcBarrelHandle;
00259   edm::Handle<edm::View<reco::CaloCluster> > bcEndcapHandle;//TODO check cluster type if BasicCluster or PFCluster
00260 
00261   iEvent.getByLabel( bcBarrelCollection_, bcBarrelHandle);
00262   if (!bcBarrelHandle.isValid()) {
00263     edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the product "<<bcBarrelCollection_;
00264   }
00265 
00266   iEvent.getByLabel( bcEndcapCollection_, bcEndcapHandle);
00267   if (! bcEndcapHandle.isValid()) {
00268     edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the product "<<bcEndcapCollection_;
00269   }
00270 
00271   edm::Handle<edm::View<reco::CaloCluster> > bcHandle = bcBarrelHandle;
00272   edm::Handle<edm::View<reco::CaloCluster> > scHandle = scBarrelHandle;
00273 
00274   if ( bcHandle.isValid()  ) {    
00275     for (unsigned jj = 0; jj < 2; ++jj ){
00276       for (unsigned ii = 0; ii < bcHandle->size(); ++ii ) {
00277         if (bcHandle->ptrAt(ii)->energy()>energyBC_)
00278           basicClusterPtrs.insert(std::make_pair(bcHandle->ptrAt(ii)->position().eta(), bcHandle->ptrAt(ii)));
00279       }
00280       bcHandle = bcEndcapHandle;
00281     }
00282   }
00283 
00284 
00285   if ( scHandle.isValid()  ) {
00286     for (unsigned jj = 0; jj < 2; ++jj ){
00287       for (unsigned ii = 0; ii < scHandle->size(); ++ii ) {
00288         if (scHandle->ptrAt(ii)->energy()>minSCEt_)
00289           superClusterPtrs.insert(std::make_pair(scHandle->ptrAt(ii)->position().eta(), scHandle->ptrAt(ii)));
00290       }
00291       scHandle = scEndcapHandle;
00292     } 
00293   }
00294 
00295 
00296 }
00297 
00298 
00299 void ConversionProducer::buildCollection(edm::Event& iEvent, const edm::EventSetup& iSetup,
00300                                          const std::multimap<float, edm::Ptr<reco::ConversionTrack> >& allTracks,
00301                                          const std::multimap<double, reco::CaloClusterPtr>& superClusterPtrs,
00302                                          const std::multimap<double, reco::CaloClusterPtr>& basicClusterPtrs,
00303                                          const reco::Vertex& the_pvtx,
00304                                          reco::ConversionCollection & outputConvPhotonCollection){
00305   
00306   
00307   edm::ESHandle<TrackerGeometry> trackerGeomHandle;
00308   edm::ESHandle<MagneticField> magFieldHandle;
00309   iSetup.get<TrackerDigiGeometryRecord>().get( trackerGeomHandle );
00310   iSetup.get<IdealMagneticFieldRecord>().get( magFieldHandle );
00311   
00312   const TrackerGeometry* trackerGeom = trackerGeomHandle.product();
00313   const MagneticField* magField = magFieldHandle.product();
00314   
00315 //   std::vector<math::XYZPointF> trackImpactPosition;
00316 //   trackImpactPosition.reserve(allTracks.size());//track impact position at ECAL
00317 //   std::vector<bool> trackValidECAL;//Does this track reach ECAL basic cluster (reach ECAL && match with BC)
00318 //   trackValidECAL.assign(allTracks.size(), false);
00319 //   
00320 //   std::vector<reco::CaloClusterPtr> trackMatchedBC;
00321 //   reco::CaloClusterPtr empty_bc;
00322 //   trackMatchedBC.assign(allTracks.size(), empty_bc);//TODO find a better way to avoid copy constructor
00323 //   
00324 //   std::vector<int> bcHandleId;//the associated BC handle id, -1 invalid, 0 barrel 1 endcap
00325 //   bcHandleId.assign(allTracks.size(), -1);
00326   
00327   // not used    std::multimap<double, int> trackInnerEta;//Track innermost state Eta map to TrackRef index, to be used in track pair sorting
00328   
00329   std::map<edm::Ptr<reco::ConversionTrack>, math::XYZPointF> trackImpactPosition;
00330   std::map<edm::Ptr<reco::ConversionTrack>, reco::CaloClusterPtr> trackMatchedBC;
00331   
00332   ConversionHitChecker hitChecker;
00333   
00334   
00335   //2 propagate all tracks into ECAL, record its eta and phi
00336  
00337   for (std::multimap<float, edm::Ptr<reco::ConversionTrack> >::const_iterator tk_ref = allTracks.begin(); tk_ref != allTracks.end(); ++tk_ref ){
00338     const reco::Track* tk = tk_ref->second->trackRef().get()  ;
00339     
00340     
00341     //check impact position then match with BC
00342     math::XYZPointF ew;
00343     if ( getTrackImpactPosition(tk, trackerGeom, magField, ew) ){
00344       trackImpactPosition[tk_ref->second] = ew;
00345       
00346       reco::CaloClusterPtr closest_bc;//the closest matching BC to track
00347       
00348       if ( getMatchedBC(basicClusterPtrs, ew, closest_bc) ){
00349         trackMatchedBC[tk_ref->second] = closest_bc;
00350       }
00351     }    
00352   }
00353   
00354   
00355   
00356   //3. pair up tracks: 
00357   //TODO it is k-Closest pair of point problem
00358   //std::cout << " allTracks.size() " <<  allTracks.size() << std::endl;
00359   for(std::multimap<float, edm::Ptr<reco::ConversionTrack> >::const_iterator ll = allTracks.begin(); ll != allTracks.end();  ++ll ) {
00360     bool track1HighPurity=true;
00361     //std::cout << " Loop on allTracks " << std::endl;
00362     const  edm::RefToBase<reco::Track> & left = ll->second->trackRef();
00363     
00364 
00365     //TODO: This is a workaround, should be fixed with a proper function in the TTBuilder
00366     //(Note that the TrackRef and GsfTrackRef versions of the constructor are needed
00367     // to properly get refit tracks in the output vertex)
00368     reco::TransientTrack ttk_l;
00369     if (dynamic_cast<const reco::GsfTrack*>(left.get())) {
00370       ttk_l = thettbuilder_->build(left.castTo<reco::GsfTrackRef>());
00371     }
00372     else {
00373       ttk_l = thettbuilder_->build(left.castTo<reco::TrackRef>());
00374     }
00375     
00377     //      if ((allowTrackBC_ && !trackValidECAL[ll-allTracks.begin()]) )//this Left leg should have valid BC
00378     // continue;
00379     
00380     
00381     if (the_pvtx.isValid()){
00382       if (!(trackD0Cut(left, the_pvtx)))   track1HighPurity=false;
00383     } else {
00384       if (!(trackD0Cut(left)))  track1HighPurity=false;
00385     }
00386     
00387     std::vector<int> right_candidates;//store all right legs passed the cut (theta/approach and ref pair)
00388     std::vector<double> right_candidate_theta, right_candidate_approach;
00389     std::vector<std::pair<bool, reco::Vertex> > vertex_candidates;
00390     
00391     //inner loop only over tracks between eta and eta + deltaEta of the first track
00392     float etasearch = ll->first + deltaEta_;
00393     std::multimap<float, edm::Ptr<reco::ConversionTrack> >::const_iterator rr = ll;
00394     ++rr;    
00395     for (; rr != allTracks.lower_bound(etasearch); ++rr ) {
00396       bool track2HighPurity = true;
00397       bool highPurityPair = true;
00398       
00399       const  edm::RefToBase<reco::Track> & right = rr->second->trackRef();
00400       
00401       
00402       //TODO: This is a workaround, should be fixed with a proper function in the TTBuilder
00403       reco::TransientTrack ttk_r;
00404       if (dynamic_cast<const reco::GsfTrack*>(right.get())) {
00405         ttk_r = thettbuilder_->build(right.castTo<reco::GsfTrackRef>());
00406       }
00407       else {
00408         ttk_r = thettbuilder_->build(right.castTo<reco::TrackRef>());
00409       }
00410       //std::cout << " This track is " <<  right->algoName() << std::endl;
00411       
00412       
00413       //all vertexing preselection should go here
00414       
00415       //check for opposite charge
00416       if (  allowOppCharge_ && (left->charge()*right->charge() > 0) )  
00417         continue; //same sign, reject pair
00418           
00420       //if ( (allowTrackBC_ && !trackValidECAL[rr-allTracks.begin()] && rightBC_) )// if right track matches ECAL
00421       //  continue;
00422           
00423           
00424       double approachDist = -999.;
00425       //apply preselection to track pair, overriding preselection for gsf+X or ecalseeded+X pairs if so configured
00426       bool preselected = preselectTrackPair(ttk_l,ttk_r, approachDist);
00427       preselected = preselected || (bypassPreselGsf_ && (left->algo()==29 || right->algo()==29));
00428       preselected = preselected || (bypassPreselEcal_ && (left->algo()==15 || right->algo()==15 || left->algo()==16 || right->algo()==16));
00429       preselected = preselected || (bypassPreselEcalEcal_ && (left->algo()==15 || left->algo()==16) && (right->algo()==15 || right->algo()==16));
00430       
00431       if (!preselected) {
00432         continue;
00433       }
00434           
00435       //do the actual vertex fit
00436       reco::Vertex theConversionVertex;//by default it is invalid          
00437       bool goodVertex = checkVertex(ttk_l, ttk_r, magField, theConversionVertex);
00438           
00439       //bail as early as possible in case the fit didn't return a good vertex
00440       if (!goodVertex) {
00441         continue;
00442       }
00443 
00444           
00445           
00446       //track pair pass the quality cut
00447       if (   !( (trackQualityFilter(left, true) && trackQualityFilter(right, false))
00448                 || (trackQualityFilter(left, false) && trackQualityFilter(right, true)) ) ) {
00449         highPurityPair=false;
00450       }
00451 
00452       if (the_pvtx.isValid()){
00453         if (!(trackD0Cut(right, the_pvtx))) track2HighPurity=false; 
00454       } else {
00455         if (!(trackD0Cut(right))) track2HighPurity=false; 
00456       }
00457         
00458 
00459       //if all cuts passed, go ahead to make conversion candidates
00460       std::vector<edm::RefToBase<reco::Track> > trackPairRef;
00461       trackPairRef.push_back(left);//left track
00462       trackPairRef.push_back(right);//right track
00463 
00464       std::vector<math::XYZVectorF> trackPin;
00465       std::vector<math::XYZVectorF> trackPout;
00466       std::vector<math::XYZPointF> trackInnPos;
00467       std::vector<uint8_t> nHitsBeforeVtx;
00468       std::vector<Measurement1DFloat> dlClosestHitToVtx;
00469 
00470       if (left->extra().isNonnull() && right->extra().isNonnull()){//only available on TrackExtra
00471         trackInnPos.push_back(  toFConverterP(left->innerPosition()));
00472         trackInnPos.push_back(  toFConverterP(right->innerPosition()));
00473         trackPin.push_back(toFConverterV(left->innerMomentum()));
00474         trackPin.push_back(toFConverterV(right->innerMomentum()));
00475         trackPout.push_back(toFConverterV(left->outerMomentum()));
00476         trackPout.push_back(toFConverterV(right->outerMomentum()));
00477       }
00478           
00479       if (ll->second->trajRef().isNonnull() && rr->second->trajRef().isNonnull()) {
00480         std::pair<uint8_t,Measurement1DFloat> leftWrongHits = hitChecker.nHitsBeforeVtx(*ll->second->trajRef().get(),theConversionVertex);
00481         std::pair<uint8_t,Measurement1DFloat> rightWrongHits = hitChecker.nHitsBeforeVtx(*rr->second->trajRef().get(),theConversionVertex);
00482         nHitsBeforeVtx.push_back(leftWrongHits.first);
00483         nHitsBeforeVtx.push_back(rightWrongHits.first);
00484         dlClosestHitToVtx.push_back(leftWrongHits.second);
00485         dlClosestHitToVtx.push_back(rightWrongHits.second);            
00486       }
00487           
00488       uint8_t nSharedHits = hitChecker.nSharedHits(*left.get(),*right.get());
00489 
00490 
00491       //if using kinematic fit, check with chi2 post cut
00492       if (theConversionVertex.isValid()){
00493         const float chi2Prob = ChiSquaredProbability(theConversionVertex.chi2(), theConversionVertex.ndof());
00494         if (chi2Prob<vtxChi2_)  highPurityPair=false;
00495       }
00496 
00497       //std::cout << "  highPurityPair after vertex cut " <<  highPurityPair << std::endl;
00498       std::vector<math::XYZPointF> trkPositionAtEcal;
00499       std::vector<reco::CaloClusterPtr> matchingBC;
00500 
00501       if (allowTrackBC_){//TODO find out the BC ptrs if not doing matching, otherwise, leave it empty
00502         //const int lbc_handle = bcHandleId[ll-allTracks.begin()],
00503         //            rbc_handle = bcHandleId[rr-allTracks.begin()];
00504 
00505         std::map<edm::Ptr<reco::ConversionTrack>, math::XYZPointF>::const_iterator trackImpactPositionLeft = trackImpactPosition.find(ll->second);
00506         std::map<edm::Ptr<reco::ConversionTrack>, math::XYZPointF>::const_iterator trackImpactPositionRight = trackImpactPosition.find(rr->second);
00507         std::map<edm::Ptr<reco::ConversionTrack>, reco::CaloClusterPtr>::const_iterator trackMatchedBCLeft = trackMatchedBC.find(ll->second);        
00508         std::map<edm::Ptr<reco::ConversionTrack>, reco::CaloClusterPtr>::const_iterator trackMatchedBCRight = trackMatchedBC.find(rr->second);        
00509         
00510         if (trackImpactPositionLeft!=trackImpactPosition.end()) {
00511           trkPositionAtEcal.push_back(trackImpactPositionLeft->second);//left track
00512         }
00513         else {
00514           trkPositionAtEcal.push_back(math::XYZPointF());//left track
00515         }
00516         if (trackImpactPositionRight!=trackImpactPosition.end()) {//second track ECAL position may be invalid
00517           trkPositionAtEcal.push_back(trackImpactPositionRight->second);
00518         }
00519 
00520         double total_e_bc = 0.;
00521         if (trackMatchedBCLeft!=trackMatchedBC.end()) {
00522           matchingBC.push_back(trackMatchedBCLeft->second);//left track
00523           total_e_bc += trackMatchedBCLeft->second->energy();
00524         }
00525         else {
00526           matchingBC.push_back( reco::CaloClusterPtr() );//left track
00527         }
00528         if (trackMatchedBCRight!=trackMatchedBC.end()) {//second track ECAL position may be invalid
00529           matchingBC.push_back(trackMatchedBCRight->second);
00530           total_e_bc += trackMatchedBCRight->second->energy();
00531         }
00532         
00533         if (total_e_bc<energyTotalBC_) {
00534           highPurityPair = false;
00535         }
00536 
00537 
00538       }
00539       //signature cuts, then check if vertex, then post-selection cuts
00540       highPurityPair = highPurityPair && track1HighPurity && track2HighPurity && goodVertex && checkPhi(left, right, trackerGeom, magField, theConversionVertex) ;
00541 
00542 
00544       /*
00545         for ( std::vector<edm::RefToBase<reco::Track> >::iterator iTk=trackPairRef.begin();  iTk!=trackPairRef.end(); iTk++) {
00546         math::XYZPointF impPos;
00547         if ( getTrackImpactPosition(*iTk, trackerGeom, magField, impPos) ) {
00548                
00549         }
00550 
00551         }
00552       */
00553 
00554       const float minAppDist = approachDist;
00555       reco::Conversion::ConversionAlgorithm algo = reco::Conversion::algoByName(algoName_);
00556       float dummy=0;
00557       reco::CaloClusterPtrVector scPtrVec;
00558       reco::Conversion  newCandidate(scPtrVec,  trackPairRef, trkPositionAtEcal, theConversionVertex, matchingBC, minAppDist,  trackInnPos, trackPin, trackPout, nHitsBeforeVtx, dlClosestHitToVtx, nSharedHits, dummy, algo );
00559       // Fill in scPtrVec with the macthing SC
00560       if ( matchingSC ( superClusterPtrs, newCandidate, scPtrVec) ) 
00561         newCandidate.setMatchingSuperCluster( scPtrVec);
00562           
00563       //std::cout << " ConversionProducer  scPtrVec.size " <<  scPtrVec.size() << std::endl;
00564           
00565       newCandidate.setQuality(reco::Conversion::highPurity,  highPurityPair);
00566       bool generalTracksOnly = ll->second->isTrackerOnly() && rr->second->isTrackerOnly() && !dynamic_cast<const reco::GsfTrack*>(ll->second->trackRef().get()) && !dynamic_cast<const reco::GsfTrack*>(rr->second->trackRef().get());
00567       bool arbitratedEcalSeeded = ll->second->isArbitratedEcalSeeded() && rr->second->isArbitratedEcalSeeded();
00568       bool arbitratedMerged = ll->second->isArbitratedMerged() && rr->second->isArbitratedMerged();
00569       bool arbitratedMergedEcalGeneral = ll->second->isArbitratedMergedEcalGeneral() && rr->second->isArbitratedMergedEcalGeneral();          
00570           
00571       newCandidate.setQuality(reco::Conversion::generalTracksOnly,  generalTracksOnly);
00572       newCandidate.setQuality(reco::Conversion::arbitratedEcalSeeded,  arbitratedEcalSeeded);
00573       newCandidate.setQuality(reco::Conversion::arbitratedMerged,  arbitratedMerged);
00574       newCandidate.setQuality(reco::Conversion::arbitratedMergedEcalGeneral,  arbitratedMergedEcalGeneral);          
00575           
00576       outputConvPhotonCollection.push_back(newCandidate);
00577 
00578     }
00579       
00580   }
00581 
00582 
00583 
00584 
00585 
00586 
00587 }
00588 
00589 
00590 
00591 
00592 
00593 //
00594 // member functions
00595 //
00596 
00597 inline bool ConversionProducer::trackQualityFilter(const  edm::RefToBase<reco::Track>&   ref, bool isLeft){
00598   bool pass = true;
00599   if (isLeft){
00600     pass = (ref->normalizedChi2() < maxChi2Left_ && ref->found() >= minHitsLeft_);
00601   } else {
00602     pass = (ref->normalizedChi2() < maxChi2Right_ && ref->found() >= minHitsRight_);
00603   }
00604 
00605   return pass;
00606 }
00607 
00608 inline bool ConversionProducer::trackD0Cut(const  edm::RefToBase<reco::Track>&  ref){
00609   //NOTE if not allow d0 cut, always true
00610   return ((!allowD0_) || !(ref->d0()*ref->charge()/ref->d0Error()<d0Cut_));
00611 }
00612 
00613 inline bool ConversionProducer::trackD0Cut(const edm::RefToBase<reco::Track>&  ref, const reco::Vertex& the_pvtx){
00614   //
00615   return ((!allowD0_) || !(-ref->dxy(the_pvtx.position())*ref->charge()/ref->dxyError()<d0Cut_));
00616 }
00617 
00618 
00619 bool ConversionProducer::getTrackImpactPosition(const reco::Track* tk_ref,
00620                                                 const TrackerGeometry* trackerGeom, const MagneticField* magField,
00621                                                 math::XYZPointF& ew){
00622 
00623   PropagatorWithMaterial propag( alongMomentum, 0.000511, magField );
00624   
00625   ReferenceCountingPointer<Surface> ecalWall(
00626                                              new  BoundCylinder( GlobalPoint(0.,0.,0.), TkRotation<float>(),
00627                                                                  SimpleCylinderBounds( 129, 129, -320.5, 320.5 ) ) );
00628   const float epsilon = 0.001;
00629   Surface::RotationType rot; // unit rotation matrix
00630   const float barrelRadius = 129.f;
00631   const float barrelHalfLength = 270.9f;
00632   const float endcapRadius = 171.1f;
00633   const float endcapZ = 320.5f;
00634   ReferenceCountingPointer<BoundCylinder>  theBarrel_(new BoundCylinder( Surface::PositionType(0,0,0), rot,
00635                                                                          SimpleCylinderBounds( barrelRadius-epsilon, barrelRadius+epsilon, -barrelHalfLength, barrelHalfLength)));
00636   ReferenceCountingPointer<BoundDisk>      theNegativeEtaEndcap_(
00637                                                                  new BoundDisk( Surface::PositionType( 0, 0, -endcapZ), rot,
00638                                                                                 SimpleDiskBounds( 0, endcapRadius, -epsilon, epsilon)));
00639   ReferenceCountingPointer<BoundDisk>      thePositiveEtaEndcap_(
00640                                                                  new BoundDisk( Surface::PositionType( 0, 0, endcapZ), rot,
00641                                                                                 SimpleDiskBounds( 0, endcapRadius, -epsilon, epsilon)));
00642 
00643   //const TrajectoryStateOnSurface myTSOS = trajectoryStateTransform::innerStateOnSurface(*(*ref), *trackerGeom, magField);
00644   const TrajectoryStateOnSurface myTSOS = trajectoryStateTransform::outerStateOnSurface(*tk_ref, *trackerGeom, magField);
00645   TrajectoryStateOnSurface  stateAtECAL;
00646   stateAtECAL = propag.propagate(myTSOS, *theBarrel_);
00647   if (!stateAtECAL.isValid() || ( stateAtECAL.isValid() && fabs(stateAtECAL.globalPosition().eta() ) >1.479 )  ) {
00648     //endcap propagator
00649     if (myTSOS.globalPosition().eta() > 0.) {
00650             stateAtECAL = propag.propagate(myTSOS, *thePositiveEtaEndcap_);
00651     } else {
00652             stateAtECAL = propag.propagate(myTSOS, *theNegativeEtaEndcap_);
00653     }
00654   }       
00655   if (stateAtECAL.isValid()){
00656     ew = stateAtECAL.globalPosition();
00657     return true;
00658   }
00659   else
00660     return false;
00661 }
00662 
00663 
00664 
00665 
00666 bool ConversionProducer::matchingSC(const std::multimap<double, reco::CaloClusterPtr>& scMap, 
00667                                     reco::Conversion& aConv,
00668                                     // reco::CaloClusterPtr& mSC){
00669                                     reco::CaloClusterPtrVector& mSC) {
00670 
00671   //  double dRMin=999.;
00672   double detaMin=999.;
00673   double dphiMin=999.;                 
00674   reco::CaloClusterPtr match;
00675   for (std::multimap<double, reco::CaloClusterPtr>::const_iterator scItr = scMap.begin();  scItr != scMap.end(); scItr++) {
00676     const reco::CaloClusterPtr& sc = scItr->second; 
00677     const double delta_phi = reco::deltaPhi( aConv.refittedPairMomentum().phi(), sc->phi());
00678     double sceta = sc->eta();
00679     double conveta = etaTransformation(aConv.refittedPairMomentum().eta(), aConv.zOfPrimaryVertexFromTracks() );
00680     const double delta_eta = fabs(conveta - sceta);
00681     if ( fabs(delta_eta) < fabs(detaMin) && fabs(delta_phi) < fabs(dphiMin) ) {
00682       detaMin=  fabs(delta_eta);
00683       dphiMin=  fabs(delta_phi);
00684       match=sc;
00685     }
00686   }
00687   
00688   if ( fabs(detaMin) < dEtacutForSCmatching_ && fabs(dphiMin) < dPhicutForSCmatching_ ) {
00689     mSC.push_back(match);
00690     return true;
00691   } else 
00692     return false;
00693 }
00694 
00695 bool ConversionProducer::getMatchedBC(const std::multimap<double, reco::CaloClusterPtr>& bcMap, 
00696                                       const math::XYZPointF& trackImpactPosition,
00697                                       reco::CaloClusterPtr& closestBC){
00698   const double track_eta = trackImpactPosition.eta();
00699   const double track_phi = trackImpactPosition.phi();
00700 
00701   double min_eta = 999., min_phi = 999.;
00702   reco::CaloClusterPtr closest_bc;
00703   for (std::multimap<double, reco::CaloClusterPtr>::const_iterator bc = bcMap.lower_bound(track_eta - halfWayEta_);
00704        bc != bcMap.upper_bound(track_eta + halfWayEta_); ++bc){//use eta map to select possible BC collection then loop in
00705     const reco::CaloClusterPtr& ebc = bc->second;
00706     const double delta_eta = track_eta-(ebc->position().eta());
00707     const double delta_phi = reco::deltaPhi(track_phi, (ebc->position().phi()));
00708     if (fabs(delta_eta)<dEtaTkBC_ && fabs(delta_phi)<dPhiTkBC_){
00709             if (fabs(min_eta)>fabs(delta_eta) && fabs(min_phi)>fabs(delta_phi)){//take the closest to track BC
00710         min_eta = delta_eta;
00711         min_phi = delta_phi;
00712         closest_bc = bc->second;
00713         //TODO check if min_eta>delta_eta but min_phi<delta_phi
00714             }
00715     }
00716   }
00717 
00718   if (min_eta < 999.){
00719     closestBC = closest_bc;
00720     return true;
00721   } else
00722     return false;
00723 }
00724 
00725 
00726 
00727 
00728 
00729 
00730 //check track open angle of phi at vertex
00731 bool ConversionProducer::checkPhi(const edm::RefToBase<reco::Track>& tk_l, const edm::RefToBase<reco::Track>& tk_r,
00732                                   const TrackerGeometry* trackerGeom, const MagneticField* magField,
00733                                   const reco::Vertex& vtx){
00734   if (!allowDeltaPhi_)
00735     return true;
00736   //if track has innermost momentum, check with innermost phi
00737   //if track also has valid vertex, propagate to vertex then calculate phi there
00738   //if track has no innermost momentum, just return true, because track->phi() makes no sense
00739   if (tk_l->extra().isNonnull() && tk_r->extra().isNonnull()){
00740     double iphi1 = tk_l->innerMomentum().phi(), iphi2 = tk_r->innerMomentum().phi();
00741     if (vtx.isValid()){
00742             PropagatorWithMaterial propag( anyDirection, 0.000511, magField );
00743             
00744             double recoPhoR = vtx.position().Rho();
00745             Surface::RotationType rot;
00746             ReferenceCountingPointer<BoundCylinder>  theBarrel_(new BoundCylinder( Surface::PositionType(0,0,0), rot,
00747                                                                              SimpleCylinderBounds( recoPhoR-0.001, recoPhoR+0.001, -fabs(vtx.position().z()), fabs(vtx.position().z()))));
00748             ReferenceCountingPointer<BoundDisk>      theDisk_(
00749                                                         new BoundDisk( Surface::PositionType( 0, 0, vtx.position().z()), rot,
00750                                                                        SimpleDiskBounds( 0, recoPhoR, -0.001, 0.001)));
00751 
00752             const TrajectoryStateOnSurface myTSOS1 = trajectoryStateTransform::innerStateOnSurface(*tk_l, *trackerGeom, magField);
00753             const TrajectoryStateOnSurface myTSOS2 = trajectoryStateTransform::innerStateOnSurface(*tk_r, *trackerGeom, magField);
00754             TrajectoryStateOnSurface  stateAtVtx1, stateAtVtx2;
00755             stateAtVtx1 = propag.propagate(myTSOS1, *theBarrel_);
00756             if (!stateAtVtx1.isValid() ) {
00757         stateAtVtx1 = propag.propagate(myTSOS1, *theDisk_);
00758             }
00759             if (stateAtVtx1.isValid()){
00760         iphi1 = stateAtVtx1.globalDirection().phi();
00761             }
00762             stateAtVtx2 = propag.propagate(myTSOS2, *theBarrel_);
00763             if (!stateAtVtx2.isValid() ) {
00764         stateAtVtx2 = propag.propagate(myTSOS2, *theDisk_);
00765             }
00766             if (stateAtVtx2.isValid()){
00767         iphi2 = stateAtVtx2.globalDirection().phi();
00768             }
00769     }
00770     const double dPhi = reco::deltaPhi(iphi1, iphi2);
00771     return (fabs(dPhi) < deltaPhi_);
00772   } else {
00773     return true;
00774   }
00775 }
00776 
00777 bool ConversionProducer::preselectTrackPair(const reco::TransientTrack &ttk_l, const reco::TransientTrack &ttk_r,
00778                                             double& appDist) {
00779   
00780 
00781   double dCotTheta =  1./tan(ttk_l.track().innerMomentum().theta()) - 1./tan(ttk_r.track().innerMomentum().theta());
00782   if (allowDeltaCot_ && (std::abs(dCotTheta) > deltaCotTheta_)) {
00783     return false;
00784   }
00785   
00786   //non-conversion hypothesis, reject prompt track pairs
00787   ClosestApproachInRPhi closest;
00788   closest.calculate(ttk_l.innermostMeasurementState(),ttk_r.innermostMeasurementState());
00789   if (!closest.status()) {
00790     return false;
00791   }
00792   
00793   if (closest.crossingPoint().perp() < r_cut) {
00794     return false;
00795   }
00796 
00797   
00798   //compute tangent point btw tracks (conversion hypothesis)
00799   TangentApproachInRPhi tangent;
00800   tangent.calculate(ttk_l.innermostMeasurementState(),ttk_r.innermostMeasurementState());
00801   if (!tangent.status()) {
00802     return false;
00803   }
00804   
00805   GlobalPoint tangentPoint = tangent.crossingPoint();
00806   double rho = tangentPoint.perp();
00807   
00808   //reject candidates well outside of tracker bounds
00809   if (rho > maxTrackRho_) {
00810     return false;
00811   }
00812   
00813   if (std::abs(tangentPoint.z()) > maxTrackZ_) {
00814     return false;
00815   }
00816   
00817   std::pair<GlobalTrajectoryParameters,GlobalTrajectoryParameters> trajs = tangent.trajectoryParameters();
00818   
00819   //very large separation in z, no hope
00820   if (std::abs(trajs.first.position().z() - trajs.second.position().z()) > dzCut_) {
00821     return false;
00822   }
00823   
00824   
00825   float minApproach = tangent.perpdist();
00826   appDist = minApproach;
00827   
00828   if (allowMinApproach_ && (minApproach < minApproachLow_ || minApproach > minApproachHigh_) ) {
00829     return false;
00830   }
00831   
00832   return true;
00833   
00834   
00835 }
00836 
00837 bool ConversionProducer::checkTrackPair(const std::pair<edm::RefToBase<reco::Track>, reco::CaloClusterPtr>& ll, 
00838                                         const std::pair<edm::RefToBase<reco::Track>, reco::CaloClusterPtr>& rr){
00839 
00840   const reco::CaloClusterPtr& bc_l = ll.second;//can be null, so check isNonnull()
00841   const reco::CaloClusterPtr& bc_r = rr.second;
00842     
00843   //The cuts should be ordered by considering if takes time and if cuts off many fakes
00844   if (allowTrackBC_){
00845     //check energy of BC
00846     double total_e_bc = 0;
00847     if (bc_l.isNonnull()) total_e_bc += bc_l->energy();
00848     if (rightBC_) 
00849             if (bc_r.isNonnull())
00850         total_e_bc += bc_r->energy();
00851 
00852     if (total_e_bc < energyTotalBC_) return false;
00853   }
00854 
00855   return true;
00856 }
00857 
00858 
00859 
00860 //because reco::vertex uses track ref, so have to keep them
00861 bool ConversionProducer::checkVertex(const reco::TransientTrack &ttk_l, const reco::TransientTrack &ttk_r, 
00862                                      const MagneticField* magField,
00863                                      reco::Vertex& the_vertex){
00864   bool found = false;
00865 
00866   std::vector<reco::TransientTrack>  pair;
00867   pair.push_back(ttk_l);
00868   pair.push_back(ttk_r);
00869    
00870   found = theVertexFinder_->run(pair, the_vertex);
00871 
00872 
00873 
00874   return found;
00875 }
00876 
00877 
00878 
00879 double ConversionProducer::etaTransformation(  float EtaParticle , float Zvertex)  {
00880 
00881   //---Definitions
00882   const float PI    = 3.1415927;
00883 
00884   //---Definitions for ECAL
00885   const float R_ECAL           = 136.5;
00886   const float Z_Endcap         = 328.0;
00887   const float etaBarrelEndcap  = 1.479; 
00888    
00889   //---ETA correction
00890 
00891   float Theta = 0.0  ; 
00892   float ZEcal = R_ECAL*sinh(EtaParticle)+Zvertex;
00893 
00894   if(ZEcal != 0.0) Theta = atan(R_ECAL/ZEcal);
00895   if(Theta<0.0) Theta = Theta+PI ;
00896   double ETA = - log(tan(0.5*Theta));
00897          
00898   if( fabs(ETA) > etaBarrelEndcap )
00899     {
00900       float Zend = Z_Endcap ;
00901       if(EtaParticle<0.0 )  Zend = -Zend ;
00902       float Zlen = Zend - Zvertex ;
00903       float RR = Zlen/sinh(EtaParticle); 
00904       Theta = atan(RR/Zend);
00905       if(Theta<0.0) Theta = Theta+PI ;
00906       ETA = - log(tan(0.5*Theta));                    
00907     } 
00908   //---Return the result
00909   return ETA;
00910   //---end
00911 }
00912