CMS 3D CMS Logo

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 edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginRun (edm::Run const &, edm::EventSetup const &) override
 
void endJob () override
 
 KineExample (const edm::ParameterSet &)
 
 ~KineExample () override
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

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

edm::ParameterSet kvfPSet
 
std::string outputFile_
 
edm::ParameterSet theConfig
 
edm::EDGetTokenT< reco::TrackCollectiontoken_tracks
 
edm::EDGetTokenT< TrackingVertexCollectiontoken_VertexTruth
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

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 39 of file KineExample.h.

Constructor & Destructor Documentation

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

Definition at line 37 of file KineExample.cc.

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), kvfPSet, outputFile_, AlCaHLTBitMon_QueryRunRegistry::string, token_tracks, and token_VertexTruth.

38  : theConfig(iConfig)
39 {
40  token_tracks = consumes<TrackCollection>(iConfig.getParameter<string>("TrackLabel"));
41  outputFile_ = iConfig.getUntrackedParameter<std::string>("outputFile");
42  kvfPSet = iConfig.getParameter<edm::ParameterSet>("KVFParameters");
43 // rootFile_ = TFile::Open(outputFile_.c_str(),"RECREATE");
44  edm::LogInfo("RecoVertex/KineExample")
45  << "Initializing KVF TEST analyser - Output file: " << outputFile_ <<"\n";
46 // token_TrackTruth = consumes<TrackingParticleCollection>(edm::InputTag("trackingtruth", "TrackTruth"));
47  token_VertexTruth = consumes<TrackingVertexCollection>(edm::InputTag("trackingtruth", "VertexTruth"));
48 }
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: KineExample.h:63
T getUntrackedParameter(std::string const &, T const &) const
edm::ParameterSet kvfPSet
Definition: KineExample.h:58
edm::ParameterSet theConfig
Definition: KineExample.h:57
edm::EDGetTokenT< TrackingVertexCollection > token_VertexTruth
Definition: KineExample.h:65
std::string outputFile_
Definition: KineExample.h:62
KineExample::~KineExample ( )
override

Definition at line 51 of file KineExample.cc.

51  {
52 // delete rootFile_;
53 }

Member Function Documentation

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

Definition at line 71 of file KineExample.cc.

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

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

Reimplemented from edm::EDAnalyzer.

Definition at line 55 of file KineExample.cc.

55  {
56 
57 // edm::ESHandle<MagneticField> magField;
58 // setup.get<IdealMagneticFieldRecord>().get(magField);
59 // tree.reset(new SimpleVertexTree("VertexFitter", magField.product()));
60 }
void KineExample::endJob ( void  )
overridevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 63 of file KineExample.cc.

63  {
64 }
TrackingVertex KineExample::getSimVertex ( const edm::Event iEvent) const
private

Definition at line 290 of file KineExample.cc.

References DEFINE_FWK_MODULE, edm::Event::getByToken(), edm::Handle< T >::product(), and token_VertexTruth.

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

Definition at line 234 of file KineExample.cc.

References gather_cfg::cout.

Referenced by analyze(), and printout().

235 {
236  if (myVertex->vertexIsValid()) {
237  cout << "Decay vertex: " << myVertex->position() <<myVertex->chiSquared()<< " "<<myVertex->degreesOfFreedom()<<endl;
238  } else cout << "Decay vertex Not valid\n";
239 }
void KineExample::printout ( const RefCountedKinematicParticle myParticle) const
private

Definition at line 241 of file KineExample.cc.

References gather_cfg::cout.

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

Definition at line 253 of file KineExample.cc.

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

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

Member Data Documentation

edm::ParameterSet KineExample::kvfPSet
private

Definition at line 58 of file KineExample.h.

Referenced by KineExample().

std::string KineExample::outputFile_
private

Definition at line 62 of file KineExample.h.

Referenced by KineExample().

edm::ParameterSet KineExample::theConfig
private

Definition at line 57 of file KineExample.h.

edm::EDGetTokenT<reco::TrackCollection> KineExample::token_tracks
private

Definition at line 63 of file KineExample.h.

Referenced by analyze(), and KineExample().

edm::EDGetTokenT<TrackingVertexCollection> KineExample::token_VertexTruth
private

Definition at line 65 of file KineExample.h.

Referenced by getSimVertex(), and KineExample().