CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Protected Member Functions | Protected Attributes
CandCommonVertexFitterBase Class Referenceabstract

#include <CandCommonVertexFitter.h>

Inheritance diagram for CandCommonVertexFitterBase:
CandCommonVertexFitter< Fitter >

Public Types

typedef
reco::Vertex::CovarianceMatrix 
CovarianceMatrix
 

Public Member Functions

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

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 More...
 
CovarianceMatrix cov_
 covariance matrix (3x3) More...
 
double ndof_
 number of degrees of freedom More...
 

Detailed Description

Definition at line 17 of file CandCommonVertexFitter.h.

Member Typedef Documentation

Definition at line 19 of file CandCommonVertexFitter.h.

Constructor & Destructor Documentation

CandCommonVertexFitterBase::CandCommonVertexFitterBase ( const edm::ParameterSet )
inline

Definition at line 20 of file CandCommonVertexFitter.h.

20 : bField_(0) { }
virtual CandCommonVertexFitterBase::~CandCommonVertexFitterBase ( )
inlinevirtual

Definition at line 21 of file CandCommonVertexFitter.h.

21 { }

Member Function Documentation

void CandCommonVertexFitterBase::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 CandCommonVertexFitter.cc.

References dtNoiseDBValidation_cfg::cerr, reco::Candidate::daughter(), edm::hlt::Exception, lumiContext::fill, reco::Candidate::get(), edm::errors::InvalidReference, j, python.rootplot.argparse::message, reco::Candidate::numberOfDaughters(), and reco::Candidate::pdgId().

61  {
62  size_t nDau = c.numberOfDaughters();
63  for(unsigned int j = 0; j < nDau ; ++j) {
64  Candidate * d = c.daughter(j);
65  if(d == 0) {
66  ostringstream message;
67  message << "Can't access in write mode candidate daughters. "
68  << "pdgId = " << c.pdgId() << ".\n";
69  const Candidate * d1 = c.daughter(j);
70  if(d1 == 0)
71  message << "Null daughter also found in read-only mode\n";
72  else
73  message << "Daughter found in read-only mode with id: " << d1->pdgId() << "\n";
74  throw edm::Exception(edm::errors::InvalidReference) << message.str();
75  }
76  if(d->numberOfDaughters() > 0)
77  fill(tracks, daughters, trackTypes, * d);
78  else {
79  const Track * trk = d->get<const Track *>();
81  if(trk != 0) {
82  tracks.push_back(TransientTrack(* trk, bField_));
83  daughters.push_back(d);
84  trackTypes.push_back(type);
85  } else {
86  cerr << ">>> warning: candidate of type " << d->pdgId()
87  << " has no track reference." << endl;
88  }
89  }
90  }
91 }
type
Definition: HCALResponse.h:21
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
TrackType
track type
Definition: RecoCandidate.h:57
virtual size_type numberOfDaughters() const =0
number of daughters
int j
Definition: DBlmapReader.cc:9
void fill(std::vector< reco::TransientTrack > &, std::vector< reco::Candidate * > &, std::vector< reco::RecoCandidate::TrackType > &, reco::Candidate &) const
virtual int pdgId() const =0
PDG identifier.
tuple tracks
Definition: testEve_cfg.py:39
T get() const
get a component
Definition: Candidate.h:219
virtual bool CandCommonVertexFitterBase::fit ( TransientVertex ,
const std::vector< reco::TransientTrack > &   
) const
protectedpure virtual
void CandCommonVertexFitterBase::set ( const MagneticField bField)
inline
void CandCommonVertexFitterBase::set ( reco::VertexCompositeCandidate c) const

Definition at line 9 of file CandCommonVertexFitter.cc.

References reco::Vertex::chi2(), relval_parameters_module::energy, reco::Candidate::energy(), edm::hlt::Exception, reco::Vertex::fill(), lumiContext::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(), pileupReCalc_HLTpaths::scale, reco::VertexCompositeCandidate::setChi2AndNdof(), reco::VertexCompositeCandidate::setCovariance(), reco::Candidate::setP4(), reco::LeafCandidate::setP4(), reco::Candidate::setVertex(), reco::LeafCandidate::setVertex(), mathSSE::sqrt(), testEve_cfg::tracks, and findQualityFiles::v.

Referenced by betterConfigParser.BetterConfigParser::getGeneral().

9  {
10  if(bField_ == 0)
12  << "B-Field was not set up CandCommonVertexFitter.\n"
13  << "the following method must be called before fitting a candidate:\n"
14  << " CandCommonVertexFitter:.set( const MagneticField * )" << endl;
15  vector<TransientTrack> tracks;
16  vector<Candidate *> daughters;
17  vector<RecoCandidate::TrackType> trackTypes;
18  fill(tracks, daughters, trackTypes, c);
19  assert(tracks.size() == daughters.size());
20  TransientVertex vertex;
21  if(fit(vertex, tracks)) {
22  tracks = vertex.refittedTracks();
23  Candidate::Point vtx(vertex.position());
24  c.setVertex(vtx);
25  vector<TransientTrack>::const_iterator trackIt = tracks.begin(), tracksEnd = tracks.end();
26  vector<Candidate *>::const_iterator daughterIt = daughters.begin();
27  vector<RecoCandidate::TrackType>::const_iterator trackTypeIt = trackTypes.begin();
28  Candidate::LorentzVector mp4(0, 0, 0, 0);
29  for(; trackIt != tracksEnd; ++ trackIt, ++ daughterIt, ++trackTypeIt) {
30  const Track & track = trackIt->track();
31  Candidate & daughter = * * daughterIt;
32  double px = track.px(), py = track.py(), pz = track.pz(), p = track.p();
33  double energy;
34  daughter.setVertex( vtx );
35  if(*trackTypeIt == RecoCandidate::recoTrackType) {
36  double mass = daughter.mass();
37  energy = sqrt(p*p + mass*mass);
38  } else {
39  energy = daughter.energy();
40  double scale = energy / p;
41  px *= scale; py *= scale; pz *= scale;
42  }
43  Candidate::LorentzVector dp4(px, py, pz, energy);
44  daughter.setP4(dp4);
45  mp4 += dp4;
46  }
47  c.setP4(mp4);
48  Vertex v = vertex;
49  c.setChi2AndNdof(chi2_ = v.chi2(), ndof_ = v.ndof());
50  v.fill(cov_);
52  } else {
53  c.setChi2AndNdof(chi2_ = -1, ndof_ = 0);
54  c.setCovariance(cov_ = CovarianceMatrix(ROOT::Math::SMatrixIdentity()));
55  }
56 }
double p() const
momentum vector magnitude
Definition: TrackBase.h:127
virtual double energy() const =0
energy
virtual float mass() const =0
mass
virtual void setP4(const LorentzVector &p4)=0
set 4-momentum
virtual void setP4(const LorentzVector &p4)
set 4-momentum
void fill(CovarianceMatrix &v) const
fill SMatrix
Definition: Vertex.cc:27
reco::Vertex::CovarianceMatrix CovarianceMatrix
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:131
CovarianceMatrix cov_
covariance matrix (3x3)
void setChi2AndNdof(double chi2, double ndof)
set chi2 and ndof
T sqrt(T t)
Definition: SSEVec.h:48
virtual void setVertex(const Point &vertex)=0
set vertex
void setCovariance(const CovarianceMatrix &m)
set covariance matrix
double chi2() const
chi-squares
Definition: Vertex.h:95
virtual void setVertex(const Point &vertex)
set vertex
void fill(std::vector< reco::TransientTrack > &, std::vector< reco::Candidate * > &, std::vector< reco::RecoCandidate::TrackType > &, reco::Candidate &) const
double ndof() const
Definition: Vertex.h:102
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:135
tuple tracks
Definition: testEve_cfg.py:39
double ndof_
number of degrees of freedom
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
math::XYZPoint Point
point in the space
Definition: Candidate.h:45
std::vector< reco::TransientTrack > const & refittedTracks() const
virtual bool fit(TransientVertex &, const std::vector< reco::TransientTrack > &) const =0
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:133

Member Data Documentation

const MagneticField* CandCommonVertexFitterBase::bField_
protected

Definition at line 26 of file CandCommonVertexFitter.h.

Referenced by set().

double CandCommonVertexFitterBase::chi2_
mutableprotected

chi-sqared

Definition at line 34 of file CandCommonVertexFitter.h.

CovarianceMatrix CandCommonVertexFitterBase::cov_
mutableprotected

covariance matrix (3x3)

Definition at line 38 of file CandCommonVertexFitter.h.

double CandCommonVertexFitterBase::ndof_
mutableprotected

number of degrees of freedom

Definition at line 36 of file CandCommonVertexFitter.h.