CMS 3D CMS Logo

KineExample.cc
Go to the documentation of this file.
2 
11 
17 //#include "MagneticField/Engine/interface/MagneticField.h"
18 //#include "MagneticField/Records/interface/IdealMagneticField.h"
19 
23 // #include "RecoVertex/KinematicFitPrimitives/interface/"
29 
30 #include <iostream>
31 
32 using namespace reco;
33 using namespace edm;
34 using namespace std;
35 
36 KineExample::KineExample(const edm::ParameterSet& iConfig) : theConfig(iConfig) {
37  token_tracks = consumes<TrackCollection>(iConfig.getParameter<string>("TrackLabel"));
38  outputFile_ = iConfig.getUntrackedParameter<std::string>("outputFile");
39  kvfPSet = iConfig.getParameter<edm::ParameterSet>("KVFParameters");
40  // rootFile_ = TFile::Open(outputFile_.c_str(),"RECREATE");
41  edm::LogInfo("RecoVertex/KineExample") << "Initializing KVF TEST analyser - Output file: " << outputFile_ << "\n";
42  // token_TrackTruth = consumes<TrackingParticleCollection>(edm::InputTag("trackingtruth", "TrackTruth"));
43  token_VertexTruth = consumes<TrackingVertexCollection>(edm::InputTag("trackingtruth", "VertexTruth"));
44 }
45 
47  // delete rootFile_;
48 }
49 
51  // edm::ESHandle<MagneticField> magField;
52  // setup.get<IdealMagneticFieldRecord>().get(magField);
53  // tree.reset(new SimpleVertexTree("VertexFitter", magField.product()));
54 }
55 
57 
58 //
59 // member functions
60 //
61 
63  try {
64  cout << "Reconstructing event number: " << iEvent.id() << "\n";
65 
66  // get RECO tracks from the event
67  // `tks` can be used as a ptr to a reco::TrackCollection
69  iEvent.getByToken(token_tracks, tks);
70  if (!tks.isValid()) {
71  cout << "Couln't find track collection: " << iEvent.id() << "\n";
72  } else {
73  edm::LogInfo("RecoVertex/KineExample") << "Found: " << (*tks).size() << " reconstructed tracks"
74  << "\n";
75  cout << "got " << (*tks).size() << " tracks " << endl;
76 
77  // Transform Track to TransientTrack
78  //get the builder:
80  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", theB);
81  //do the conversion:
82  vector<TransientTrack> t_tks = (*theB).build(tks);
83 
84  cout << "Found: " << t_tks.size() << " reconstructed tracks"
85  << "\n";
86 
87  // Do a KindFit, if >= 4 tracks.
88  if (t_tks.size() > 3) {
89  // For a first test, suppose that the first four tracks are the 2 muons,
90  // then the 2 kaons. Since this will not be true, the result of the fit
91  // will not be meaningfull, but at least you will get the idea of how to
92  // do such a fit.
93 
94  //First, to get started, a simple vertex fit:
95 
96  vector<TransientTrack> ttv;
97  ttv.push_back(t_tks[0]);
98  ttv.push_back(t_tks[1]);
99  ttv.push_back(t_tks[2]);
100  ttv.push_back(t_tks[3]);
101  KalmanVertexFitter kvf(false);
102  TransientVertex tv = kvf.vertex(ttv);
103  if (!tv.isValid())
104  cout << "KVF failed\n";
105  else
106  std::cout << "KVF fit Position: " << Vertex::Point(tv.position()) << "\n";
107 
108  TransientTrack ttMuPlus = t_tks[0];
109  TransientTrack ttMuMinus = t_tks[1];
110  TransientTrack ttKPlus = t_tks[2];
111  TransientTrack ttKMinus = t_tks[3];
112 
113  //the final state muons and kaons from the Bs->J/PsiPhi->mumuKK decay
114  //Creating a KinematicParticleFactory
116 
117  //The mass of a muon and the insignificant mass sigma to avoid singularities in the covariance matrix.
118  ParticleMass muon_mass = 0.1056583;
119  ParticleMass kaon_mass = 0.493677;
120  float muon_sigma = 0.0000001;
121  float kaon_sigma = 0.000016;
122 
123  //initial chi2 and ndf before kinematic fits. The chi2 of the reconstruction is not considered
124  float chi = 0.;
125  float ndf = 0.;
126 
127  //making particles
128  vector<RefCountedKinematicParticle> muonParticles;
129  vector<RefCountedKinematicParticle> phiParticles;
130  vector<RefCountedKinematicParticle> allParticles;
131  muonParticles.push_back(pFactory.particle(ttMuPlus, muon_mass, chi, ndf, muon_sigma));
132  muonParticles.push_back(pFactory.particle(ttMuMinus, muon_mass, chi, ndf, muon_sigma));
133  allParticles.push_back(pFactory.particle(ttMuPlus, muon_mass, chi, ndf, muon_sigma));
134  allParticles.push_back(pFactory.particle(ttMuMinus, muon_mass, chi, ndf, muon_sigma));
135 
136  phiParticles.push_back(pFactory.particle(ttKPlus, kaon_mass, chi, ndf, kaon_sigma));
137  phiParticles.push_back(pFactory.particle(ttKMinus, kaon_mass, chi, ndf, kaon_sigma));
138  allParticles.push_back(pFactory.particle(ttKPlus, kaon_mass, chi, ndf, kaon_sigma));
139  allParticles.push_back(pFactory.particle(ttKMinus, kaon_mass, chi, ndf, kaon_sigma));
140 
141  /* Example of a simple vertex fit, without other constraints
142  * The reconstructed decay tree is a result of the kinematic fit
143  * The KinematicParticleVertexFitter fits the final state particles to their vertex and
144  * reconstructs the decayed state
145  */
147  cout << "Simple vertex fit with KinematicParticleVertexFitter:\n";
148  RefCountedKinematicTree vertexFitTree = fitter.fit(allParticles);
149 
150  printout(vertexFitTree);
151 
153 
154  //creating the constraint for the J/Psi mass
155  ParticleMass jpsi = 3.09687;
156 
157  //creating the two track mass constraint
159 
160  //creating the fitter
162 
163  //obtaining the resulting tree
164  RefCountedKinematicTree myTree = kcvFitter.fit(allParticles, j_psi_c);
165 
166  cout << "\nGlobal fit done:\n";
167  printout(myTree);
168 
169  //creating the vertex fitter
171 
172  //reconstructing a J/Psi decay
173  RefCountedKinematicTree jpTree = kpvFitter.fit(muonParticles);
174 
175  //creating the particle fitter
176  KinematicParticleFitter csFitter;
177 
178  // creating the constraint
179  float jp_m_sigma = 0.00004;
180  KinematicConstraint* jpsi_c2 = new MassKinematicConstraint(jpsi, jp_m_sigma);
181 
182  //the constrained fit:
183  jpTree = csFitter.fit(jpsi_c2, jpTree);
184 
185  //getting the J/Psi KinematicParticle and putting it together with the kaons.
186  //The J/Psi KinematicParticle has a pointer to the tree it belongs to
187  jpTree->movePointerToTheTop();
188  RefCountedKinematicParticle jpsi_part = jpTree->currentParticle();
189  phiParticles.push_back(jpsi_part);
190 
191  //making a vertex fit and thus reconstructing the Bs parameters
192  // the resulting tree includes all the final state tracks, the J/Psi meson,
193  // its decay vertex, the Bs meson and its decay vertex.
194  RefCountedKinematicTree bsTree = kpvFitter.fit(phiParticles);
195  cout << "Sequential fit done:\n";
196  printout(bsTree);
197 
198  // // For the analysis: compare to your SimVertex
199  // TrackingVertex sv = getSimVertex(iEvent);
200  // edm::Handle<TrackingParticleCollection> TPCollectionH ;
201  // iEvent.getByToken(token_TrackTruth, TPCollectionH);
202  // const TrackingParticleCollection tPC = *(TPCollectionH.product());
203  // reco::RecoToSimCollection recSimColl=associatorForParamAtPca->associateRecoToSim(tks,
204  // TPCollectionH,
205  // &iEvent,
206  // &iSetup);
207  //
208  // tree->fill(tv, &sv, &recSimColl);
209  // }
210  }
211  }
212 
213  } catch (std::exception& err) {
214  cout << "Exception during event number: " << iEvent.id() << "\n" << err.what() << "\n";
215  }
216 }
217 
219  if (myVertex->vertexIsValid()) {
220  cout << "Decay vertex: " << myVertex->position() << myVertex->chiSquared() << " " << myVertex->degreesOfFreedom()
221  << endl;
222  } else
223  cout << "Decay vertex Not valid\n";
224 }
225 
226 void KineExample::printout(const RefCountedKinematicParticle& myParticle) const {
227  cout << "Particle: \n";
228  //accessing the reconstructed Bs meson parameters:
229  //SK: uncomment if needed AlgebraicVector7 bs_par = myParticle->currentState().kinematicParameters().vector();
230 
231  //and their joint covariance matrix:
232  //SK:uncomment if needed AlgebraicSymMatrix77 bs_er = myParticle->currentState().kinematicParametersError().matrix();
233  cout << "Momentum at vertex: " << myParticle->currentState().globalMomentum() << endl;
234  cout << "Parameters at vertex: " << myParticle->currentState().kinematicParameters().vector() << endl;
235 }
236 
238  if (!myTree->isValid()) {
239  cout << "Tree is invalid. Fit failed.\n";
240  return;
241  }
242 
243  //accessing the tree components, move pointer to top
244  myTree->movePointerToTheTop();
245 
246  //We are now at the top of the decay tree getting the B_s reconstructed KinematicPartlcle
247  RefCountedKinematicParticle b_s = myTree->currentParticle();
248  printout(b_s);
249 
250  // The B_s decay vertex
251  RefCountedKinematicVertex b_dec_vertex = myTree->currentDecayVertex();
252  printout(b_dec_vertex);
253 
254  // Get all the children of Bs:
255  //In this way, the pointer is not moved
256  vector<RefCountedKinematicParticle> bs_children = myTree->finalStateParticles();
257 
258  for (unsigned int i = 0; i < bs_children.size(); ++i) {
259  printout(bs_children[i]);
260  }
261 
262  //Now navigating down the tree , pointer is moved:
263  bool child = myTree->movePointerToTheFirstChild();
264 
265  if (child)
266  while (myTree->movePointerToTheNextChild()) {
267  RefCountedKinematicParticle aChild = myTree->currentParticle();
268  printout(aChild);
269  }
270 }
271 
272 //Returns the first vertex in the list.
273 
275  // get the simulated vertices
277  iEvent.getByToken(token_VertexTruth, TVCollectionH);
278  const TrackingVertexCollection tPC = *(TVCollectionH.product());
279 
280  // Handle<edm::SimVertexContainer> simVtcs;
281  // iEvent.getByLabel("g4SimHits", simVtcs);
282  // std::cout << "SimVertex " << simVtcs->size() << std::endl;
283  // for(edm::SimVertexContainer::const_iterator v=simVtcs->begin();
284  // v!=simVtcs->end(); ++v){
285  // std::cout << "simvtx "
286  // << v->position().x() << " "
287  // << v->position().y() << " "
288  // << v->position().z() << " "
289  // << v->parentIndex() << " "
290  // << v->noParent() << " "
291  // << std::endl;
292  // }
293  return *(tPC.begin());
294 }
KineExample::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: KineExample.cc:62
KalmanVertexFitter::vertex
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const override
Definition: KalmanVertexFitter.h:49
Handle.h
mps_fire.i
i
Definition: mps_fire.py:428
MessageLogger.h
edm::Handle::product
T const * product() const
Definition: Handle.h:70
MultiTrackKinematicConstraint.h
KalmanVertexFitter.h
ESHandle.h
KineExample::beginRun
void beginRun(edm::Run const &, edm::EventSetup const &) override
Definition: KineExample.cc:50
TransientVertex::isValid
bool isValid() const
Definition: TransientVertex.h:195
KinematicParticleFactoryFromTransientTrack
Definition: KinematicParticleFactoryFromTransientTrack.h:16
ParticleMass
double ParticleMass
Definition: ParticleMass.h:4
edm::Run
Definition: Run.h:45
edm
HLT enums.
Definition: AlignableModifier.h:19
gather_cfg.cout
cout
Definition: gather_cfg.py:144
TransientVertex::position
GlobalPoint position() const
Definition: TransientVertex.h:169
MultiTrackKinematicConstraint
Definition: MultiTrackKinematicConstraint.h:23
TrackingVertexCollection
std::vector< TrackingVertex > TrackingVertexCollection
Definition: TrackingVertexContainer.h:8
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
MassKinematicConstraint
Definition: MassKinematicConstraint.h:17
ReferenceCountingPointer< KinematicTree >
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
TransientTrack.h
edm::Handle< reco::TrackCollection >
KinematicParticleVertexFitter.h
singleTopDQM_cfi.setup
setup
Definition: singleTopDQM_cfi.py:37
TwoTrackMassKinematicConstraint
Definition: TwoTrackMassKinematicConstraint.h:21
KinematicConstrainedVertexFitter::fit
RefCountedKinematicTree fit(const std::vector< RefCountedKinematicParticle > &part)
Definition: KinematicConstrainedVertexFitter.h:46
KineExample
Definition: KineExample.h:38
MakerMacros.h
Point
math::XYZPoint Point
Definition: TrackerDpgAnalysis.cc:106
Track.h
edm::EventSetup::get
T get() const
Definition: EventSetup.h:80
TrackFwd.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
TransientTrackRecord
Definition: TransientTrackRecord.h:11
edm::ESHandle< TransientTrackBuilder >
KineExample::kvfPSet
edm::ParameterSet kvfPSet
Definition: KineExample.h:56
KinematicConstraint
Definition: KinematicConstraint.h:21
KineExample::KineExample
KineExample(const edm::ParameterSet &)
Definition: KineExample.cc:36
KinematicParticleFactoryFromTransientTrack::particle
RefCountedKinematicParticle particle(const reco::TransientTrack &initialTrack, const ParticleMass &massGuess, float chiSquared, float degreesOfFr, float &m_sigma) const
Definition: KinematicParticleFactoryFromTransientTrack.cc:15
KinematicParticleVertexFitter::fit
RefCountedKinematicTree fit(const std::vector< RefCountedKinematicParticle > &particles) const
Definition: KinematicParticleVertexFitter.cc:47
KineExample.h
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
Vertex.h
KineExample::getSimVertex
TrackingVertex getSimVertex(const edm::Event &iEvent) const
Definition: KineExample.cc:274
KineExample::endJob
void endJob() override
Definition: KineExample.cc:56
cppFunctionSkipper.exception
exception
Definition: cppFunctionSkipper.py:10
TransientTrackBuilder.h
edm::ParameterSet
Definition: ParameterSet.h:47
KinematicConstrainedVertexFitter
Definition: KinematicConstrainedVertexFitter.h:21
TrackingVertex
Definition: TrackingVertex.h:22
KineExample::token_VertexTruth
edm::EDGetTokenT< TrackingVertexCollection > token_VertexTruth
Definition: KineExample.h:63
KineExample::printout
void printout(const RefCountedKinematicVertex &myVertex) const
Definition: KineExample.cc:218
iEvent
int iEvent
Definition: GenABIO.cc:224
jpsi
#define jpsi
Definition: CascadeWrapper.h:46
TransientVertex
Definition: TransientVertex.h:18
edm::EventSetup
Definition: EventSetup.h:57
TransientTrackRecord.h
submitPVResolutionJobs.err
err
Definition: submitPVResolutionJobs.py:85
get
#define get
VertexFwd.h
TransientVertex.h
KinematicParticleFitter::fit
std::vector< RefCountedKinematicTree > fit(KinematicConstraint *cs, const std::vector< RefCountedKinematicTree > &trees) const
Definition: KinematicParticleFitter.cc:20
std
Definition: JetResolutionObject.h:76
writedatasetfile.run
run
Definition: writedatasetfile.py:27
KinematicConstrainedVertexFitter.h
reco::TransientTrack
Definition: TransientTrack.h:19
KinematicParticleVertexFitter
Definition: KinematicParticleVertexFitter.h:25
KineExample::~KineExample
~KineExample() override
Definition: KineExample.cc:46
MassKinematicConstraint.h
KinematicParticleFactoryFromTransientTrack.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
KineExample::token_tracks
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: KineExample.h:61
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
KinematicParticleFitter.h
edm::Event
Definition: Event.h:73
child
Definition: simpleInheritance.h:11
KineExample::outputFile_
std::string outputFile_
Definition: KineExample.h:60
ParticleMass.h
edm::InputTag
Definition: InputTag.h:15
TwoTrackMassKinematicConstraint.h
KinematicParticleFitter
Definition: KinematicParticleFitter.h:23
KalmanVertexFitter
Definition: KalmanVertexFitter.h:22