CMS 3D CMS Logo

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