CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
KineExample Class Reference

#include <RecoVertex/KineExample/src/KineExample.cc>

Inheritance diagram for KineExample:
edm::EDAnalyzer

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 
virtual void beginRun (edm::Run const &, edm::EventSetup const &)
 
virtual void endJob ()
 
 KineExample (const edm::ParameterSet &)
 
 ~KineExample ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Member Functions

TrackingVertex getSimVertex (const edm::Event &iEvent) const
 
void printout (const RefCountedKinematicVertex &myVertex) const
 
void printout (const RefCountedKinematicParticle &myParticle) const
 
void printout (const RefCountedKinematicTree &myTree) const
 

Private Attributes

TrackAssociatorByChi2associatorForParamAtPca
 
edm::ParameterSet kvfPSet
 
std::string outputFile_
 
edm::ParameterSet theConfig
 
std::string trackLabel_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

This is a very simple test analyzer mean to test the KalmanVertexFitter

Description: steers tracker primary vertex reconstruction and storage

Implementation: <Notes on="" implementation>="">

Definition at line 46 of file KineExample.h.

Constructor & Destructor Documentation

KineExample::KineExample ( const edm::ParameterSet iConfig)
explicit

Definition at line 35 of file KineExample.cc.

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), kvfPSet, outputFile_, and trackLabel_.

36  : theConfig(iConfig)
37 {
38  trackLabel_ = iConfig.getParameter<std::string>("TrackLabel");
39  outputFile_ = iConfig.getUntrackedParameter<std::string>("outputFile");
40  kvfPSet = iConfig.getParameter<edm::ParameterSet>("KVFParameters");
41 // rootFile_ = TFile::Open(outputFile_.c_str(),"RECREATE");
42  edm::LogInfo("RecoVertex/KineExample")
43  << "Initializing KVF TEST analyser - Output file: " << outputFile_ <<"\n";
44 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::ParameterSet kvfPSet
Definition: KineExample.h:65
std::string trackLabel_
Definition: KineExample.h:71
edm::ParameterSet theConfig
Definition: KineExample.h:64
std::string outputFile_
Definition: KineExample.h:70
KineExample::~KineExample ( )

Definition at line 47 of file KineExample.cc.

47  {
48 // delete rootFile_;
49 }

Member Function Documentation

void KineExample::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 69 of file KineExample.cc.

References gather_cfg::cout, cppFunctionSkipper::exception, KinematicConstrainedVertexFitter::fit(), KinematicParticleFitter::fit(), KinematicParticleVertexFitter::fit(), edm::EventSetup::get(), edm::Event::getByLabel(), edm::EventBase::id(), edm::HandleBase::isValid(), TransientVertex::isValid(), jpsi, KinematicParticleFactoryFromTransientTrack::particle(), TransientVertex::position(), printout(), trackLabel_, and KalmanVertexFitter::vertex().

70 {
71  try {
72 
73  cout << "Reconstructing event number: " << iEvent.id() << "\n";
74 
75  // get RECO tracks from the event
76  // `tks` can be used as a ptr to a reco::TrackCollection
78  iEvent.getByLabel(trackLabel_, tks);
79  if (!tks.isValid()) {
80  cout
81  << "Couln't find track collection: " << iEvent.id()
82  << "\n";
83  } else {
84 
85  edm::LogInfo("RecoVertex/KineExample")
86  << "Found: " << (*tks).size() << " reconstructed tracks" << "\n";
87  cout << "got " << (*tks).size() << " tracks " << endl;
88 
89  // Transform Track to TransientTrack
90  //get the builder:
92  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
93  //do the conversion:
94  vector<TransientTrack> t_tks = (*theB).build(tks);
95 
96  cout << "Found: " << t_tks.size() << " reconstructed tracks" << "\n";
97 
98  // Do a KindFit, if >= 4 tracks.
99  if (t_tks.size() > 3) {
100 
101  // For a first test, suppose that the first four tracks are the 2 muons,
102  // then the 2 kaons. Since this will not be true, the result of the fit
103  // will not be meaningfull, but at least you will get the idea of how to
104  // do such a fit.
105 
106  //First, to get started, a simple vertex fit:
107 
108  vector<TransientTrack> ttv;
109  ttv.push_back(t_tks[0]); ttv.push_back(t_tks[1]); ttv.push_back(t_tks[2]);ttv.push_back(t_tks[3]);
110  KalmanVertexFitter kvf(false);
111  TransientVertex tv = kvf.vertex(ttv);
112  if (!tv.isValid()) cout << "KVF failed\n";
113  else std::cout << "KVF fit Position: " << Vertex::Point(tv.position()) << "\n";
114 
115 
116  TransientTrack ttMuPlus = t_tks[0];
117  TransientTrack ttMuMinus = t_tks[1];
118  TransientTrack ttKPlus = t_tks[2];
119  TransientTrack ttKMinus = t_tks[3];
120 
121  //the final state muons and kaons from the Bs->J/PsiPhi->mumuKK decay
122  //Creating a KinematicParticleFactory
124 
125  //The mass of a muon and the insignificant mass sigma to avoid singularities in the covariance matrix.
126  ParticleMass muon_mass = 0.1056583;
127  ParticleMass kaon_mass = 0.493677;
128  float muon_sigma = 0.0000001;
129  float kaon_sigma = 0.000016;
130 
131  //initial chi2 and ndf before kinematic fits. The chi2 of the reconstruction is not considered
132  float chi = 0.;
133  float ndf = 0.;
134 
135  //making particles
136  vector<RefCountedKinematicParticle> muonParticles;
137  vector<RefCountedKinematicParticle> phiParticles;
138  vector<RefCountedKinematicParticle> allParticles;
139  muonParticles.push_back(pFactory.particle (ttMuPlus,muon_mass,chi,ndf,muon_sigma));
140  muonParticles.push_back(pFactory.particle (ttMuMinus,muon_mass,chi,ndf,muon_sigma));
141  allParticles.push_back(pFactory.particle (ttMuPlus,muon_mass,chi,ndf,muon_sigma));
142  allParticles.push_back(pFactory.particle (ttMuMinus,muon_mass,chi,ndf,muon_sigma));
143 
144  phiParticles.push_back(pFactory.particle (ttKPlus,kaon_mass,chi,ndf,kaon_sigma));
145  phiParticles.push_back(pFactory.particle (ttKMinus,kaon_mass,chi,ndf,kaon_sigma));
146  allParticles.push_back(pFactory.particle (ttKPlus,kaon_mass,chi,ndf,kaon_sigma));
147  allParticles.push_back(pFactory.particle (ttKMinus,kaon_mass,chi,ndf,kaon_sigma));
148 
149  /* Example of a simple vertex fit, without other constraints
150  * The reconstructed decay tree is a result of the kinematic fit
151  * The KinematicParticleVertexFitter fits the final state particles to their vertex and
152  * reconstructs the decayed state
153  */
155  cout <<"Simple vertex fit with KinematicParticleVertexFitter:\n";
156  RefCountedKinematicTree vertexFitTree = fitter.fit(allParticles);
157 
158  printout(vertexFitTree);
159 
161 
162  //creating the constraint for the J/Psi mass
163  ParticleMass jpsi = 3.09687;
164 
165  //creating the two track mass constraint
167 
168  //creating the fitter
170 
171  //obtaining the resulting tree
172  RefCountedKinematicTree myTree = kcvFitter.fit(allParticles, j_psi_c);
173 
174  cout << "\nGlobal fit done:\n";
175  printout(myTree);
176 
177  //creating the vertex fitter
179 
180  //reconstructing a J/Psi decay
181  RefCountedKinematicTree jpTree = kpvFitter.fit(muonParticles);
182 
183  //creating the particle fitter
184  KinematicParticleFitter csFitter;
185 
186  // creating the constraint
187  float jp_m_sigma = 0.00004;
188  KinematicConstraint * jpsi_c2 = new MassKinematicConstraint(jpsi,jp_m_sigma);
189 
190  //the constrained fit:
191  jpTree = csFitter.fit(jpsi_c2,jpTree);
192 
193  //getting the J/Psi KinematicParticle and putting it together with the kaons.
194  //The J/Psi KinematicParticle has a pointer to the tree it belongs to
195  jpTree->movePointerToTheTop();
196  RefCountedKinematicParticle jpsi_part = jpTree->currentParticle();
197  phiParticles.push_back(jpsi_part);
198 
199  //making a vertex fit and thus reconstructing the Bs parameters
200  // the resulting tree includes all the final state tracks, the J/Psi meson,
201  // its decay vertex, the Bs meson and its decay vertex.
202  RefCountedKinematicTree bsTree = kpvFitter.fit(phiParticles);
203  cout << "Sequential fit done:\n";
204  printout(bsTree);
205 
206 
207 
208 // // For the analysis: compare to your SimVertex
209 // TrackingVertex sv = getSimVertex(iEvent);
210 // edm::Handle<TrackingParticleCollection> TPCollectionH ;
211 // iEvent.getByLabel("trackingtruth","TrackTruth",TPCollectionH);
212 // const TrackingParticleCollection tPC = *(TPCollectionH.product());
213 // reco::RecoToSimCollection recSimColl=associatorForParamAtPca->associateRecoToSim(tks,
214 // TPCollectionH,
215 // &iEvent);
216 //
217 // tree->fill(tv, &sv, &recSimColl);
218 // }
219 
220  }
221  }
222 
223  }
224  catch (std::exception & err) {
225  cout << "Exception during event number: " << iEvent.id()
226  << "\n" << err.what() << "\n";
227  }
228 
229 }
std::vector< RefCountedKinematicTree > fit(KinematicConstraint *cs, std::vector< RefCountedKinematicTree > trees) const
RefCountedKinematicTree fit(std::vector< RefCountedKinematicParticle > part)
void printout(const RefCountedKinematicVertex &myVertex) const
Definition: KineExample.cc:231
double ParticleMass
Definition: ParticleMass.h:5
std::string trackLabel_
Definition: KineExample.h:71
std::pair< double, double > Point
Definition: CaloEllipse.h:18
GlobalPoint position() const
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
RefCountedKinematicTree fit(std::vector< RefCountedKinematicParticle > particles) const
const T & get() const
Definition: EventSetup.h:55
#define jpsi
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:121
bool isValid() const
void KineExample::beginRun ( edm::Run const &  run,
edm::EventSetup const &  setup 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 51 of file KineExample.cc.

References associatorForParamAtPca, edm::EventSetup::get(), and edm::ESHandle< class >::product().

51  {
52  edm::ESHandle<TrackAssociatorBase> theAssociatorForParamAtPca;
53  setup.get<TrackAssociatorRecord>().get("TrackAssociatorByChi2",theAssociatorForParamAtPca);
54  associatorForParamAtPca = (TrackAssociatorByChi2 *) theAssociatorForParamAtPca.product();
55 
56 // tree = new SimpleVertexTree("VertexFitter", associatorForParamAtPca);
57 }
TrackAssociatorByChi2 * associatorForParamAtPca
Definition: KineExample.h:66
T const * product() const
Definition: ESHandle.h:62
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
void KineExample::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 60 of file KineExample.cc.

60  {
61 // delete tree;
62 }
TrackingVertex KineExample::getSimVertex ( const edm::Event iEvent) const
private

Definition at line 287 of file KineExample.cc.

References edm::Event::getByLabel(), and edm::Handle< T >::product().

288 {
289  // get the simulated vertices
291  iEvent.getByLabel("trackingtruth","VertexTruth",TVCollectionH);
292  const TrackingVertexCollection tPC = *(TVCollectionH.product());
293 
294 // Handle<edm::SimVertexContainer> simVtcs;
295 // iEvent.getByLabel("g4SimHits", simVtcs);
296 // std::cout << "SimVertex " << simVtcs->size() << std::endl;
297 // for(edm::SimVertexContainer::const_iterator v=simVtcs->begin();
298 // v!=simVtcs->end(); ++v){
299 // std::cout << "simvtx "
300 // << v->position().x() << " "
301 // << v->position().y() << " "
302 // << v->position().z() << " "
303 // << v->parentIndex() << " "
304 // << v->noParent() << " "
305 // << std::endl;
306 // }
307  return *(tPC.begin());
308 }
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
std::vector< TrackingVertex > TrackingVertexCollection
T const * product() const
Definition: Handle.h:74
void KineExample::printout ( const RefCountedKinematicVertex myVertex) const
private

Definition at line 231 of file KineExample.cc.

References gather_cfg::cout.

Referenced by analyze(), and printout().

232 {
233  if (myVertex->vertexIsValid()) {
234  cout << "Decay vertex: " << myVertex->position() <<myVertex->chiSquared()<< " "<<myVertex->degreesOfFreedom()<<endl;
235  } else cout << "Decay vertex Not valid\n";
236 }
tuple cout
Definition: gather_cfg.py:121
void KineExample::printout ( const RefCountedKinematicParticle myParticle) const
private

Definition at line 238 of file KineExample.cc.

References gather_cfg::cout.

239 {
240  cout << "Particle: \n";
241 //accessing the reconstructed Bs meson parameters:
242 //SK: uncomment if needed AlgebraicVector7 bs_par = myParticle->currentState().kinematicParameters().vector();
243 
244 //and their joint covariance matrix:
245 //SK:uncomment if needed AlgebraicSymMatrix77 bs_er = myParticle->currentState().kinematicParametersError().matrix();
246  cout << "Momentum at vertex: " << myParticle->currentState().globalMomentum ()<<endl;
247  cout << "Parameters at vertex: " << myParticle->currentState().kinematicParameters().vector()<<endl;
248 }
tuple cout
Definition: gather_cfg.py:121
void KineExample::printout ( const RefCountedKinematicTree myTree) const
private

Definition at line 250 of file KineExample.cc.

References gather_cfg::cout, i, and printout().

251 {
252  if (!myTree->isValid()) {
253  cout <<"Tree is invalid. Fit failed.\n";
254  return;
255  }
256 
257 //accessing the tree components, move pointer to top
258  myTree->movePointerToTheTop();
259 
260 //We are now at the top of the decay tree getting the B_s reconstructed KinematicPartlcle
261  RefCountedKinematicParticle b_s = myTree->currentParticle();
262  printout(b_s);
263 
264 // The B_s decay vertex
265  RefCountedKinematicVertex b_dec_vertex = myTree->currentDecayVertex();
266  printout(b_dec_vertex);
267 
268  // Get all the children of Bs:
269  //In this way, the pointer is not moved
270  vector< RefCountedKinematicParticle > bs_children = myTree->finalStateParticles();
271 
272  for (unsigned int i=0;i< bs_children.size();++i) {
273  printout(bs_children[i]);
274  }
275 
276 //Now navigating down the tree , pointer is moved:
277  bool child = myTree->movePointerToTheFirstChild();
278 
279  if(child) while (myTree->movePointerToTheNextChild()) {
280  RefCountedKinematicParticle aChild = myTree->currentParticle();
281  printout(aChild);
282  }
283 }
int i
Definition: DBlmapReader.cc:9
void printout(const RefCountedKinematicVertex &myVertex) const
Definition: KineExample.cc:231
tuple cout
Definition: gather_cfg.py:121

Member Data Documentation

TrackAssociatorByChi2* KineExample::associatorForParamAtPca
private

Definition at line 66 of file KineExample.h.

Referenced by beginRun().

edm::ParameterSet KineExample::kvfPSet
private

Definition at line 65 of file KineExample.h.

Referenced by KineExample().

std::string KineExample::outputFile_
private

Definition at line 70 of file KineExample.h.

Referenced by KineExample().

edm::ParameterSet KineExample::theConfig
private

Definition at line 64 of file KineExample.h.

std::string KineExample::trackLabel_
private

Definition at line 71 of file KineExample.h.

Referenced by analyze(), and KineExample().