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::one::EDAnalyzer< edm::one::WatchRuns > edm::one::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginRun (edm::Run const &, edm::EventSetup const &) override
 
void endRun (edm::Run const &, edm::EventSetup const &) override
 
 KineExample (const edm::ParameterSet &)
 
 ~KineExample () override
 
- Public Member Functions inherited from edm::one::EDAnalyzer< edm::one::WatchRuns >
 EDAnalyzer ()=default
 
 EDAnalyzer (const EDAnalyzer &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
const EDAnalyzeroperator= (const EDAnalyzer &)=delete
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsInputProcessBlocks () const final
 
bool wantsProcessBlocks () const final
 
- Public Member Functions inherited from edm::one::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDAnalyzerBase () 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
 
ESResolverIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESResolverIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex > const & esGetTokenRecordIndicesVector (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::array< std::vector< ModuleDescription const *> *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, 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
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProductResolverIndices 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

const edm::ESGetToken< MagneticField, IdealMagneticFieldRecordestoken_mf
 
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecordestoken_ttk
 
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::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
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)
 
void resetItemsToGetFrom (BranchType iType)
 

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 59 of file KineExample.cc.

Constructor & Destructor Documentation

◆ KineExample()

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

Definition at line 92 of file KineExample.cc.

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

93  : estoken_ttk(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
94  estoken_mf(esConsumes<edm::Transition::BeginRun>()) {
95  token_tracks = consumes<TrackCollection>(iConfig.getParameter<string>("TrackLabel"));
96  outputFile_ = iConfig.getUntrackedParameter<std::string>("outputFile");
97  kvfPSet = iConfig.getParameter<edm::ParameterSet>("KVFParameters");
98  // rootFile_ = TFile::Open(outputFile_.c_str(),"RECREATE");
99  edm::LogInfo("RecoVertex/KineExample") << "Initializing KVF TEST analyser - Output file: " << outputFile_ << "\n";
100  // token_TrackTruth = consumes<TrackingParticleCollection>(edm::InputTag("trackingtruth", "TrackTruth"));
101  token_VertexTruth = consumes<TrackingVertexCollection>(edm::InputTag("trackingtruth", "VertexTruth"));
102 }
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: KineExample.cc:83
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::ParameterSet kvfPSet
Definition: KineExample.cc:79
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > estoken_mf
Definition: KineExample.cc:76
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< TrackingVertexCollection > token_VertexTruth
Definition: KineExample.cc:85
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > estoken_ttk
Definition: KineExample.cc:75
Log< level::Info, false > LogInfo
std::string outputFile_
Definition: KineExample.cc:82

◆ ~KineExample()

KineExample::~KineExample ( )
overridedefault

Member Function Documentation

◆ analyze()

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

Implements edm::one::EDAnalyzerBase.

Definition at line 115 of file KineExample.cc.

References submitPVResolutionJobs::err, estoken_ttk, KinematicConstrainedVertexFitter::fit(), KinematicParticleFitter::fit(), KinematicParticleVertexFitter::fit(), edm::EventSetup::getHandle(), iEvent, edm::HandleBase::isValid(), TransientVertex::isValid(), KinematicParticleFactoryFromTransientTrack::particle(), TransientVertex::position(), printout(), token_tracks, and KalmanVertexFitter::vertex().

115  {
116  try {
117  edm::LogPrint("KineExample") << "Reconstructing event number: " << iEvent.id() << "\n";
118 
119  // get RECO tracks from the event
120  // `tks` can be used as a ptr to a reco::TrackCollection
122  iEvent.getByToken(token_tracks, tks);
123  if (!tks.isValid()) {
124  edm::LogPrint("KineExample") << "Couln't find track collection: " << iEvent.id() << "\n";
125  } else {
126  edm::LogInfo("RecoVertex/KineExample") << "Found: " << (*tks).size() << " reconstructed tracks"
127  << "\n";
128  edm::LogPrint("KineExample") << "got " << (*tks).size() << " tracks " << endl;
129 
130  // Transform Track to TransientTrack
131  //get the builder:
133  //do the conversion:
134  vector<TransientTrack> t_tks = (*theB).build(tks);
135 
136  edm::LogPrint("KineExample") << "Found: " << t_tks.size() << " reconstructed tracks"
137  << "\n";
138 
139  // Do a KindFit, if >= 4 tracks.
140  if (t_tks.size() > 3) {
141  // For a first test, suppose that the first four tracks are the 2 muons,
142  // then the 2 kaons. Since this will not be true, the result of the fit
143  // will not be meaningfull, but at least you will get the idea of how to
144  // do such a fit.
145 
146  //First, to get started, a simple vertex fit:
147 
148  vector<TransientTrack> ttv;
149  ttv.push_back(t_tks[0]);
150  ttv.push_back(t_tks[1]);
151  ttv.push_back(t_tks[2]);
152  ttv.push_back(t_tks[3]);
153  KalmanVertexFitter kvf(false);
154  TransientVertex tv = kvf.vertex(ttv);
155  if (!tv.isValid())
156  edm::LogPrint("KineExample") << "KVF failed\n";
157  else
158  edm::LogPrint("KineExample") << "KVF fit Position: " << Vertex::Point(tv.position()) << "\n";
159 
160  TransientTrack ttMuPlus = t_tks[0];
161  TransientTrack ttMuMinus = t_tks[1];
162  TransientTrack ttKPlus = t_tks[2];
163  TransientTrack ttKMinus = t_tks[3];
164 
165  //the final state muons and kaons from the Bs->J/PsiPhi->mumuKK decay
166  //Creating a KinematicParticleFactory
168 
169  //The mass of a muon and the insignificant mass sigma to avoid singularities in the covariance matrix.
170  ParticleMass muon_mass = 0.1056583;
171  ParticleMass kaon_mass = 0.493677;
172  float muon_sigma = 0.0000001;
173  float kaon_sigma = 0.000016;
174 
175  //initial chi2 and ndf before kinematic fits. The chi2 of the reconstruction is not considered
176  float chi = 0.;
177  float ndf = 0.;
178 
179  //making particles
180  vector<RefCountedKinematicParticle> muonParticles;
181  vector<RefCountedKinematicParticle> phiParticles;
182  vector<RefCountedKinematicParticle> allParticles;
183  muonParticles.push_back(pFactory.particle(ttMuPlus, muon_mass, chi, ndf, muon_sigma));
184  muonParticles.push_back(pFactory.particle(ttMuMinus, muon_mass, chi, ndf, muon_sigma));
185  allParticles.push_back(pFactory.particle(ttMuPlus, muon_mass, chi, ndf, muon_sigma));
186  allParticles.push_back(pFactory.particle(ttMuMinus, muon_mass, chi, ndf, muon_sigma));
187 
188  phiParticles.push_back(pFactory.particle(ttKPlus, kaon_mass, chi, ndf, kaon_sigma));
189  phiParticles.push_back(pFactory.particle(ttKMinus, kaon_mass, chi, ndf, kaon_sigma));
190  allParticles.push_back(pFactory.particle(ttKPlus, kaon_mass, chi, ndf, kaon_sigma));
191  allParticles.push_back(pFactory.particle(ttKMinus, kaon_mass, chi, ndf, kaon_sigma));
192 
193  /* Example of a simple vertex fit, without other constraints
194  * The reconstructed decay tree is a result of the kinematic fit
195  * The KinematicParticleVertexFitter fits the final state particles to their vertex and
196  * reconstructs the decayed state
197  */
199  edm::LogPrint("KineExample") << "Simple vertex fit with KinematicParticleVertexFitter:\n";
200  RefCountedKinematicTree vertexFitTree = fitter.fit(allParticles);
201 
202  printout(vertexFitTree);
203 
205 
206  //creating the constraint for the J/Psi mass
207  ParticleMass jpsi = 3.09687;
208 
209  //creating the two track mass constraint
211 
212  //creating the fitter
214 
215  //obtaining the resulting tree
216  RefCountedKinematicTree myTree = kcvFitter.fit(allParticles, j_psi_c);
217 
218  edm::LogPrint("KineExample") << "\nGlobal fit done:\n";
219  printout(myTree);
220 
221  //creating the vertex fitter
223 
224  //reconstructing a J/Psi decay
225  RefCountedKinematicTree jpTree = kpvFitter.fit(muonParticles);
226 
227  //creating the particle fitter
228  KinematicParticleFitter csFitter;
229 
230  // creating the constraint
231  float jp_m_sigma = 0.00004;
232  KinematicConstraint* jpsi_c2 = new MassKinematicConstraint(jpsi, jp_m_sigma);
233 
234  //the constrained fit:
235  jpTree = csFitter.fit(jpsi_c2, jpTree);
236 
237  //getting the J/Psi KinematicParticle and putting it together with the kaons.
238  //The J/Psi KinematicParticle has a pointer to the tree it belongs to
239  jpTree->movePointerToTheTop();
240  RefCountedKinematicParticle jpsi_part = jpTree->currentParticle();
241  phiParticles.push_back(jpsi_part);
242 
243  //making a vertex fit and thus reconstructing the Bs parameters
244  // the resulting tree includes all the final state tracks, the J/Psi meson,
245  // its decay vertex, the Bs meson and its decay vertex.
246  RefCountedKinematicTree bsTree = kpvFitter.fit(phiParticles);
247  edm::LogPrint("KineExample") << "Sequential fit done:\n";
248  printout(bsTree);
249 
250  // // For the analysis: compare to your SimVertex
251  // TrackingVertex sv = getSimVertex(iEvent);
252  // edm::Handle<TrackingParticleCollection> TPCollectionH ;
253  // iEvent.getByToken(token_TrackTruth, TPCollectionH);
254  // const TrackingParticleCollection tPC = *(TPCollectionH.product());
255  // reco::RecoToSimCollection recSimColl=associatorForParamAtPca->associateRecoToSim(tks,
256  // TPCollectionH,
257  // &iEvent,
258  // &iSetup);
259  //
260  // tree->fill(tv, &sv, &recSimColl);
261  // }
262  }
263  }
264  } catch (cms::Exception& err) {
265  edm::LogError("KineExample") << "Exception during event number: " << iEvent.id() << "\n" << err.what() << "\n";
266  }
267 }
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: KineExample.cc:83
RefCountedKinematicTree fit(const std::vector< RefCountedKinematicParticle > &part)
GlobalPoint position() const
RefCountedKinematicTree fit(const std::vector< RefCountedKinematicParticle > &particles) const
double ParticleMass
Definition: ParticleMass.h:4
Log< level::Error, false > LogError
int iEvent
Definition: GenABIO.cc:224
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > estoken_ttk
Definition: KineExample.cc:75
void printout(const RefCountedKinematicVertex &myVertex) const
Definition: KineExample.cc:269
bool isValid() const
math::XYZPoint Point
RefCountedKinematicParticle particle(const reco::TransientTrack &initialTrack, const ParticleMass &massGuess, float chiSquared, float degreesOfFr, float &m_sigma) const
Log< level::Warning, true > LogPrint
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
Log< level::Info, false > LogInfo
std::vector< RefCountedKinematicTree > fit(KinematicConstraint *cs, const std::vector< RefCountedKinematicTree > &trees) const
bool isValid() const
Definition: HandleBase.h:70

◆ beginRun()

void KineExample::beginRun ( edm::Run const &  run,
edm::EventSetup const &  setup 
)
override

Definition at line 106 of file KineExample.cc.

106  {
107  //const MagneticField* magField = &setup.getData(estoken_mf);
108  //tree.reset(new SimpleVertexTree("VertexFitter", magField.product()));
109 }

◆ endRun()

void KineExample::endRun ( edm::Run const &  ,
edm::EventSetup const &   
)
inlineoverride

Definition at line 66 of file KineExample.cc.

66 {};

◆ getSimVertex()

TrackingVertex KineExample::getSimVertex ( const edm::Event iEvent) const
private

Definition at line 326 of file KineExample.cc.

References iEvent, edm::Handle< T >::product(), and token_VertexTruth.

326  {
327  // get the simulated vertices
329  iEvent.getByToken(token_VertexTruth, TVCollectionH);
330  const TrackingVertexCollection tPC = *(TVCollectionH.product());
331 
332  // Handle<edm::SimVertexContainer> simVtcs;
333  // iEvent.getByLabel("g4SimHits", simVtcs);
334  // edm::LogPrint("KineExample") << "SimVertex " << simVtcs->size() << std::endl;
335  // for(edm::SimVertexContainer::const_iterator v=simVtcs->begin();
336  // v!=simVtcs->end(); ++v){
337  // edm::LogPrint("KineExample") << "simvtx "
338  // << v->position().x() << " "
339  // << v->position().y() << " "
340  // << v->position().z() << " "
341  // << v->parentIndex() << " "
342  // << v->noParent() << " "
343  // << std::endl;
344  // }
345  return *(tPC.begin());
346 }
T const * product() const
Definition: Handle.h:70
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< TrackingVertexCollection > token_VertexTruth
Definition: KineExample.cc:85
std::vector< TrackingVertex > TrackingVertexCollection

◆ printout() [1/3]

void KineExample::printout ( const RefCountedKinematicVertex myVertex) const
private

Definition at line 269 of file KineExample.cc.

Referenced by analyze(), and printout().

269  {
270  if (myVertex->vertexIsValid()) {
271  edm::LogPrint("KineExample") << "Decay vertex: " << myVertex->position() << myVertex->chiSquared() << " "
272  << myVertex->degreesOfFreedom() << endl;
273  } else
274  edm::LogPrint("KineExample") << "Decay vertex Not valid\n";
275 }
Log< level::Warning, true > LogPrint

◆ printout() [2/3]

void KineExample::printout ( const RefCountedKinematicParticle myParticle) const
private

Definition at line 277 of file KineExample.cc.

277  {
278  edm::LogPrint("KineExample") << "Particle: \n";
279  //accessing the reconstructed Bs meson parameters:
280  //SK: uncomment if needed AlgebraicVector7 bs_par = myParticle->currentState().kinematicParameters().vector();
281 
282  //and their joint covariance matrix:
283  //SK:uncomment if needed AlgebraicSymMatrix77 bs_er = myParticle->currentState().kinematicParametersError().matrix();
284  edm::LogPrint("KineExample") << "Momentum at vertex: " << myParticle->currentState().globalMomentum() << endl;
285  edm::LogPrint("KineExample") << "Parameters at vertex: " << myParticle->currentState().kinematicParameters().vector()
286  << endl;
287 }
Log< level::Warning, true > LogPrint

◆ printout() [3/3]

void KineExample::printout ( const RefCountedKinematicTree myTree) const
private

Definition at line 289 of file KineExample.cc.

References mps_fire::i, and printout().

289  {
290  if (!myTree->isValid()) {
291  edm::LogPrint("KineExample") << "Tree is invalid. Fit failed.\n";
292  return;
293  }
294 
295  //accessing the tree components, move pointer to top
296  myTree->movePointerToTheTop();
297 
298  //We are now at the top of the decay tree getting the B_s reconstructed KinematicPartlcle
299  RefCountedKinematicParticle b_s = myTree->currentParticle();
300  printout(b_s);
301 
302  // The B_s decay vertex
303  RefCountedKinematicVertex b_dec_vertex = myTree->currentDecayVertex();
304  printout(b_dec_vertex);
305 
306  // Get all the children of Bs:
307  //In this way, the pointer is not moved
308  vector<RefCountedKinematicParticle> bs_children = myTree->finalStateParticles();
309 
310  for (unsigned int i = 0; i < bs_children.size(); ++i) {
311  printout(bs_children[i]);
312  }
313 
314  //Now navigating down the tree , pointer is moved:
315  bool child = myTree->movePointerToTheFirstChild();
316 
317  if (child)
318  while (myTree->movePointerToTheNextChild()) {
319  RefCountedKinematicParticle aChild = myTree->currentParticle();
320  printout(aChild);
321  }
322 }
void printout(const RefCountedKinematicVertex &myVertex) const
Definition: KineExample.cc:269
Log< level::Warning, true > LogPrint

Member Data Documentation

◆ estoken_mf

const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> KineExample::estoken_mf
private

Definition at line 76 of file KineExample.cc.

◆ estoken_ttk

const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> KineExample::estoken_ttk
private

Definition at line 75 of file KineExample.cc.

Referenced by analyze().

◆ kvfPSet

edm::ParameterSet KineExample::kvfPSet
private

Definition at line 79 of file KineExample.cc.

Referenced by KineExample().

◆ outputFile_

std::string KineExample::outputFile_
private

Definition at line 82 of file KineExample.cc.

Referenced by KineExample().

◆ theConfig

edm::ParameterSet KineExample::theConfig
private

Definition at line 78 of file KineExample.cc.

◆ token_tracks

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

Definition at line 83 of file KineExample.cc.

Referenced by analyze(), and KineExample().

◆ token_VertexTruth

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

Definition at line 85 of file KineExample.cc.

Referenced by getSimVertex(), and KineExample().