CMS 3D CMS Logo

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