CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ConversionVertexFinder.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <memory>
4 // Framework
7 //
8 //Kinematic constraint vertex fitter
10 
12 
13 // new templated one
16 
17 
18 
19 //
20 
21 
22 
23 
24 
26  conf_(config)
27 {
28  LogDebug("ConversionVertexFinder") << "ConversionVertexFinder CTOR " << "\n";
29  maxDelta_ = conf_.getParameter<double>("maxDelta");
30  maxReducedChiSq_ = conf_.getParameter<double>("maxReducedChiSq");
31  minChiSqImprovement_ = conf_.getParameter<double>("minChiSqImprovement");
32  maxNbrOfIterations_ = conf_.getParameter<int>("maxNbrOfIterations");
34  kcvFitter_->setParameters(conf_);
35 
36 }
37 
39 
40  LogDebug("ConversionVertexFinder") << "ConversionVertexFinder DTOR " << "\n";
41  delete kcvFitter_;
42 
43 }
44 
45 bool ConversionVertexFinder::run(std::vector<reco::TransientTrack> pair, reco::Vertex& the_vertex) {
46  bool found= false;
47 
48  if ( pair.size() < 2) return found;
49 
50  float sigma = 0.00000001;
51  float chi = 0.;
52  float ndf = 0.;
53  float mass = 0.000511;
54 
55 
57 
58  std::vector<RefCountedKinematicParticle> particles;
59 
60  particles.push_back(pFactory.particle (pair[0],mass,chi,ndf,sigma,*pair[0].innermostMeasurementState().freeState()));
61  particles.push_back(pFactory.particle (pair[1],mass,chi,ndf,sigma,*pair[1].innermostMeasurementState().freeState()));
62 
63 
64 #ifdef OldKineFit
66 
67  RefCountedKinematicTree myTree = kcvFitter_->fit(particles, &constr);
68 #else
69 
70  // bizzare way to the get the field...
71  const MagneticField* mf = pair[0].field();
72 
75  kcvFitter.setParameters(conf_);
76  RefCountedKinematicTree myTree = kcvFitter.fit(particles, &constr);
77 
78 #ifdef KineFitDebug
79 
81 
82  RefCountedKinematicTree oldTree = kcvFitter_->fit(particles, &oldconstr);
83 
84 
85  if( oldTree->isValid() ) {
86  std::cout << "old " << kcvFitter_->getNit() << std::endl;
87  RefCountedKinematicVertex gamma_dec_vertex = oldTree->currentDecayVertex();
88  if( gamma_dec_vertex->vertexIsValid())
89  std::cout << gamma_dec_vertex->chiSquared() << " " << gamma_dec_vertex->degreesOfFreedom() << std::endl;
90  std::cout << oldTree->currentParticle()->currentState().globalMomentum() <<
91  oldTree->currentParticle()->currentState().globalPosition()<< std::endl;
92  std::vector<RefCountedKinematicParticle> fStates=oldTree->finalStateParticles();
93  for (unsigned int kk=0; kk<fStates.size(); kk++) {
94  std::cout << fStates[kk]->currentState().globalMomentum() <<
95  fStates[kk]->currentState().globalPosition() << std::endl;
96  std::cout << fStates[kk]->currentState().kinematicParametersError().matrix()<<std::endl;
97  }
98  } else std::cout << "old invalid " << kcvFitter_->getNit() << std::endl;
99 
100  if( myTree->isValid() ) {
101  std::cout << "new " << kcvFitter.getNit() << std::endl;
102  RefCountedKinematicVertex gamma_dec_vertex = myTree->currentDecayVertex();
103  if( gamma_dec_vertex->vertexIsValid())
104  std::cout << gamma_dec_vertex->chiSquared() << " " << gamma_dec_vertex->degreesOfFreedom() << std::endl;
105  std::cout << myTree->currentParticle()->currentState().globalMomentum() <<
106  myTree->currentParticle()->currentState().globalPosition()<< std::endl;
107  std::vector<RefCountedKinematicParticle> fStates=myTree->finalStateParticles();
108  for (unsigned int kk=0; kk<fStates.size(); kk++) {
109  std::cout << fStates[kk]->currentState().globalMomentum() <<
110  fStates[kk]->currentState().globalPosition() << std::endl;
111  std::cout << fStates[kk]->currentState().kinematicParametersError().matrix()<<std::endl;
112  }
113  } else std::cout << "new invalid " << kcvFitter.getNit() << std::endl;
114 
115 #endif // TemplateKineFitDebug
116 
117 #endif
118 
119  if( myTree->isValid() ) {
120  myTree->movePointerToTheTop();
121  RefCountedKinematicParticle the_photon = myTree->currentParticle();
122  if (the_photon->currentState().isValid()){
123  //const ParticleMass photon_mass = the_photon->currentState().mass();
124  RefCountedKinematicVertex gamma_dec_vertex;
125  gamma_dec_vertex = myTree->currentDecayVertex();
126  if( gamma_dec_vertex->vertexIsValid() ){
127  const float chi2Prob = ChiSquaredProbability(gamma_dec_vertex->chiSquared(), gamma_dec_vertex->degreesOfFreedom());
128  if (chi2Prob>0.){// no longer cut here, only ask positive probability here
129  //const math::XYZPoint vtxPos(gamma_dec_vertex->position());
130  the_vertex = *gamma_dec_vertex;
131  found = true;
132  }
133  }
134  }
135  }
136 
137  return found;
138 }
139 
140 
141 TransientVertex ConversionVertexFinder::run(std::vector<reco::TransientTrack> pair) {
142  LogDebug("ConversionVertexFinder") << "ConversionVertexFinder run pair size " << pair.size() << "\n";
143 
144  //for ( std::vector<reco::TransientTrack>::const_iterator iTk=pair.begin(); iTk!=pair.end(); ++iTk) {
145  // LogDebug("ConversionVertexFinder") << " ConversionVertexFinder Tracks in the pair charge " << iTk->charge() << " Num of RecHits " << iTk->recHitsSize() << " inner momentum " << iTk->track().innerMomentum() << "\n";
146  //}
147 
148 
149  reco::Vertex theVertex;
150  KalmanVertexFitter fitter(true);
151  TransientVertex transientVtx;
152 
153  const std::string metname = "ConversionVertexFinder| ConversionVertexFinder";
154  try{
155 
156  transientVtx = fitter.vertex(pair);
157 
158  } catch ( cms::Exception& e ) {
159 
160 
161  edm::LogWarning(metname) << "cms::Exception caught in ConversionVertexFinder::run\n"
162  << e.explainSelf();
163 
164  }
165 
166 
167  return transientVtx;
168 
169 
170 
171 }
172 
#define LogDebug(id)
ConversionVertexFinder(const edm::ParameterSet &config)
T getParameter(std::string const &) const
RefCountedKinematicTree fit(std::vector< RefCountedKinematicParticle > part)
const std::string metname
virtual std::string explainSelf() const
Definition: Exception.cc:146
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const
void setParameters(const edm::ParameterSet &pSet)
TransientVertex run(std::vector< reco::TransientTrack > pair)
float ChiSquaredProbability(double chiSquared, double nrDOF)
RefCountedKinematicTree fit(std::vector< RefCountedKinematicParticle > part)
tuple mass
Definition: scaleCards.py:27
KinematicConstrainedVertexFitter * kcvFitter_
RefCountedKinematicParticle particle(const reco::TransientTrack &initialTrack, const ParticleMass &massGuess, float chiSquared, float degreesOfFr, float &m_sigma) const
tuple cout
Definition: gather_cfg.py:121