CMS 3D CMS Logo

KinematicParticleVertexFitter.cc
Go to the documentation of this file.
2 // #include "Vertex/LinearizationPointFinders/interface/LMSLinearizationPointFinder.h"
14 
16 {
18  setup(pSet);
19 }
20 
22 {
23  setup(pSet);
24 }
25 
26 void
28 {
29 
32 
36  SequentialVertexSmoother<6> smoother(vtu, vse, covCalc);
37  fitter
40 }
41 
43 {
44  delete vFactory;
45  delete pointFinder;
46  delete fitter;
47 }
48 
50 {
51  edm::ParameterSet pSet;
52  pSet.addParameter<double>("maxDistance", 0.01);
53  pSet.addParameter<int>("maxNbrOfIterations", 100); //10
54  return pSet;
55 }
56 
57 RefCountedKinematicTree KinematicParticleVertexFitter::fit(const std::vector<RefCountedKinematicParticle> &particles) const
58 {
59  typedef ReferenceCountingPointer<VertexTrack<6> > RefCountedVertexTrack;
60 //sorting the input
61  if(particles.size()<2) throw VertexException("KinematicParticleVertexFitter::input states are less than 2");
62  InputSort iSort;
63  std::pair<std::vector<RefCountedKinematicParticle>, std::vector<FreeTrajectoryState> > input = iSort.sort(particles);
64  std::vector<RefCountedKinematicParticle> & newPart = input.first;
65  std::vector<FreeTrajectoryState> & freeStates = input.second;
66 
67  GlobalPoint linPoint = pointFinder->getLinearizationPoint(freeStates);
68 
69 // cout<<"Linearization point found"<<endl;
70 
71 //making initial veretx seed with lin point as position and a fake error
73  we(0,0)=we(1,1)=we(2,2) = 10000.;
74  GlobalError error(we);
75  VertexState state(linPoint, error);
76 
77 //vector of Vertex Tracks to fit
78  std::vector<RefCountedVertexTrack> ttf;
80  for(auto const & i : newPart){
81  if( !(i)->currentState().isValid() || !propagator_.willPropagateToTheTransversePCA((i)->currentState(), linPoint) ) {
82  // std::cout << "Here's the bad state." << std::endl;
83  return ReferenceCountingPointer<KinematicTree>(new KinematicTree()); // return invalid vertex
84  }
85  ttf.push_back(vFactory->vertexTrack((i)->particleLinearizedTrackState(linPoint),state,1.));
86  }
87 
88 // //debugging code to check neutrals:
89 // for(std::vector<RefCountedVertexTrack>::const_iterator i = ttf.begin(); i!=ttf.end(); i++)
90 // {
91 // // cout<<"predicted state momentum error"<<(*i)->linearizedTrack()->predictedStateMomentumError()<<endl;
92 // // cout<<"Momentum jacobian"<<(*i)->linearizedTrack()->momentumJacobian() <<endl;
93 // // cout<<"predicted state momentum "<<(*i)->linearizedTrack()->predictedStateMomentum()<<endl;
94 // // cout<<"constant term"<<(*i)->linearizedTrack()->constantTerm()<<endl;
95 //
96 // }
97 
99  if (!vtx.isValid()) {
100  LogDebug("RecoVertex/KinematicParticleVertexFitter")
101  << "Fitted position is invalid. Returned Tree is invalid\n";
102  return ReferenceCountingPointer<KinematicTree>(new KinematicTree()); // return invalid vertex
103  }
104  FinalTreeBuilder tBuilder;
105  return tBuilder.buildTree(vtx, newPart);
106 }
#define LogDebug(id)
edm::ParameterSet defaultParameters() const
void setup(const edm::ParameterSet &pSet)
RefCountedVertexTrack vertexTrack(const RefCountedLinearizedTrackState lt, const VertexState vs, float weight=1.0) const
Common base class.
static std::string const input
Definition: EdmProvDump.cc:44
bool willPropagateToTheTransversePCA(const KinematicState &state, const GlobalPoint &point) const override
std::pair< std::vector< RefCountedKinematicParticle >, std::vector< FreeTrajectoryState > > sort(const std::vector< RefCountedKinematicParticle > &particles) const
Definition: InputSort.cc:6
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:144
virtual CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const =0
RefCountedKinematicTree fit(const std::vector< RefCountedKinematicParticle > &particles) const
virtual GlobalPoint getLinearizationPoint(const std::vector< reco::TransientTrack > &) const =0
bool isValid() const