CMS 3D CMS Logo

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 
20  LogDebug("ConversionVertexFinder") << "ConversionVertexFinder CTOR "
21  << "\n";
22  maxDelta_ = conf_.getParameter<double>("maxDelta");
23  maxReducedChiSq_ = conf_.getParameter<double>("maxReducedChiSq");
24  minChiSqImprovement_ = conf_.getParameter<double>("minChiSqImprovement");
25  maxNbrOfIterations_ = conf_.getParameter<int>("maxNbrOfIterations");
28 }
29 
31  LogDebug("ConversionVertexFinder") << "ConversionVertexFinder DTOR "
32  << "\n";
33  delete kcvFitter_;
34 }
35 
36 bool ConversionVertexFinder::run(const std::vector<reco::TransientTrack>& _pair, reco::Vertex& the_vertex) {
37  std::vector<reco::TransientTrack> pair = _pair;
38  bool found = false;
39 
40  if (pair.size() < 2)
41  return found;
42 
43  float sigma = 0.00000001;
44  float chi = 0.;
45  float ndf = 0.;
46  float mass = 0.000511;
47 
49 
50  std::vector<RefCountedKinematicParticle> particles;
51 
52  particles.push_back(
53  pFactory.particle(pair[0], mass, chi, ndf, sigma, *pair[0].innermostMeasurementState().freeState()));
54  particles.push_back(
55  pFactory.particle(pair[1], mass, chi, ndf, sigma, *pair[1].innermostMeasurementState().freeState()));
56 
57 #ifdef OldKineFit
59 
60  RefCountedKinematicTree myTree = kcvFitter_->fit(particles, &constr);
61 #else
62 
63  // bizzare way to the get the field...
64  const MagneticField* mf = pair[0].field();
65 
68  kcvFitter.setParameters(conf_);
69  RefCountedKinematicTree myTree = kcvFitter.fit(particles, &constr);
70 
71 #ifdef KineFitDebug
72 
74 
75  RefCountedKinematicTree oldTree = kcvFitter_->fit(particles, &oldconstr);
76 
77  if (oldTree->isValid()) {
78  std::cout << "old " << kcvFitter_->getNit() << std::endl;
79  RefCountedKinematicVertex gamma_dec_vertex = oldTree->currentDecayVertex();
80  if (gamma_dec_vertex->vertexIsValid())
81  std::cout << gamma_dec_vertex->chiSquared() << " " << gamma_dec_vertex->degreesOfFreedom() << std::endl;
82  std::cout << oldTree->currentParticle()->currentState().globalMomentum()
83  << oldTree->currentParticle()->currentState().globalPosition() << std::endl;
84  std::vector<RefCountedKinematicParticle> fStates = oldTree->finalStateParticles();
85  for (unsigned int kk = 0; kk < fStates.size(); kk++) {
86  std::cout << fStates[kk]->currentState().globalMomentum() << fStates[kk]->currentState().globalPosition()
87  << std::endl;
88  std::cout << fStates[kk]->currentState().kinematicParametersError().matrix() << std::endl;
89  }
90  } else
91  std::cout << "old invalid " << kcvFitter_->getNit() << std::endl;
92 
93  if (myTree->isValid()) {
94  std::cout << "new " << kcvFitter.getNit() << std::endl;
95  RefCountedKinematicVertex gamma_dec_vertex = myTree->currentDecayVertex();
96  if (gamma_dec_vertex->vertexIsValid())
97  std::cout << gamma_dec_vertex->chiSquared() << " " << gamma_dec_vertex->degreesOfFreedom() << std::endl;
98  std::cout << myTree->currentParticle()->currentState().globalMomentum()
99  << myTree->currentParticle()->currentState().globalPosition() << std::endl;
100  std::vector<RefCountedKinematicParticle> fStates = myTree->finalStateParticles();
101  for (unsigned int kk = 0; kk < fStates.size(); kk++) {
102  std::cout << fStates[kk]->currentState().globalMomentum() << fStates[kk]->currentState().globalPosition()
103  << std::endl;
104  std::cout << fStates[kk]->currentState().kinematicParametersError().matrix() << std::endl;
105  }
106  } else
107  std::cout << "new invalid " << kcvFitter.getNit() << std::endl;
108 
109 #endif // TemplateKineFitDebug
110 
111 #endif
112 
113  if (myTree->isValid()) {
114  myTree->movePointerToTheTop();
115  RefCountedKinematicParticle the_photon = myTree->currentParticle();
116  if (the_photon->currentState().isValid()) {
117  //const ParticleMass photon_mass = the_photon->currentState().mass();
118  RefCountedKinematicVertex gamma_dec_vertex;
119  gamma_dec_vertex = myTree->currentDecayVertex();
120  if (gamma_dec_vertex->vertexIsValid()) {
121  const float chi2Prob =
122  ChiSquaredProbability(gamma_dec_vertex->chiSquared(), gamma_dec_vertex->degreesOfFreedom());
123  if (chi2Prob > 0.) { // no longer cut here, only ask positive probability here
124  //const math::XYZPoint vtxPos(gamma_dec_vertex->position());
125  the_vertex = *gamma_dec_vertex;
126  found = true;
127  }
128  }
129  }
130  }
131 
132  return found;
133 }
134 
135 TransientVertex ConversionVertexFinder::run(const std::vector<reco::TransientTrack>& pair) {
136  LogDebug("ConversionVertexFinder") << "ConversionVertexFinder run pair size " << pair.size() << "\n";
137 
138  //for ( std::vector<reco::TransientTrack>::const_iterator iTk=pair.begin(); iTk!=pair.end(); ++iTk) {
139  // LogDebug("ConversionVertexFinder") << " ConversionVertexFinder Tracks in the pair charge " << iTk->charge() << " Num of RecHits " << iTk->recHitsSize() << " inner momentum " << iTk->track().innerMomentum() << "\n";
140  //}
141 
142  reco::Vertex theVertex;
143  KalmanVertexFitter fitter(true);
144  TransientVertex transientVtx;
145 
146  const std::string metname = "ConversionVertexFinder| ConversionVertexFinder";
147  try {
148  transientVtx = fitter.vertex(pair);
149 
150  } catch (cms::Exception& e) {
151  edm::LogWarning(metname) << "cms::Exception caught in ConversionVertexFinder::run\n" << e.explainSelf();
152  }
153 
154  return transientVtx;
155 }
KalmanVertexFitter::vertex
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const override
Definition: KalmanVertexFitter.h:49
MessageLogger.h
ConversionVertexFinder::~ConversionVertexFinder
~ConversionVertexFinder()
Definition: ConversionVertexFinder.cc:30
KinematicConstrainedVertexFitter::setParameters
void setParameters(const edm::ParameterSet &pSet)
Definition: KinematicConstrainedVertexFitter.cc:34
KinematicParticleFactoryFromTransientTrack
Definition: KinematicParticleFactoryFromTransientTrack.h:16
ChiSquaredProbability
float ChiSquaredProbability(double chiSquared, double nrDOF)
Definition: ChiSquaredProbability.cc:13
gather_cfg.cout
cout
Definition: gather_cfg.py:144
KinematicConstrainedVertexFitterT::fit
RefCountedKinematicTree fit(const std::vector< RefCountedKinematicParticle > &part)
Definition: KinematicConstrainedVertexFitterT.h:50
ConversionVertexFinder.h
ChiSquaredProbability.h
ReferenceCountingPointer< KinematicTree >
newFWLiteAna.found
found
Definition: newFWLiteAna.py:118
ConversionVertexFinder::kcvFitter_
KinematicConstrainedVertexFitter * kcvFitter_
Definition: ConversionVertexFinder.h:45
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
ecalTrigSettings_cff.particles
particles
Definition: ecalTrigSettings_cff.py:11
KinematicConstrainedVertexFitterT::getNit
int getNit() const
Definition: KinematicConstrainedVertexFitterT.h:275
config
Definition: config.py:1
KinematicConstrainedVertexFitter::fit
RefCountedKinematicTree fit(const std::vector< RefCountedKinematicParticle > &part)
Definition: KinematicConstrainedVertexFitter.h:46
GetRecoTauVFromDQM_MC_cff.kk
kk
Definition: GetRecoTauVFromDQM_MC_cff.py:84
KinematicParticleFactoryFromTransientTrack::particle
RefCountedKinematicParticle particle(const reco::TransientTrack &initialTrack, const ParticleMass &massGuess, float chiSquared, float degreesOfFr, float &m_sigma) const
Definition: KinematicParticleFactoryFromTransientTrack.cc:15
ConversionVertexFinder::ConversionVertexFinder
ConversionVertexFinder(const edm::ParameterSet &config)
Definition: ConversionVertexFinder.cc:19
ConversionVertexFinder::maxDelta_
double maxDelta_
Definition: ConversionVertexFinder.h:41
ConversionVertexFinder::run
TransientVertex run(const std::vector< reco::TransientTrack > &pair)
Definition: ConversionVertexFinder.cc:135
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
edm::ParameterSet
Definition: ParameterSet.h:47
KinematicConstrainedVertexFitter
Definition: KinematicConstrainedVertexFitter.h:21
ConversionVertexFinder::maxReducedChiSq_
double maxReducedChiSq_
Definition: ConversionVertexFinder.h:42
ConversionVertexFinder::minChiSqImprovement_
double minChiSqImprovement_
Definition: ConversionVertexFinder.h:43
KinematicConstrainedVertexFitterT.h
ColinearityKinematicConstraint::PhiTheta
Definition: ColinearityKinematicConstraint.h:18
ColinearityKinematicConstraintT.h
TransientVertex
Definition: TransientVertex.h:18
ConversionVertexFinder::maxNbrOfIterations_
int maxNbrOfIterations_
Definition: ConversionVertexFinder.h:44
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
KinematicConstrainedVertexFitter::getNit
int getNit() const
Definition: KinematicConstrainedVertexFitter.cc:195
KinematicConstrainedVertexFitterT::setParameters
void setParameters(const edm::ParameterSet &pSet)
Definition: KinematicConstrainedVertexFitterT.h:129
ColinearityKinematicConstraint
Definition: ColinearityKinematicConstraint.h:16
ConversionVertexFinder::conf_
edm::ParameterSet conf_
Definition: ConversionVertexFinder.h:40
EgHLTOffHistBins_cfi.mass
mass
Definition: EgHLTOffHistBins_cfi.py:34
KinematicParticleFactoryFromTransientTrack.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
ColinearityKinematicConstraintT
Definition: ColinearityKinematicConstraintT.h:23
cms::Exception
Definition: Exception.h:70
ParameterSet.h
MagneticField
Definition: MagneticField.h:19
KinematicConstrainedVertexFitterT
Definition: KinematicConstrainedVertexFitterT.h:25
reco::Vertex
Definition: Vertex.h:35
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
metname
const std::string metname
Definition: MuonSeedOrcaPatternRecognition.cc:40
KalmanVertexFitter
Definition: KalmanVertexFitter.h:22