CMS 3D CMS Logo

Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Attributes

reco::Conversion Class Reference

#include <DataFormats/EgammaCandidates/interface/Conversion.h>

List of all members.

Public Types

enum  ConversionAlgorithm {
  undefined = 0, ecalSeeded = 1, trackerOnly = 2, mixed = 3,
  pflow = 4, algoSize = 5
}
enum  ConversionQuality {
  generalTracksOnly = 0, arbitratedEcalSeeded = 1, arbitratedMerged = 2, arbitratedMergedEcalGeneral = 3,
  highPurity = 8, highEfficiency = 9, ecalMatched1Track = 10, ecalMatched2Track = 11
}

Public Member Functions

ConversionAlgorithm algo () const
std::string algoName () const
const std::vector
< reco::CaloClusterPtr > & 
bcMatchingWithTracks () const
reco::CaloClusterPtrVector caloCluster () const
 Pointer to CaloCluster (foe Egamma Conversions it points to a SuperCluster)
Conversionclone () const
 returns a clone of the candidate
 Conversion (const reco::CaloClusterPtrVector clu, const std::vector< edm::RefToBase< reco::Track > > tr, const std::vector< math::XYZPointF > trackPositionAtEcal, const reco::Vertex &convVtx, const std::vector< reco::CaloClusterPtr > &matchingBC, const float DCA, const std::vector< math::XYZPointF > &innPoint, const std::vector< math::XYZVectorF > &trackPin, const std::vector< math::XYZVectorF > &trackPout, const std::vector< uint8_t > nHitsBeforeVtx, const std::vector< Measurement1DFloat > &dlClosestHitToVtx, uint8_t nSharedHits, const float mva, ConversionAlgorithm=undefined)
 Conversion (const reco::CaloClusterPtrVector clu, const std::vector< reco::TrackRef > tr, const reco::Vertex &convVtx, ConversionAlgorithm=undefined)
 Conversion (const reco::CaloClusterPtrVector clu, const std::vector< edm::RefToBase< reco::Track > > tr, const reco::Vertex &convVtx, ConversionAlgorithm=undefined)
 Conversion ()
 Conversion (const reco::CaloClusterPtrVector clu, const std::vector< reco::TrackRef > tr, const std::vector< math::XYZPointF > trackPositionAtEcal, const reco::Vertex &convVtx, const std::vector< reco::CaloClusterPtr > &matchingBC, const float DCA, const std::vector< math::XYZPointF > &innPoint, const std::vector< math::XYZVectorF > &trackPin, const std::vector< math::XYZVectorF > &trackPout, const float mva, ConversionAlgorithm=undefined)
const reco::VertexconversionVertex () const
 returns the reco conversion vertex
double dEtaTracksAtEcal () const
double distOfMinimumApproach () const
const std::vector
< Measurement1DFloat > & 
dlClosestHitToVtx () const
 Vector of signed decay length with uncertainty from nearest hit on track to the conversion vtx positions.
double dPhiTracksAtEcal () const
double dPhiTracksAtVtx () const
double dxy (const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
double dz (const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
const std::vector
< math::XYZPointF > & 
ecalImpactPosition () const
double EoverP () const
double EoverPrefittedTracks () const
bool isConverted () const
 Bool flagging objects having track size >0.
double lxy (const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
double lz (const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
double MVAout () const
 get the value of the TMVA output
const std::vector< uint8_t > & nHitsBeforeVtx () const
 Vector of the number of hits before the vertex along each track trajector.
uint8_t nSharedHits () const
 number of shared hits btw the two track
unsigned int nTracks () const
 Number of tracks= 0,1,2.
std::vector< float > const oneLegMVA ()
 get the MVS output from PF for one leg conversions
double pairCotThetaSeparation () const
 Delta cot(Theta) where Theta is the angle in the (y,z) plane between the two tracks. Original tracks are used.
double pairInvariantMass () const
 if nTracks=2 returns the pair invariant mass. Original tracks are used here
math::XYZVectorF pairMomentum () const
 Conversion tracks momentum from the tracks inner momentum.
bool quality (ConversionQuality q) const
math::XYZTLorentzVectorF refittedPair4Momentum () const
 Conversion track pair 4-momentum from the tracks refitted with vertex constraint.
math::XYZVectorF refittedPairMomentum () const
 Conversion tracks momentum from the tracks refitted with vertex constraint.
void setConversionAlgorithm (const ConversionAlgorithm a, bool set=true)
 Conversion Track algorithm/provenance.
void setMatchingSuperCluster (const reco::CaloClusterPtrVector &sc)
void setMVAout (const float &mva)
 set the value of the TMVA output
void setOneLegMVA (const std::vector< float > &mva)
 set the MVS output from PF for one leg conversions
void setQuality (ConversionQuality q, bool b)
std::vector< edm::RefToBase
< reco::Track > > 
tracks () const
 vector of track to base references
const std::vector
< math::XYZPointF > & 
tracksInnerPosition () const
 Vector containing the position of the innermost hit of each track.
const std::vector
< math::XYZVectorF > & 
tracksPin () const
 Vector of track momentum measured at the innermost hit.
const std::vector
< math::XYZVectorF > & 
tracksPout () const
 Vector of track momentum measured at the outermost hit.
std::vector< double > tracksSigned_d0 () const
 signed transverse impact parameter for each track
double zOfPrimaryVertexFromTracks (const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
virtual ~Conversion ()
 destructor

Static Public Member Functions

static ConversionAlgorithm algoByName (const std::string &name)
static std::string algoName (ConversionAlgorithm)

Static Public Attributes

static const std::string algoNames [] = { "undefined","ecalSeeded","trackerOnly","mixed","pflow"}

Private Attributes

uint8_t algorithm_
 conversion algorithm/provenance
reco::CaloClusterPtrVector caloCluster_
 vector pointer to a/multiple seed CaloCluster(s)
std::vector< Measurement1DFloatdlClosestHitToVtx_
 signed decay length and uncertainty from nearest hit on track to conversion vertex
std::vector< uint8_t > nHitsBeforeVtx_
 number of hits before the vertex on each trackerOnly
uint8_t nSharedHits_
 number of shared hits between tracks
uint16_t qualityMask_
reco::Vertex theConversionVertex_
 Fitted Kalman conversion vertex.
std::vector< reco::CaloClusterPtrtheMatchingBCs_
 Clusters mathing the tracks (these are not the seeds)
float theMinDistOfApproach_
 Distance of min approach of the two tracks.
float theMVAout_
 TMVA output.
std::vector< float > theOneLegMVA_
 vectors of TMVA outputs from pflow for one leg conversions
std::vector< math::XYZPointFthePositionAtEcal_
 position at the ECAl surface of the track extrapolation
std::vector< math::XYZPointFtheTrackInnerPosition_
 P_in of tracks.
std::vector< math::XYZVectorFtheTrackPin_
 P_in of tracks.
std::vector< math::XYZVectorFtheTrackPout_
 P_out of tracks.
std::vector< reco::TrackReftracks_
 vector of Track references
std::vector< edm::RefToBase
< reco::Track > > 
trackToBaseRefs_
 vector Track RefToBase

Detailed Description

Author:
N.Marinelli University of Notre Dame, US
Version:
Id:
Conversion.h,v 1.24 2011/11/09 11:26:54 nancy Exp

Definition at line 25 of file Conversion.h.


Member Enumeration Documentation

Enumerator:
undefined 
ecalSeeded 
trackerOnly 
mixed 
pflow 
algoSize 

Definition at line 28 of file Conversion.h.

Enumerator:
generalTracksOnly 
arbitratedEcalSeeded 
arbitratedMerged 
arbitratedMergedEcalGeneral 
highPurity 
highEfficiency 
ecalMatched1Track 
ecalMatched2Track 

Definition at line 35 of file Conversion.h.


Constructor & Destructor Documentation

Conversion::Conversion ( )
Conversion::Conversion ( const reco::CaloClusterPtrVector  clu,
const std::vector< edm::RefToBase< reco::Track > >  tr,
const std::vector< math::XYZPointF trackPositionAtEcal,
const reco::Vertex convVtx,
const std::vector< reco::CaloClusterPtr > &  matchingBC,
const float  DCA,
const std::vector< math::XYZPointF > &  innPoint,
const std::vector< math::XYZVectorF > &  trackPin,
const std::vector< math::XYZVectorF > &  trackPout,
const std::vector< uint8_t >  nHitsBeforeVtx,
const std::vector< Measurement1DFloat > &  dlClosestHitToVtx,
uint8_t  nSharedHits,
const float  mva,
ConversionAlgorithm  algo = undefined 
)
Conversion::Conversion ( const reco::CaloClusterPtrVector  clu,
const std::vector< reco::TrackRef tr,
const std::vector< math::XYZPointF trackPositionAtEcal,
const reco::Vertex convVtx,
const std::vector< reco::CaloClusterPtr > &  matchingBC,
const float  DCA,
const std::vector< math::XYZPointF > &  innPoint,
const std::vector< math::XYZVectorF > &  trackPin,
const std::vector< math::XYZVectorF > &  trackPout,
const float  mva,
ConversionAlgorithm  algo = undefined 
)

Definition at line 9 of file Conversion.cc.

                                                  :  
                         

  caloCluster_(sc), tracks_(tr), 
  thePositionAtEcal_(trackPositionAtEcal), 
  theConversionVertex_(convVtx), 
  theMatchingBCs_(matchingBC), 
  theMinDistOfApproach_(DCA),
  theTrackInnerPosition_(innPoint),
  theTrackPin_(trackPin),
  theTrackPout_(trackPout),
  nSharedHits_(0),  
  theMVAout_(mva),
  algorithm_(algo),
  qualityMask_(0)
 { 
   
 }
Conversion::Conversion ( const reco::CaloClusterPtrVector  clu,
const std::vector< reco::TrackRef tr,
const reco::Vertex convVtx,
ConversionAlgorithm  algo = undefined 
)
Conversion::Conversion ( const reco::CaloClusterPtrVector  clu,
const std::vector< edm::RefToBase< reco::Track > >  tr,
const reco::Vertex convVtx,
ConversionAlgorithm  algo = undefined 
)
Conversion::~Conversion ( ) [virtual]

destructor

Definition at line 152 of file Conversion.cc.

{ }

Member Function Documentation

Conversion::ConversionAlgorithm reco::Conversion::algo ( ) const [inline]

Definition at line 229 of file Conversion.h.

References algorithm_.

Conversion::ConversionAlgorithm Conversion::algoByName ( const std::string &  name) [static]
std::string reco::Conversion::algoName ( ) const [inline]

Definition at line 234 of file Conversion.h.

References algorithm_, ecalSeeded, mixed, pflow, trackerOnly, and undefined.

                                               {
            
      switch(algorithm_)
        {
        case undefined: return "undefined";
        case ecalSeeded: return "ecalSeeded";
        case trackerOnly: return "trackerOnly";
        case mixed: return "mixed";
        case pflow: return "pflow";

        }
      return "undefined";
    }
std::string reco::Conversion::algoName ( ConversionAlgorithm  a) [inline, static]

Definition at line 248 of file Conversion.h.

References algoNames, and algoSize.

                                                              {
      if(int(a) < int(algoSize) && int(a)>0) return algoNames[int(a)];
      return "undefined";
    }
const std::vector<reco::CaloClusterPtr>& reco::Conversion::bcMatchingWithTracks ( ) const [inline]

Definition at line 155 of file Conversion.h.

References theMatchingBCs_.

Referenced by dEtaTracksAtEcal(), and dPhiTracksAtEcal().

{ return theMatchingBCs_;}
reco::CaloClusterPtrVector reco::Conversion::caloCluster ( ) const [inline]

Pointer to CaloCluster (foe Egamma Conversions it points to a SuperCluster)

Definition at line 97 of file Conversion.h.

References caloCluster_.

Referenced by EoverP(), and EoverPrefittedTracks().

{return caloCluster_ ;}
Conversion * Conversion::clone ( void  ) const

returns a clone of the candidate

Definition at line 165 of file Conversion.cc.

References Conversion().

Referenced by ConvertedPhotonProducer::cleanCollections().

                                     { 
  return new Conversion( * this ); 
}
const reco::Vertex& reco::Conversion::conversionVertex ( ) const [inline]
double Conversion::dEtaTracksAtEcal ( ) const

Definition at line 356 of file Conversion.cc.

References bcMatchingWithTracks(), ecalImpactPosition(), nTracks(), and query::result.

                                            {
  double result=-99.;


  if ( nTracks()==2 && bcMatchingWithTracks()[0].isNonnull() && bcMatchingWithTracks()[1].isNonnull() ) {

    result =ecalImpactPosition()[0].eta() - ecalImpactPosition()[1].eta();

  }



  return result;


}
double reco::Conversion::distOfMinimumApproach ( ) const [inline]
const std::vector<Measurement1DFloat>& reco::Conversion::dlClosestHitToVtx ( ) const [inline]

Vector of signed decay length with uncertainty from nearest hit on track to the conversion vtx positions.

Definition at line 167 of file Conversion.h.

References dlClosestHitToVtx_.

Referenced by TkConvValidator::analyze().

{ return dlClosestHitToVtx_; }
double Conversion::dPhiTracksAtEcal ( ) const

Definition at line 331 of file Conversion.cc.

References bcMatchingWithTracks(), ecalImpactPosition(), nTracks(), pi, and query::result.

                                            {
  double result =-99.;
  
  if (  nTracks()==2  && bcMatchingWithTracks()[0].isNonnull() && bcMatchingWithTracks()[1].isNonnull() ) {
    
    float recoPhi1 = ecalImpactPosition()[0].phi();
    if( recoPhi1   > pi)  { recoPhi1 = recoPhi1 - twopi;}
    if( recoPhi1  < -pi)  { recoPhi1 = recoPhi1 + twopi;}

    float recoPhi2 = ecalImpactPosition()[1].phi();
    if( recoPhi2   > pi)  { recoPhi2 = recoPhi2 - twopi;}
    if( recoPhi2  < -pi)  { recoPhi2 = recoPhi2 + twopi;}

    result = recoPhi1 -recoPhi2;

    if( result   > pi)  { result = result - twopi;}
    if( result  < -pi)  { result = result + twopi;}

  }

  return result;


}
double Conversion::dPhiTracksAtVtx ( ) const

Definition at line 318 of file Conversion.cc.

References nTracks(), pi, query::result, and tracksPin().

Referenced by TkConvValidator::analyze(), and FWConvTrackHitsDetailView::setTextInfo().

                                           {
  double result=-99;
  if  ( nTracks()==2 ) {
    result = tracksPin()[0].phi() - tracksPin()[1].phi();
    if( result   > pi)  { result = result - twopi;}
    if( result  < -pi)  { result = result + twopi;}
  }

  return result;


}
double Conversion::dxy ( const math::XYZPoint myBeamSpot = math::XYZPoint()) const

Definition at line 373 of file Conversion.cc.

References conversionVertex(), reco::Vertex::isValid(), refittedPairMomentum(), reco::Vertex::x(), and reco::Vertex::y().

                                                           {

  const reco::Vertex &vtx = conversionVertex();
  if (!vtx.isValid()) return -9999.;

  math::XYZVectorF mom = refittedPairMomentum();
  
  double dxy = (-(vtx.x() - myBeamSpot.x())*mom.y() + (vtx.y() - myBeamSpot.y())*mom.x())/mom.rho();
  return dxy;  
  
}
double Conversion::dz ( const math::XYZPoint myBeamSpot = math::XYZPoint()) const

Definition at line 385 of file Conversion.cc.

References conversionVertex(), reco::Vertex::isValid(), refittedPairMomentum(), reco::Vertex::x(), reco::Vertex::y(), and reco::Vertex::z().

Referenced by zOfPrimaryVertexFromTracks().

                                                          {

  const reco::Vertex &vtx = conversionVertex();
  if (!vtx.isValid()) return -9999.;

  math::XYZVectorF mom = refittedPairMomentum();
  
  double dz = (vtx.z()-myBeamSpot.z()) - ((vtx.x()-myBeamSpot.x())*mom.x()+(vtx.y()-myBeamSpot.y())*mom.y())/mom.rho() * mom.z()/mom.rho();
  return dz;  
  
}
const std::vector<math::XYZPointF>& reco::Conversion::ecalImpactPosition ( ) const [inline]

The following are variables provided per each track positions of the track extrapolation at the ECAL front face

Definition at line 153 of file Conversion.h.

References thePositionAtEcal_.

Referenced by dEtaTracksAtEcal(), and dPhiTracksAtEcal().

double Conversion::EoverP ( ) const

Super Cluster energy divided by track pair momentum if Standard seeding method. If a pointer to two (or more clusters) is stored in the conversion, this method returns the energy sum of clusters divided by the track pair momentum Track innermost momentum is used here

Definition at line 263 of file Conversion.cc.

References caloCluster(), i, nTracks(), pairMomentum(), findQualityFiles::size, edm::PtrVectorBase::size(), and mathSSE::sqrt().

Referenced by ConversionLikelihoodCalculator::calculateLikelihood().

                                  {


  double ep=-99.;

  if ( nTracks() > 0  ) {
    unsigned int size= this->caloCluster().size();
    float etot=0.;
    for ( unsigned int i=0; i<size; i++) {
      etot+= caloCluster()[i]->energy();
    }
    if (this->pairMomentum().Mag2() !=0) ep= etot/sqrt(this->pairMomentum().Mag2());
  }



  return ep;  

}
double Conversion::EoverPrefittedTracks ( ) const

Super Cluster energy divided by track pair momentum if Standard seeing method. If a pointer to two (or more clusters) is stored in the conversion, this method returns the energy sum of clusters divided by the track pair momentum Track momentum refitted with vertex constraint is used

Definition at line 285 of file Conversion.cc.

References caloCluster(), i, nTracks(), refittedPairMomentum(), findQualityFiles::size, edm::PtrVectorBase::size(), and mathSSE::sqrt().

                                                {


  double ep=-99.;
 
  if ( nTracks() > 0  ) {
    unsigned int size= this->caloCluster().size();
    float etot=0.;
    for ( unsigned int i=0; i<size; i++) {
      etot+= caloCluster()[i]->energy();
    }
    if (this->refittedPairMomentum().Mag2() !=0) ep= etot/sqrt(this->refittedPairMomentum().Mag2());
  }



  return ep;  

}
bool Conversion::isConverted ( ) const

Bool flagging objects having track size >0.

Definition at line 186 of file Conversion.cc.

References nTracks().

                                   {
  
  if ( this->nTracks() == 2 ) 
    return true;
  else
    return false;
}
double Conversion::lxy ( const math::XYZPoint myBeamSpot = math::XYZPoint()) const

Definition at line 397 of file Conversion.cc.

References conversionVertex(), reco::Vertex::isValid(), refittedPairMomentum(), reco::Vertex::x(), and reco::Vertex::y().

                                                           {

  const reco::Vertex &vtx = conversionVertex();
  if (!vtx.isValid()) return -9999.;

  math::XYZVectorF mom = refittedPairMomentum();
  
  double dbsx = vtx.x() - myBeamSpot.x();
  double dbsy = vtx.y() - myBeamSpot.y();
  double lxy = (mom.x()*dbsx + mom.y()*dbsy)/mom.rho();
  return lxy;  
  
}
double Conversion::lz ( const math::XYZPoint myBeamSpot = math::XYZPoint()) const

Definition at line 411 of file Conversion.cc.

References abs, conversionVertex(), reco::Vertex::isValid(), refittedPairMomentum(), and reco::Vertex::z().

                                                          {

  const reco::Vertex &vtx = conversionVertex();
  if (!vtx.isValid()) return -9999.;

  math::XYZVectorF mom = refittedPairMomentum();
  
  double lz = (vtx.z() - myBeamSpot.z())*mom.z()/std::abs(mom.z());
  return lz;  
  
}
double reco::Conversion::MVAout ( ) const [inline]

get the value of the TMVA output

Definition at line 107 of file Conversion.h.

References theMVAout_.

{ return theMVAout_;}
const std::vector<uint8_t>& reco::Conversion::nHitsBeforeVtx ( ) const [inline]

Vector of the number of hits before the vertex along each track trajector.

Definition at line 165 of file Conversion.h.

References nHitsBeforeVtx_.

Referenced by TkConvValidator::analyze(), and ConversionTools::isGoodConversion().

{ return nHitsBeforeVtx_; }
uint8_t reco::Conversion::nSharedHits ( ) const [inline]

number of shared hits btw the two track

Definition at line 169 of file Conversion.h.

References nSharedHits_.

Referenced by TkConvValidator::analyze().

{ return nSharedHits_; }
unsigned int reco::Conversion::nTracks ( ) const [inline]
std::vector<float> const reco::Conversion::oneLegMVA ( ) [inline]

get the MVS output from PF for one leg conversions

Definition at line 109 of file Conversion.h.

References theOneLegMVA_.

{return theOneLegMVA_;}
double Conversion::pairCotThetaSeparation ( ) const

Delta cot(Theta) where Theta is the angle in the (y,z) plane between the two tracks. Original tracks are used.

Definition at line 213 of file Conversion.cc.

References nTracks(), funct::tan(), and tracksPin().

Referenced by TkConvValidator::analyze(), ConversionLikelihoodCalculator::calculateLikelihood(), and FWConvTrackHitsDetailView::setTextInfo().

                                                  {
  double dCotTheta=-99.;
  
  if ( nTracks()==2 ) {
    double theta1=this->tracksPin()[0].Theta();
    double theta2=this->tracksPin()[1].Theta();
    dCotTheta =  1./tan(theta1) - 1./tan(theta2) ;
  }

  return dCotTheta;

}
double Conversion::pairInvariantMass ( ) const

if nTracks=2 returns the pair invariant mass. Original tracks are used here

Definition at line 194 of file Conversion.cc.

References alignCSCRings::e, nTracks(), mathSSE::sqrt(), and tracksPin().

Referenced by TkConvValidator::analyze().

                                          {
  double invMass=-99.;
  const float mElec= 0.000511;
  if ( nTracks()==2 ) {
    double px= tracksPin()[0].x() +  tracksPin()[1].x();
    double py= tracksPin()[0].y() +  tracksPin()[1].y();
    double pz= tracksPin()[0].z() +  tracksPin()[1].z();
    double mom1= tracksPin()[0].Mag2();
    double mom2= tracksPin()[1].Mag2();
    double e = sqrt( mom1+ mElec*mElec ) + sqrt( mom2 + mElec*mElec );
    invMass= ( e*e - px*px -py*py - pz*pz);
    if ( invMass>0) invMass = sqrt(invMass);
    else 
      invMass = -1;
  }
  
  return invMass;
}
math::XYZVectorF Conversion::pairMomentum ( ) const

Conversion tracks momentum from the tracks inner momentum.

Definition at line 227 of file Conversion.cc.

References nTracks(), and tracksPin().

Referenced by EoverP(), FWConversionProxyBuilder::requestCommon(), and FWConvTrackHitsDetailView::setTextInfo().

                                                {
  
  if ( nTracks()==2 ) {
    return this->tracksPin()[0] +  this->tracksPin()[1];
  }
  return math::XYZVectorF(0.,0.,0.);



}
bool reco::Conversion::quality ( ConversionQuality  q) const [inline]

Definition at line 185 of file Conversion.h.

References lumiQueryAPI::q, and qualityMask_.

Referenced by TkConvValidator::analyze().

{ return  (qualityMask_ & (1<<q))>>q; }
math::XYZTLorentzVectorF Conversion::refittedPair4Momentum ( ) const

Conversion track pair 4-momentum from the tracks refitted with vertex constraint.

Definition at line 239 of file Conversion.cc.

References conversionVertex(), reco::Vertex::p4(), and p4.

Referenced by PF_PU_AssoMapAlgos::FindConversionVertex(), and refittedPairMomentum().

                                                                {

  math::XYZTLorentzVectorF p4;
  if ( this->conversionVertex().isValid() ) 
    p4 = this->conversionVertex().p4( 0.000511, 0.5);

  return p4;


}
math::XYZVectorF Conversion::refittedPairMomentum ( ) const

Conversion tracks momentum from the tracks refitted with vertex constraint.

Definition at line 252 of file Conversion.cc.

References conversionVertex(), and refittedPair4Momentum().

Referenced by TkConvValidator::analyze(), dxy(), dz(), EoverPrefittedTracks(), ConversionTools::isGoodConversion(), lxy(), lz(), ConversionTools::matchesConversion(), and ConversionProducer::matchingSC().

                                                        {

  if (  this->conversionVertex().isValid() ) {
    return this->refittedPair4Momentum().Vect();
  }
  return math::XYZVectorF(0.,0.,0.);

}
void reco::Conversion::setConversionAlgorithm ( const ConversionAlgorithm  a,
bool  set = true 
) [inline]

Conversion Track algorithm/provenance.

Definition at line 179 of file Conversion.h.

References a, algorithm_, and undefined.

{ if (set) algorithm_=a; else algorithm_=undefined;}
void reco::Conversion::setMatchingSuperCluster ( const reco::CaloClusterPtrVector sc) [inline]

Definition at line 177 of file Conversion.h.

References caloCluster_.

Referenced by ConversionProducer::buildCollection().

{ caloCluster_= sc;}
void reco::Conversion::setMVAout ( const float &  mva) [inline]

set the value of the TMVA output

Definition at line 173 of file Conversion.h.

References theMVAout_.

Referenced by ConvertedPhotonProducer::buildCollections().

{ theMVAout_=mva;}
void reco::Conversion::setOneLegMVA ( const std::vector< float > &  mva) [inline]

set the MVS output from PF for one leg conversions

Definition at line 175 of file Conversion.h.

References theOneLegMVA_.

Referenced by PFPhotonTranslator::createOneLegConversions().

{ theOneLegMVA_=mva;}
void reco::Conversion::setQuality ( ConversionQuality  q,
bool  b 
) [inline]

Definition at line 253 of file Conversion.h.

References lumiQueryAPI::q, and qualityMask_.

Referenced by ConversionProducer::buildCollection().

                                                                 {
      if (b)//regular OR if setting value to true
        qualityMask_ |= (1<<q) ;
      else // doing "half-XOR" if unsetting value
        qualityMask_ &= (~(1<<q));

    }
std::vector< edm::RefToBase< reco::Track > > Conversion::tracks ( void  ) const

vector of track to base references

Definition at line 170 of file Conversion.cc.

References tracks_, trackToBaseRefs_, and groupFilesInBlocks::tt.

Referenced by TkConvValidator::analyze(), FWConvTrackHitsDetailView::build(), ConversionLikelihoodCalculator::calculateLikelihood(), ConversionTools::matchesConversion(), nTracks(), FWConvTrackHitsDetailView::setTextInfo(), and tracksSigned_d0().

                                                              { 
  if (trackToBaseRefs_.size() ==0 ) {
 
    for (std::vector<reco::TrackRef>::const_iterator ref=tracks_.begin(); ref!=tracks_.end(); ref++ ) 
      {
        edm::RefToBase<reco::Track> tt(*ref);
        trackToBaseRefs_.push_back(tt);
        
      }  
  }

  return trackToBaseRefs_;
}
const std::vector<math::XYZPointF>& reco::Conversion::tracksInnerPosition ( ) const [inline]

Vector containing the position of the innermost hit of each track.

Definition at line 159 of file Conversion.h.

References theTrackInnerPosition_.

const std::vector<math::XYZVectorF>& reco::Conversion::tracksPin ( ) const [inline]

Vector of track momentum measured at the innermost hit.

Definition at line 163 of file Conversion.h.

References theTrackPin_.

Referenced by ConversionLikelihoodCalculator::calculateLikelihood(), dPhiTracksAtVtx(), pairCotThetaSeparation(), pairInvariantMass(), and pairMomentum().

{return theTrackPin_;}
const std::vector<math::XYZVectorF>& reco::Conversion::tracksPout ( ) const [inline]

Vector of track momentum measured at the outermost hit.

Definition at line 161 of file Conversion.h.

References theTrackPout_.

{return theTrackPout_;}
std::vector< double > Conversion::tracksSigned_d0 ( ) const

signed transverse impact parameter for each track

Definition at line 307 of file Conversion.cc.

References DeDxDiscriminatorTools::charge(), i, nTracks(), query::result, and tracks().

                                                      {
  std::vector<double> result;

  for (unsigned int i=0; i< nTracks(); i++)   
    result.push_back(tracks()[i]->d0()* tracks()[i]->charge()) ;

  return result;


}
double reco::Conversion::zOfPrimaryVertexFromTracks ( const math::XYZPoint myBeamSpot = math::XYZPoint()) const [inline]

Definition at line 149 of file Conversion.h.

References dz().

Referenced by TkConvValidator::analyze(), ConversionProducer::matchingSC(), and FWConversionProxyBuilder::requestCommon().

{ return dz(myBeamSpot) + myBeamSpot.z(); }

Member Data Documentation

std::string const Conversion::algoNames = { "undefined","ecalSeeded","trackerOnly","mixed","pflow"} [static]

Definition at line 44 of file Conversion.h.

Referenced by algoByName(), and algoName().

uint8_t reco::Conversion::algorithm_ [private]

conversion algorithm/provenance

Definition at line 223 of file Conversion.h.

Referenced by algo(), algoName(), Conversion(), and setConversionAlgorithm().

vector pointer to a/multiple seed CaloCluster(s)

Definition at line 193 of file Conversion.h.

Referenced by caloCluster(), and setMatchingSuperCluster().

signed decay length and uncertainty from nearest hit on track to conversion vertex

Definition at line 215 of file Conversion.h.

Referenced by dlClosestHitToVtx().

std::vector<uint8_t> reco::Conversion::nHitsBeforeVtx_ [private]

number of hits before the vertex on each trackerOnly

Definition at line 213 of file Conversion.h.

Referenced by nHitsBeforeVtx().

uint8_t reco::Conversion::nSharedHits_ [private]

number of shared hits between tracks

Definition at line 217 of file Conversion.h.

Referenced by Conversion(), and nSharedHits().

uint16_t reco::Conversion::qualityMask_ [private]

Definition at line 224 of file Conversion.h.

Referenced by Conversion(), quality(), and setQuality().

Fitted Kalman conversion vertex.

Definition at line 201 of file Conversion.h.

Referenced by conversionVertex().

Clusters mathing the tracks (these are not the seeds)

Definition at line 203 of file Conversion.h.

Referenced by bcMatchingWithTracks().

Distance of min approach of the two tracks.

Definition at line 205 of file Conversion.h.

Referenced by Conversion(), and distOfMinimumApproach().

TMVA output.

Definition at line 219 of file Conversion.h.

Referenced by Conversion(), MVAout(), and setMVAout().

std::vector<float> reco::Conversion::theOneLegMVA_ [private]

vectors of TMVA outputs from pflow for one leg conversions

Definition at line 221 of file Conversion.h.

Referenced by oneLegMVA(), and setOneLegMVA().

position at the ECAl surface of the track extrapolation

Definition at line 199 of file Conversion.h.

Referenced by Conversion(), and ecalImpactPosition().

P_in of tracks.

Definition at line 207 of file Conversion.h.

Referenced by Conversion(), and tracksInnerPosition().

P_in of tracks.

Definition at line 209 of file Conversion.h.

Referenced by Conversion(), and tracksPin().

P_out of tracks.

Definition at line 211 of file Conversion.h.

Referenced by Conversion(), and tracksPout().

std::vector<reco::TrackRef> reco::Conversion::tracks_ [private]

vector of Track references

Definition at line 195 of file Conversion.h.

Referenced by tracks().

vector Track RefToBase

Definition at line 197 of file Conversion.h.

Referenced by tracks().