CMS 3D CMS Logo

Public Member Functions

CosmicParametersDefinerForTP Class Reference

#include <CosmicParametersDefinerForTP.h>

Inheritance diagram for CosmicParametersDefinerForTP:
ParametersDefinerForTP

List of all members.

Public Member Functions

 CosmicParametersDefinerForTP ()
virtual ParticleBase::Vector momentum (const edm::Event &iEvent, const edm::EventSetup &iSetup, const TrackingParticle &tp) const
virtual ParticleBase::Point vertex (const edm::Event &iEvent, const edm::EventSetup &iSetup, const TrackingParticle &tp) const

Detailed Description

Author:
Boris Mangano (UCSD) 5/7/2009

Definition at line 13 of file CosmicParametersDefinerForTP.h.


Constructor & Destructor Documentation

CosmicParametersDefinerForTP::CosmicParametersDefinerForTP ( ) [inline]

Definition at line 16 of file CosmicParametersDefinerForTP.h.

{};

Member Function Documentation

ParticleBase::Vector CosmicParametersDefinerForTP::momentum ( const edm::Event iEvent,
const edm::EventSetup iSetup,
const TrackingParticle tp 
) const [virtual]

Reimplemented from ParametersDefinerForTP.

Definition at line 22 of file CosmicParametersDefinerForTP.cc.

References ParticleBase::charge(), newFWLiteAna::found, edm::EventSetup::get(), edm::Event::getByLabel(), TrajectoryStateClosestToBeamLine::isValid(), FreeTrajectoryState::momentum(), L1TEmulatorMonitor_cff::p, PV3DBase< T, PVType, FrameType >::perp(), CosmicsPD_Skims::radius, dt_offlineAnalysis_common_cff::reco, trackerHits::simHits, GeomDet::surface(), Surface::toGlobal(), patCandidatesForDimuonsSequences_cff::tracker, DetId::Tracker, TrackingParticle::trackPSimHit(), TrajectoryStateClosestToBeamLine::trackStateAtPCA(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

                                                                                                                            {
  // to add a new implementation for cosmic. For the moment, it is just as for the base class:
  using namespace edm;
  using namespace std;
  using namespace reco;

  ESHandle<TrackerGeometry> tracker;
  iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
  
  edm::ESHandle<MagneticField> theMF;
  iSetup.get<IdealMagneticFieldRecord>().get(theMF);
  
  edm::Handle<reco::BeamSpot> bs;
  iEvent.getByLabel(InputTag("offlineBeamSpot"),bs);

  // cout<<"TrackingParticle pdgId = "<<tp.pdgId()<<endl;
  // cout<<"with tp.vertex(): ("<<tp.vertex().x()<<", "<<tp.vertex().y()<<", "<<tp.vertex().z()<<")"<<endl;
  // cout<<"with tp.momentum(): ("<<tp.momentum().x()<<", "<<tp.momentum().y()<<", "<<tp.momentum().z()<<")"<<endl;
  
  GlobalVector finalGV;
  GlobalPoint finalGP;
  double radius(9999);
  bool found(0);
  ParticleBase::Vector momentum(0,0,0);
  
  const vector<PSimHit> & simHits = tp.trackPSimHit(DetId::Tracker);
  for(vector<PSimHit>::const_iterator it=simHits.begin(); it!=simHits.end(); ++it){
    const GeomDet* tmpDet  = tracker->idToDet( DetId(it->detUnitId()) ) ;
    LocalVector  lv = it->momentumAtEntry();
    Local3DPoint lp = it->localPosition ();
    GlobalVector gv = tmpDet->surface().toGlobal( lv );
    GlobalPoint  gp = tmpDet->surface().toGlobal( lp );
    if(gp.perp()<radius){
      found=true;
      radius = gp.perp();
      finalGV = gv;
      finalGP = gp;
    }
  }

  //cout<<"found = "<<found<<endl;
  // cout<<"Closest Hit Position: ("<<finalGP.x()<<", "<<finalGP.y()<<", "<<finalGP.z()<<")"<<endl;
  //cout<<"Momentum at Closest Hit to BL: ("<<finalGV.x()<<", "<<finalGV.y()<<", "<<finalGV.z()<<")"<<endl;

  if(found) 
    {
      FreeTrajectoryState ftsAtProduction(finalGP,finalGV,TrackCharge(tp.charge()),theMF.product());
      TSCBLBuilderNoMaterial tscblBuilder;
      TrajectoryStateClosestToBeamLine tsAtClosestApproach = tscblBuilder(ftsAtProduction,*bs);//as in TrackProducerAlgorithm
      if(tsAtClosestApproach.isValid()){
        GlobalVector p = tsAtClosestApproach.trackStateAtPCA().momentum();
        momentum = ParticleBase::Vector(p.x(), p.y(), p.z());
      }
      return momentum;
    }
  return momentum;
}
ParticleBase::Point CosmicParametersDefinerForTP::vertex ( const edm::Event iEvent,
const edm::EventSetup iSetup,
const TrackingParticle tp 
) const [virtual]

Reimplemented from ParametersDefinerForTP.

Definition at line 80 of file CosmicParametersDefinerForTP.cc.

References ParticleBase::charge(), newFWLiteAna::found, edm::EventSetup::get(), edm::Event::getByLabel(), TrajectoryStateClosestToBeamLine::isValid(), PV3DBase< T, PVType, FrameType >::perp(), FreeTrajectoryState::position(), CosmicsPD_Skims::radius, dt_offlineAnalysis_common_cff::reco, trackerHits::simHits, GeomDet::surface(), Surface::toGlobal(), patCandidatesForDimuonsSequences_cff::tracker, DetId::Tracker, TrackingParticle::trackPSimHit(), TrajectoryStateClosestToBeamLine::trackStateAtPCA(), v, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

                                                                                                                                           {
  
  using namespace edm;
  using namespace std;
  using namespace reco;

  ESHandle<TrackerGeometry> tracker;
  iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
  
  edm::ESHandle<MagneticField> theMF;
  iSetup.get<IdealMagneticFieldRecord>().get(theMF);
  
  edm::Handle<reco::BeamSpot> bs;
  iEvent.getByLabel(InputTag("offlineBeamSpot"),bs);

  GlobalVector finalGV;
  GlobalPoint finalGP;
  double radius(9999);
  bool found(0);
  ParticleBase::Point vertex(0,0,0);

  const vector<PSimHit> & simHits = tp.trackPSimHit(DetId::Tracker);
  for(vector<PSimHit>::const_iterator it=simHits.begin(); it!=simHits.end(); ++it){
    const GeomDet* tmpDet  = tracker->idToDet( DetId(it->detUnitId()) ) ;
    LocalVector  lv = it->momentumAtEntry();
    Local3DPoint lp = it->localPosition ();
    GlobalVector gv = tmpDet->surface().toGlobal( lv );
    GlobalPoint  gp = tmpDet->surface().toGlobal( lp );
    if(gp.perp()<radius){
      found=true;
      radius = gp.perp();
      finalGV = gv;
      finalGP = gp;
    }
  }
  if(found)
    {
      FreeTrajectoryState ftsAtProduction(finalGP,finalGV,TrackCharge(tp.charge()),theMF.product());
      TSCBLBuilderNoMaterial tscblBuilder;
      TrajectoryStateClosestToBeamLine tsAtClosestApproach = tscblBuilder(ftsAtProduction,*bs);//as in TrackProducerAlgorithm
      if(tsAtClosestApproach.isValid()){
        GlobalPoint v = tsAtClosestApproach.trackStateAtPCA().position();
        vertex = ParticleBase::Point(v.x()-bs->x0(),v.y()-bs->y0(),v.z()-bs->z0());
      }
      return vertex;
    }
  return vertex;
}