CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/CommonTools/RecoAlgos/interface/CosmicTrackingParticleSelector.h

Go to the documentation of this file.
00001 #ifndef RecoSelectors_CosmicTrackingParticleSelector_h
00002 #define RecoSelectors_CosmicTrackingParticleSelector_h
00003 /* \class CosmicTrackingParticleSelector
00004  *
00005  * \author Yanyan Gao, FNAL
00006  *
00007  *  $Date: 2010/02/11 00:10:48 $
00008  *  $Revision: 1.2 $
00009  *
00010  */
00011 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00012 #include "FWCore/Utilities/interface/InputTag.h"
00013 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/Framework/interface/ESHandle.h"
00016 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00017 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
00018 #include "TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h"
00019 #include "MagneticField/Engine/interface/MagneticField.h" 
00020 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00021 
00022 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00023 #include <Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h>
00024 #include <DataFormats/GeometrySurface/interface/Surface.h>
00025 #include <DataFormats/GeometrySurface/interface/GloballyPositioned.h>
00026 #include <Geometry/CommonDetUnit/interface/GeomDet.h>
00027 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00028 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00029 
00030 class TrajectoryStateClosestToBeamLineBuilder;
00031 
00032 class CosmicTrackingParticleSelector {
00033 
00034 public:
00035   typedef TrackingParticleCollection collection;
00036   typedef std::vector<const TrackingParticle *> container;
00037   typedef container::const_iterator const_iterator;
00038   
00039   CosmicTrackingParticleSelector(){}
00040 
00041   CosmicTrackingParticleSelector ( double ptMin,double minRapidity,double maxRapidity,
00042                                    double tip,double lip,int minHit, bool chargedOnly, 
00043                                    std::vector<int> pdgId = std::vector<int>()) :
00044     ptMin_( ptMin ), minRapidity_( minRapidity ), maxRapidity_( maxRapidity ),
00045     tip_( tip ), lip_( lip ), minHit_( minHit ), chargedOnly_(chargedOnly), pdgId_( pdgId ) { }
00046 
00047 
00048     CosmicTrackingParticleSelector  ( const edm::ParameterSet & cfg ) :
00049       ptMin_(cfg.getParameter<double>("ptMin")),
00050       minRapidity_(cfg.getParameter<double>("minRapidity")),
00051       maxRapidity_(cfg.getParameter<double>("maxRapidity")),
00052       tip_(cfg.getParameter<double>("tip")),
00053       lip_(cfg.getParameter<double>("lip")),
00054       minHit_(cfg.getParameter<int>("minHit")),
00055       chargedOnly_(cfg.getParameter<bool>("chargedOnly")),
00056       pdgId_(cfg.getParameter<std::vector<int> >("pdgId")) { }
00057       
00058       
00059       void select( const edm::Handle<collection>& c, const edm::Event & event, const edm::EventSetup& setup) {
00060         selected_.clear();
00061         edm::Handle<reco::BeamSpot> beamSpot;
00062         event.getByLabel(edm::InputTag("offlineBeamSpot"), beamSpot); 
00063         for( TrackingParticleCollection::const_iterator itp = c->begin(); 
00064              itp != c->end(); ++ itp )
00065           if ( operator()(*itp,beamSpot.product(),event,setup) ) {
00066             selected_.push_back( & * itp );
00067           }
00068       }
00069       
00070     const_iterator begin() const { return selected_.begin(); }
00071     const_iterator end() const { return selected_.end(); }
00072     
00073     // Operator() performs the selection: e.g. if (tPSelector(tp, bs, event, evtsetup)) {...
00074     bool operator()( const TrackingParticle & tp, const reco::BeamSpot* bs, const edm::Event &iEvent, const edm::EventSetup& iSetup ) const { 
00075       if (chargedOnly_ && tp.charge()==0) return false;//select only if charge!=0
00076       bool testId = false;
00077       unsigned int idSize = pdgId_.size();
00078       if (idSize==0) testId = true;
00079       else for (unsigned int it=0;it!=idSize;++it){
00080         if (tp.pdgId()==pdgId_[it]) testId = true;
00081       }
00082       
00083       edm::ESHandle<TrackerGeometry> tracker;
00084       iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
00085       
00086       edm::ESHandle<MagneticField> theMF;
00087       iSetup.get<IdealMagneticFieldRecord>().get(theMF);
00088       
00089       GlobalVector finalGV (0,0,0);
00090       GlobalPoint finalGP(0,0,0);      
00091       GlobalVector momentum(0,0,0);//At the PCA
00092       GlobalPoint vertex(0,0,0);//At the PCA
00093       double radius(9999);
00094       bool found(0);
00095       
00096       const std::vector<PSimHit> & simHits = tp.trackPSimHit(DetId::Tracker);
00097       for(std::vector<PSimHit>::const_iterator it=simHits.begin(); it!=simHits.end(); ++it){
00098         const GeomDet* tmpDet  = tracker->idToDet( DetId(it->detUnitId()) ) ;
00099         LocalVector  lv = it->momentumAtEntry();
00100         Local3DPoint lp = it->localPosition ();
00101         GlobalVector gv = tmpDet->surface().toGlobal( lv );
00102         GlobalPoint  gp = tmpDet->surface().toGlobal( lp );
00103         if(gp.perp()<radius){
00104           found=true;
00105           radius = gp.perp();
00106           finalGV = gv;
00107           finalGP = gp;
00108         }
00109       }
00110       if(!found) return 0;
00111       else
00112         {
00113           FreeTrajectoryState ftsAtProduction(finalGP,finalGV,TrackCharge(tp.charge()),theMF.product());
00114           TSCBLBuilderNoMaterial tscblBuilder;
00115           TrajectoryStateClosestToBeamLine tsAtClosestApproach = tscblBuilder(ftsAtProduction,*bs);//as in TrackProducerAlgorithm
00116           if(!tsAtClosestApproach.isValid()){
00117             std::cout << "WARNING: tsAtClosestApproach is not valid" << std::endl;
00118             return 0;
00119           }
00120           else
00121             {
00122               momentum = tsAtClosestApproach.trackStateAtPCA().momentum();
00123               vertex = tsAtClosestApproach.trackStateAtPCA().position();
00124               return (
00125                       tp.matchedHit() >= minHit_ &&
00126                       sqrt(momentum.perp2()) >= ptMin_ && 
00127                       momentum.eta() >= minRapidity_ && momentum.eta() <= maxRapidity_ && 
00128                       sqrt(vertex.perp2()) <= tip_ &&
00129                       fabs(vertex.z()) <= lip_ 
00130                       );
00131             }
00132         }
00133     }
00134     
00135     size_t size() const { return selected_.size(); }
00136     
00137  private:
00138     
00139     double ptMin_;
00140     double minRapidity_;
00141     double maxRapidity_;
00142     double tip_;
00143     double lip_;
00144     int    minHit_;
00145     bool chargedOnly_;
00146       std::vector<int> pdgId_;
00147       container selected_;
00148       
00149       
00150 };
00151 
00152 #endif