CMS 3D CMS Logo

SoftConversionProducer.cc

Go to the documentation of this file.
00001 #include <iostream>
00002 #include <vector>
00003 #include <memory>
00004 
00005 // Framework
00006 #include "DataFormats/Common/interface/Handle.h"
00007 #include "FWCore/Framework/interface/ESHandle.h"
00008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00009 #include "FWCore/Utilities/interface/Exception.h"
00010 //
00011 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00012 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00013 #include "RecoTracker/Record/interface/CkfComponentsRecord.h"
00014 #include "RecoTracker/Record/interface/TrackerRecoGeometryRecord.h"
00015 //
00016 #include "DataFormats/TrackReco/interface/Track.h"
00017 #include "DataFormats/TrackReco/interface/TrackBase.h"
00018 #include "DataFormats/TrackReco/interface/TrackExtra.h"
00019 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00020 #include "TrackingTools/TransientTrack/interface/TrackTransientTrack.h"
00021 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00022 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00023 
00024 #include "DataFormats/EgammaCandidates/interface/Conversion.h"
00025 #include "DataFormats/EgammaCandidates/interface/ConversionFwd.h"
00026 #include "RecoEgamma/EgammaPhotonAlgos/interface/ConversionTrackEcalImpactPoint.h"
00027 #include "RecoEgamma/EgammaPhotonAlgos/interface/ConversionTrackPairFinder.h"
00028 #include "RecoEgamma/EgammaPhotonAlgos/interface/ConversionVertexFinder.h"
00029 
00030 #include "DataFormats/EgammaTrackReco/interface/TrackCaloClusterAssociation.h"
00031 #include "RecoEgamma/EgammaPhotonProducers/interface/SoftConversionProducer.h"
00032 
00033 #include "Math/GenVector/VectorUtil.h"
00034 
00035 
00036 SoftConversionProducer::SoftConversionProducer(const edm::ParameterSet& config) : conf_(config) {
00037 
00038   LogDebug("SoftConversionProducer") << " SoftConversionProducer CTOR " << "\n";
00039   
00040   conversionOITrackProducer_              = conf_.getParameter<std::string>("conversionOITrackProducer");
00041   conversionIOTrackProducer_              = conf_.getParameter<std::string>("conversionIOTrackProducer");
00042   outInTrackClusterAssociationCollection_ = conf_.getParameter<std::string>("outInTrackClusterAssociationCollection");
00043   inOutTrackClusterAssociationCollection_ = conf_.getParameter<std::string>("inOutTrackClusterAssociationCollection");
00044   clusterType_                            = conf_.getParameter<std::string>("clusterType");
00045   clusterBarrelCollection_                = conf_.getParameter<edm::InputTag>("clusterBarrelCollection");
00046   clusterEndcapCollection_                = conf_.getParameter<edm::InputTag>("clusterEndcapCollection");
00047   softConversionCollection_               = conf_.getParameter<std::string>("softConversionCollection");
00048   trackMaxChi2_                           = conf_.getParameter<double>("trackMaxChi2");
00049   trackMinHits_                           = conf_.getParameter<double>("trackMinHits");
00050   clustersMaxDeltaEta_                    = conf_.getParameter<double>("clustersMaxDeltaEta");
00051   clustersMaxDeltaPhi_                    = conf_.getParameter<double>("clustersMaxDeltaPhi");
00052   
00053   theTrackPairFinder_ = 0;
00054   theVertexFinder_ = 0;
00055   theEcalImpactPositionFinder_ = 0;
00056   
00057   // Register the product
00058   produces< reco::ConversionCollection >(softConversionCollection_);
00059 
00060 }
00061 
00062 
00063 SoftConversionProducer::~SoftConversionProducer() {
00064 
00065   if(theTrackPairFinder_) delete theTrackPairFinder_;
00066   if(theVertexFinder_) delete theVertexFinder_;
00067   if(theEcalImpactPositionFinder_) delete theEcalImpactPositionFinder_; 
00068 
00069 }
00070 
00071 
00072 void  SoftConversionProducer::beginJob (edm::EventSetup const & theEventSetup) {
00073   
00074   //get magnetic field
00075   theEventSetup.get<IdealMagneticFieldRecord>().get(theMF_);  
00076   
00077   // instantiate the Track Pair Finder algorithm
00078   theTrackPairFinder_ = new ConversionTrackPairFinder ();
00079 
00080   // instantiate the Vertex Finder algorithm
00081   theVertexFinder_ = new ConversionVertexFinder ();
00082   
00083   // instantiate the algorithm for finding the position of the track extrapolation at the Ecal front face
00084   theEcalImpactPositionFinder_ = new   ConversionTrackEcalImpactPoint ( &(*theMF_) );
00085 
00086 }
00087 
00088 
00089 void  SoftConversionProducer::endJob () {}
00090 
00091 
00092 void SoftConversionProducer::produce(edm::Event& theEvent, const edm::EventSetup& theEventSetup) {
00093   
00094   edm::Handle<reco::TrackCollection> outInTrkHandle;
00095   theEvent.getByLabel(conversionOITrackProducer_,  outInTrkHandle);
00096    
00097   edm::Handle<reco::TrackCaloClusterPtrAssociation> outInTrkClusterAssocHandle;
00098   theEvent.getByLabel( conversionOITrackProducer_, outInTrackClusterAssociationCollection_, outInTrkClusterAssocHandle);
00099 
00100   edm::Handle<reco::TrackCollection> inOutTrkHandle;
00101   theEvent.getByLabel(conversionIOTrackProducer_, inOutTrkHandle);
00102   
00103   edm::Handle<reco::TrackCaloClusterPtrAssociation> inOutTrkClusterAssocHandle;
00104   theEvent.getByLabel( conversionIOTrackProducer_, inOutTrackClusterAssociationCollection_, inOutTrkClusterAssocHandle);
00105 
00106   edm::Handle<edm::View<reco::CaloCluster> > clusterBarrelHandle;
00107   theEvent.getByLabel(clusterBarrelCollection_, clusterBarrelHandle);
00108 
00109   edm::Handle<edm::View<reco::CaloCluster> > clusterEndcapHandle;
00110   if(clusterType_ == "BasicCluster"){
00111     //theEvent.getByLabel(clusterEndcapProducer_, clusterEndcapCollection_, clusterEndcapHandle);
00112     theEvent.getByLabel(clusterEndcapCollection_, clusterEndcapHandle);
00113   }
00114 
00115   // create temporary map to loop over tracks conveniently
00116   TrackClusterMap trackClusterMap;
00117 
00118   int nTracksOI = (int) outInTrkHandle->size();
00119   //  std::cout << "  nTracksOI " <<  nTracksOI << std::endl;
00120   for(int itrk=0; itrk<nTracksOI; itrk++){
00121     reco::TrackRef tRef(outInTrkHandle,itrk);
00122     if(!trackQualityCut(tRef)) continue;
00123     reco::CaloClusterPtr cRef = (*outInTrkClusterAssocHandle)[tRef];
00124     trackClusterMap.push_back(make_pair(tRef,cRef));
00125   }
00126 
00127   int nTracksIO = (int) inOutTrkHandle->size();
00128   // std::cout << "  nTracksIO " <<  nTracksIO << std::endl;
00129   for(int itrk=0; itrk<nTracksIO; itrk++){
00130     reco::TrackRef tRef(inOutTrkHandle,itrk);
00131     if(!trackQualityCut(tRef)) continue;
00132     reco::CaloClusterPtr cRef = (*inOutTrkClusterAssocHandle)[tRef];
00133     trackClusterMap.push_back(make_pair(tRef,cRef));
00134   }
00135 
00136   // Transform Track into TransientTrack (needed by the Vertex fitter)
00137   edm::ESHandle<TransientTrackBuilder> theTransientTrackBuilder;
00138   theEventSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theTransientTrackBuilder);
00139 
00140   // the output collection to be produced from this producer
00141   std::auto_ptr< reco::ConversionCollection > outputColl(new reco::ConversionCollection);
00142 
00143   // prepare iterator
00144   TrackClusterMap::iterator iter1    = trackClusterMap.begin();
00145   TrackClusterMap::iterator iter2    = trackClusterMap.begin();
00146   TrackClusterMap::iterator iter_end = trackClusterMap.end();
00147 
00148 
00149 
00150   // double-loop to make pairs
00151   for( ; iter1 != iter_end; iter1++) {
00152     //    std::cout << " back to start of loop 1 " << std::endl;
00153     const reco::TrackRef trk1 = iter1->first;
00154 
00155     for(iter2 = iter1+1; iter2 != iter_end; iter2++) {
00156       //      std::cout << " back to start of loop 2 " << std::endl;
00157       const reco::TrackRef trk2 = iter2->first;
00158       if(trk1 == trk2) continue;
00159 
00160 
00161       const reco::CaloClusterPtr cls1 = iter1->second;
00162       reco::TransientTrack tsk1 = theTransientTrackBuilder->build(*trk1);
00163 
00164       const reco::CaloClusterPtr cls2 = iter2->second;
00165       reco::TransientTrack tsk2 = theTransientTrackBuilder->build(*trk2);
00166 
00167       //      std::cout << " after transient " << std::endl;
00168       // std::cout << " eta1 " << cls1->position().Eta() << " eta2 " << cls2->position().Eta() << std::endl;
00169 
00170       if ( ( tsk1.innermostMeasurementState().globalParameters().position() - tsk2.innermostMeasurementState().globalParameters().position() ).mag() < 1e-7 &&
00171            ( tsk1.innermostMeasurementState().globalParameters().momentum() - tsk2.innermostMeasurementState().globalParameters().momentum() ).mag() < 1e-7 ) continue;
00172 
00173 
00174       double dEta = std::abs(cls1->position().Eta() - cls2->position().Eta());
00175       if(dEta > clustersMaxDeltaEta_) continue;
00176       double dPhi = std::abs(ROOT::Math::VectorUtil::DeltaPhi(cls1->position(),cls2->position()));
00177       if(dPhi > clustersMaxDeltaPhi_) continue;
00178 
00179       std::vector<reco::TransientTrack> toBeFitted;
00180       toBeFitted.push_back(tsk1);
00181       toBeFitted.push_back(tsk2);
00182 
00183       //      std::cout << " Before vertex " << std::endl;
00184       reco::Vertex theConversionVertex = (reco::Vertex) theVertexFinder_->run(toBeFitted);
00185       //std::cout << " After vertex " << std::endl;
00186 
00187       if(theConversionVertex.isValid()){
00188         //      std::cout << " valid vertex " << std::endl;
00189         reco::CaloClusterPtrVector scRefs;
00190         scRefs.push_back(cls1);
00191         scRefs.push_back(cls2);
00192         //std::cout << " after push back scref " << std::endl; 
00193         std::vector<reco::CaloClusterPtr> clusterRefs;
00194         clusterRefs.push_back(cls1);
00195         clusterRefs.push_back(cls2);
00196         //std::cout << " after push back slusterref " << std::endl; 
00197         std::vector<reco::TrackRef> trkRefs;
00198         trkRefs.push_back(trk1);
00199         trkRefs.push_back(trk2);
00200 
00201 
00202         //      std::cout << " Before impact finder " << std::endl;
00203         std::vector<math::XYZPoint> trkPositionAtEcal = theEcalImpactPositionFinder_->find(toBeFitted,clusterBarrelHandle);
00204         //      std::cout << " After impact finder " << std::endl;
00205         if((clusterType_ == "BasicCluster") && (std::abs(cls2->position().Eta()) > 1.5)){
00206           trkPositionAtEcal.clear();
00207           trkPositionAtEcal = theEcalImpactPositionFinder_->find(toBeFitted,clusterEndcapHandle);
00208         }
00209 
00210         reco::Conversion newCandidate(scRefs, trkRefs, trkPositionAtEcal, theConversionVertex, clusterRefs);
00211 
00212         // Check this candidate is already in the collection.
00213         // This is checking that two tracks in a conversion candidate are identicial.
00214 
00215         if(NotAlreadyIn(newCandidate,outputColl)) outputColl->push_back(newCandidate);
00216 
00217         //      printf("=====> run(%d), event(%d) <=====\n",theEvent.id().run(),theEvent.id().event());
00218         //      printf("Found a softConverion with vtxR(%f), vtxEta(%f), pt(%f), pt1(%f), pt2(%f)\n",
00219         //      newCandidate.conversionVertex().position().rho(),newCandidate.conversionVertex().position().eta(),
00220         //      newCandidate.pairMomentum().perp(),trk1->momentum().rho(),trk2->momentum().rho());
00221 
00222         clusterRefs.clear();
00223         trkRefs.clear();
00224         trkPositionAtEcal.clear();
00225       }// if(theConversionVertex.isValid()
00226 
00227       toBeFitted.clear();
00228 
00229     }// end of iter2
00230   }// end of iter1
00231 
00232   //  std::cout << " Putting stuff in the event " << std::endl;
00233   // put the collection into the event
00234   theEvent.put( outputColl, softConversionCollection_);
00235   
00236 }
00237 
00238 
00239 bool SoftConversionProducer::trackQualityCut(const reco::TrackRef& trk){
00240   return (trk->numberOfValidHits() >= trackMinHits_ && trk->normalizedChi2() < trackMaxChi2_);
00241 }
00242 
00243 
00244 bool SoftConversionProducer::NotAlreadyIn(const reco::Conversion& thisConv,
00245                                           const std::auto_ptr<reco::ConversionCollection>& outputColl) const {
00246 
00247   if(outputColl->size() == 0) return true;
00248 
00249   reco::ConversionCollection::const_iterator it = outputColl->begin();
00250   reco::ConversionCollection::const_iterator it_end = outputColl->end();
00251   for( ; it != it_end; it++){
00252     const reco::Conversion& conv = *it;
00253     if((thisConv.tracks()[0] == conv.tracks()[0] && thisConv.tracks()[1] == conv.tracks()[1]) ||
00254        (thisConv.tracks()[0] == conv.tracks()[1] && thisConv.tracks()[1] == conv.tracks()[0]))
00255       return false;
00256   }// for
00257   
00258   return true;
00259 }

Generated on Tue Jun 9 17:43:27 2009 for CMSSW by  doxygen 1.5.4