CMS 3D CMS Logo

ImpactParameterAlgorithm Class Reference

#include <RecoTauTag/ImpactParameter/interface/ImpactParameterAlgorithm.h>

List of all members.

Public Member Functions

 ImpactParameterAlgorithm ()
 ImpactParameterAlgorithm (const ParameterSet &parameters)
void setPrimaryVertex (Vertex *pv)
void setTransientTrackBuilder (const TransientTrackBuilder *)
pair< float,
TauImpactParameterInfo
tag (const IsolatedTauTagInfoRef &, const Vertex &)
 ~ImpactParameterAlgorithm ()

Private Attributes

double ip_max
double ip_min
VertexprimaryVertex
double sip_min
const TransientTrackBuildertransientTrackBuilder
bool use3D
bool use_sign


Detailed Description

Definition at line 23 of file ImpactParameterAlgorithm.h.


Constructor & Destructor Documentation

ImpactParameterAlgorithm::ImpactParameterAlgorithm ( const ParameterSet parameters  ) 

Definition at line 16 of file ImpactParameterAlgorithm.cc.

References edm::ParameterSet::getParameter(), ip_max, ip_min, sip_min, use3D, and use_sign.

00016                                                                                  {
00017         ip_min       = parameters.getParameter<double>("TauImpactParameterMin");
00018         ip_max       = parameters.getParameter<double>("TauImpactParameterMax");
00019         sip_min      = parameters.getParameter<double>("TauImpactParameterSignificanceMin");
00020         use_sign     = parameters.getParameter<bool>("UseTauImpactParameterSign");
00021         use3D        = parameters.getParameter<bool>("UseTau3DImpactParameter");
00022 
00023 }

ImpactParameterAlgorithm::ImpactParameterAlgorithm (  ) 

Definition at line 8 of file ImpactParameterAlgorithm.cc.

References ip_max, ip_min, sip_min, use3D, and use_sign.

00008                                                   {
00009         ip_min   = -9999;
00010         ip_max   = 9999;
00011         sip_min  = 0;
00012         use_sign = false;
00013         use3D    = false; 
00014 }

ImpactParameterAlgorithm::~ImpactParameterAlgorithm (  )  [inline]

Definition at line 34 of file ImpactParameterAlgorithm.h.

00034 {}


Member Function Documentation

void ImpactParameterAlgorithm::setPrimaryVertex ( Vertex pv  )  [inline]

Definition at line 36 of file ImpactParameterAlgorithm.h.

References primaryVertex.

00036 {primaryVertex = pv;}

void ImpactParameterAlgorithm::setTransientTrackBuilder ( const TransientTrackBuilder builder  ) 

Definition at line 25 of file ImpactParameterAlgorithm.cc.

References transientTrackBuilder.

Referenced by ImpactParameter::produce().

00025                                                                                              { 
00026         transientTrackBuilder = builder; 
00027 }

std::pair< float, TauImpactParameterInfo > ImpactParameterAlgorithm::tag ( const IsolatedTauTagInfoRef &  tauRef,
const Vertex pv 
)

Definition at line 30 of file ImpactParameterAlgorithm.cc.

References SignedImpactParameter3D::apply(), SignedTransverseImpactParameter::apply(), edm::RefVector< C, T, F >::begin(), TransientTrackBuilder::build(), GenMuonPlsPt100GeV_cfg::cout, direction, reco::TauImpactParameterInfo::discriminator(), edm::RefVector< C, T, F >::end(), lat::endl(), Measurement1D::error(), reco::TauImpactParameterTrackData::ip3D, ip_max, ip_min, metsig::jet, reco::Particle::px(), reco::Particle::py(), reco::Particle::pz(), reco::TauImpactParameterInfo::setIsolatedTauTag(), sip_min, reco::TauImpactParameterInfo::storeTrackData(), tracks, transientTrackBuilder, reco::TauImpactParameterTrackData::transverseIp, use3D, use_sign, and Measurement1D::value().

Referenced by ImpactParameter::produce().

00030                                                                                                                            {
00031 
00032         if(transientTrackBuilder == 0){
00033              cout << "Transient track builder is 0. abort!" << endl;
00034              abort(); //FIXME: trow an exception here
00035         }
00036 
00037         TauImpactParameterInfo resultExtended;
00038         resultExtended.setIsolatedTauTag(tauRef);
00039 
00040         const Jet* jet = tauRef->jet().get();
00041         GlobalVector direction(jet->px(),jet->py(),jet->pz());
00042 
00043         const TrackRefVector& tracks = tauRef->selectedTracks();
00044 
00045         RefVector<TrackCollection>::const_iterator iTrack;
00046         for(iTrack = tracks.begin(); iTrack!= tracks.end(); iTrack++){
00047 
00048           const TransientTrack transientTrack = (transientTrackBuilder->build(&(**iTrack)));
00049 
00050           SignedTransverseImpactParameter stip;
00051           Measurement1D ip = stip.apply(transientTrack,direction,pv).second;
00052 
00053           SignedImpactParameter3D signed_ip3D;
00054           Measurement1D ip3D = signed_ip3D.apply(transientTrack,direction,pv).second;
00055           //cout << "check pv,ip3d,track z " << pv.z() << " " << ip3D.value() << " " << transientTrack->dz() << endl;
00056           if(!use_sign){
00057             Measurement1D tmp2D(fabs(ip.value()),ip.error());
00058             ip = tmp2D;
00059 
00060             Measurement1D tmp3D(fabs(ip3D.value()),ip3D.error());
00061             ip3D = tmp3D;
00062           }
00063 
00064           reco::TauImpactParameterTrackData theData;
00065 
00066           theData.transverseIp = ip;
00067           theData.ip3D = ip3D;
00068           resultExtended.storeTrackData(*iTrack,theData);
00069 
00070         }
00071 
00072         float discriminator = resultExtended.discriminator(ip_min,ip_max,sip_min,use_sign,use3D);
00073 
00074         return std::make_pair( discriminator, resultExtended );
00075 }


Member Data Documentation

double ImpactParameterAlgorithm::ip_max [private]

Definition at line 47 of file ImpactParameterAlgorithm.h.

Referenced by ImpactParameterAlgorithm(), and tag().

double ImpactParameterAlgorithm::ip_min [private]

Definition at line 47 of file ImpactParameterAlgorithm.h.

Referenced by ImpactParameterAlgorithm(), and tag().

Vertex* ImpactParameterAlgorithm::primaryVertex [private]

Definition at line 44 of file ImpactParameterAlgorithm.h.

Referenced by setPrimaryVertex().

double ImpactParameterAlgorithm::sip_min [private]

Definition at line 47 of file ImpactParameterAlgorithm.h.

Referenced by ImpactParameterAlgorithm(), and tag().

const TransientTrackBuilder* ImpactParameterAlgorithm::transientTrackBuilder [private]

Definition at line 53 of file ImpactParameterAlgorithm.h.

Referenced by setTransientTrackBuilder(), and tag().

bool ImpactParameterAlgorithm::use3D [private]

Definition at line 50 of file ImpactParameterAlgorithm.h.

Referenced by ImpactParameterAlgorithm(), and tag().

bool ImpactParameterAlgorithm::use_sign [private]

Definition at line 50 of file ImpactParameterAlgorithm.h.

Referenced by ImpactParameterAlgorithm(), and tag().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:25:49 2009 for CMSSW by  doxygen 1.5.4