CMS 3D CMS Logo

Public Types | Public Member Functions | Protected Member Functions | Protected Attributes

PFCandCommonVertexFitterBase Class Reference

#include <PFCandCommonVertexFitter.h>

Inheritance diagram for PFCandCommonVertexFitterBase:
PFCandCommonVertexFitter< Fitter >

List of all members.

Public Types

typedef
reco::Vertex::CovarianceMatrix 
CovarianceMatrix

Public Member Functions

 PFCandCommonVertexFitterBase (const edm::ParameterSet &)
void set (reco::VertexCompositeCandidate &) const
void set (const MagneticField *bField)
virtual ~PFCandCommonVertexFitterBase ()

Protected Member Functions

void fill (std::vector< reco::TransientTrack > &, std::vector< reco::Candidate * > &, std::vector< reco::RecoCandidate::TrackType > &, reco::Candidate &) const
virtual bool fit (TransientVertex &, const std::vector< reco::TransientTrack > &) const =0

Protected Attributes

const MagneticFieldbField_
double chi2_
 chi-sqared
CovarianceMatrix cov_
 covariance matrix (3x3)
double ndof_
 number of degrees of freedom

Detailed Description

Definition at line 23 of file PFCandCommonVertexFitter.h.


Member Typedef Documentation

Definition at line 25 of file PFCandCommonVertexFitter.h.


Constructor & Destructor Documentation

PFCandCommonVertexFitterBase::PFCandCommonVertexFitterBase ( const edm::ParameterSet ) [inline]

Definition at line 26 of file PFCandCommonVertexFitter.h.

: bField_(0) { }
virtual PFCandCommonVertexFitterBase::~PFCandCommonVertexFitterBase ( ) [inline, virtual]

Definition at line 27 of file PFCandCommonVertexFitter.h.

{ }

Member Function Documentation

void PFCandCommonVertexFitterBase::fill ( std::vector< reco::TransientTrack > &  tracks,
std::vector< reco::Candidate * > &  daughters,
std::vector< reco::RecoCandidate::TrackType > &  trackTypes,
reco::Candidate c 
) const [protected]

Definition at line 58 of file PFCandCommonVertexFitter.cc.

References benchmark_cfg::cerr, debug_cff::d1, reco::Candidate::daughter(), Exception, edm::RefToBase< T >::get(), edm::Ref< C, T, F >::get(), reco::Candidate::hasMasterClone(), edm::errors::InvalidReference, j, reco::Candidate::masterClone(), python::rootplot::argparse::message, NULL, reco::Candidate::numberOfDaughters(), reco::Candidate::pdgId(), reco::RecoCandidate::recoTrackType, and reco::PFCandidate::trackRef().

                                                           {
  size_t nDau = c.numberOfDaughters();
  for(unsigned int j = 0; j < nDau ; ++j) {
    Candidate * d = c.daughter(j);
    if(d == 0) {
      ostringstream message;
      message << "Can't access in write mode candidate daughters. "
              << "pdgId = " << c.pdgId() << ".\n";
      const Candidate * d1 = c.daughter(j);
      if(d1 == 0)
        message << "Null daughter also found in read-only mode\n";
      else
        message << "Daughter found in read-only mode with id: " << d1->pdgId() << "\n";
      throw edm::Exception(edm::errors::InvalidReference) << message.str();
    }
    if(d->numberOfDaughters() > 0)
      fill(tracks, daughters, trackTypes, * d);
    else {
       const Track * trk = NULL;
       RecoCandidate::TrackType type = RecoCandidate::recoTrackType;
       if (d->hasMasterClone())
       {
          //get the PFCandidate
          const PFCandidate* myPFCand = dynamic_cast<const PFCandidate*>(d->masterClone().get());
          trk = myPFCand->trackRef().get();
       }
      if(trk != 0) {
        tracks.push_back(TransientTrack(* trk, bField_));
        daughters.push_back(d);
        trackTypes.push_back(type);
      } else {
        cerr << ">>> warning: candidate of type " << d->pdgId() 
             << " has no track reference." << endl;
      }
    }
  }
}
virtual bool PFCandCommonVertexFitterBase::fit ( TransientVertex ,
const std::vector< reco::TransientTrack > &   
) const [protected, pure virtual]
void PFCandCommonVertexFitterBase::set ( const MagneticField bField) [inline]
void PFCandCommonVertexFitterBase::set ( reco::VertexCompositeCandidate c) const

Definition at line 9 of file PFCandCommonVertexFitter.cc.

References reco::Vertex::chi2(), reco::Candidate::energy(), relval_parameters_module::energy, Exception, reco::Vertex::fill(), edm::errors::InvalidReference, reco::Candidate::mass(), reco::Vertex::ndof(), reco::TrackBase::p(), AlCaHLTBitMon_ParallelJobs::p, reco::TrackBase::px(), reco::TrackBase::py(), reco::TrackBase::pz(), reco::RecoCandidate::recoTrackType, TransientVertex::refittedTracks(), reco::VertexCompositeCandidate::setChi2AndNdof(), reco::VertexCompositeCandidate::setCovariance(), reco::LeafCandidate::setP4(), reco::Candidate::setP4(), reco::LeafCandidate::setVertex(), reco::Candidate::setVertex(), mathSSE::sqrt(), testEve_cfg::tracks, and v.

                                                                         {
  if(bField_ == 0)
    throw edm::Exception(edm::errors::InvalidReference)
      << "B-Field was not set up PFCandCommonVertexFitter.\n"
      << "the following method must be called before fitting a candidate:\n"
      << " PFCandCommonVertexFitter:.set( const MagneticField * )" << endl;
  std::vector<TransientTrack> tracks;
  std::vector<Candidate *> daughters;
  std::vector<RecoCandidate::TrackType> trackTypes;
  fill(tracks, daughters, trackTypes, c);
  assert(tracks.size() == daughters.size());
  TransientVertex vertex;
  if(fit(vertex, tracks)) {
    tracks = vertex.refittedTracks();    
    Candidate::Point vtx(vertex.position());
    c.setVertex(vtx);
    std::vector<TransientTrack>::const_iterator trackIt = tracks.begin(), tracksEnd = tracks.end();
    std::vector<Candidate *>::const_iterator daughterIt = daughters.begin();
    std::vector<RecoCandidate::TrackType>::const_iterator trackTypeIt = trackTypes.begin();
    Candidate::LorentzVector mp4(0, 0, 0, 0);
    for(; trackIt != tracksEnd; ++ trackIt, ++ daughterIt, ++trackTypeIt) {
      const Track & track = trackIt->track();
      Candidate & daughter = * * daughterIt;
      double px = track.px(), py = track.py(), pz = track.pz(), p = track.p();
      double energy;
      daughter.setVertex( vtx );
      if(*trackTypeIt == RecoCandidate::recoTrackType) {
        double mass = daughter.mass();
        energy = sqrt(p*p + mass*mass);
      } else {
        energy = daughter.energy();
        double scale = energy / p;
        px *= scale; py *= scale; pz *= scale; 
      }
      Candidate::LorentzVector dp4(px, py, pz, energy);
      daughter.setP4(dp4);
      mp4 += dp4;
    }
    c.setP4(mp4);
    Vertex v = vertex;
    c.setChi2AndNdof(chi2_ = v.chi2(), ndof_ = v.ndof());
    v.fill(cov_);
    c.setCovariance(cov_);
  } else {
    c.setChi2AndNdof(chi2_ = -1, ndof_ = 0);
    c.setCovariance(cov_ = CovarianceMatrix(ROOT::Math::SMatrixIdentity())); 
  }
}

Member Data Documentation

Definition at line 32 of file PFCandCommonVertexFitter.h.

Referenced by set().

double PFCandCommonVertexFitterBase::chi2_ [mutable, protected]

chi-sqared

Definition at line 40 of file PFCandCommonVertexFitter.h.

covariance matrix (3x3)

Definition at line 44 of file PFCandCommonVertexFitter.h.

double PFCandCommonVertexFitterBase::ndof_ [mutable, protected]

number of degrees of freedom

Definition at line 42 of file PFCandCommonVertexFitter.h.