CMS 3D CMS Logo

Public Member Functions

ParametersDefinerForTP Class Reference

#include <ParametersDefinerForTP.h>

Inheritance diagram for ParametersDefinerForTP:
CosmicParametersDefinerForTP

List of all members.

Public Member Functions

virtual ParticleBase::Vector momentum (const edm::Event &iEvent, const edm::EventSetup &iSetup, const TrackingParticle &tp) const
 ParametersDefinerForTP ()
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 15 of file ParametersDefinerForTP.h.


Constructor & Destructor Documentation

ParametersDefinerForTP::ParametersDefinerForTP ( ) [inline]

Definition at line 18 of file ParametersDefinerForTP.h.

{};

Member Function Documentation

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

Reimplemented in CosmicParametersDefinerForTP.

Definition at line 15 of file ParametersDefinerForTP.cc.

References ParticleBase::charge(), edm::EventSetup::get(), edm::Event::getByLabel(), ParticleBase::momentum(), AlCaHLTBitMon_ParallelJobs::p, ParticleBase::vertex(), 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;

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

  ParticleBase::Vector momentum(0, 0, 0); 

  FreeTrajectoryState ftsAtProduction(GlobalPoint(tp.vertex().x(),tp.vertex().y(),tp.vertex().z()),
                                      GlobalVector(tp.momentum().x(),tp.momentum().y(),tp.momentum().z()),
                                      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;
}
ParticleBase::Point ParametersDefinerForTP::vertex ( const edm::Event iEvent,
const edm::EventSetup iSetup,
const TrackingParticle tp 
) const [virtual]

Reimplemented in CosmicParametersDefinerForTP.

Definition at line 42 of file ParametersDefinerForTP.cc.

References ParticleBase::charge(), edm::EventSetup::get(), edm::Event::getByLabel(), ParticleBase::momentum(), v, ParticleBase::vertex(), 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;

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

  ParticleBase::Point vertex(0, 0, 0);
  
  FreeTrajectoryState ftsAtProduction(GlobalPoint(tp.vertex().x(),tp.vertex().y(),tp.vertex().z()),
                                      GlobalVector(tp.momentum().x(),tp.momentum().y(),tp.momentum().z()),
                                      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;
}