ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitEdm/Analyzers/src/ConversionRejectionExample.cc
Revision: 1.1
Committed: Tue Jun 8 23:30:47 2010 UTC (14 years, 11 months ago) by bendavid
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_018, Mit_017, Mit_017pre3, Mit_017pre2, Mit_017pre1, V07-05-00, Mit_016, Mit_015b, Mit_015a, Mit_015, Mit_014e, Mit_014d, Mit_014c, Mit_014b, ConvRejection-10-06-09
Log Message:
Add example analyzer

File Contents

# Content
1 //--------------------------------------------------------------------------------------------------
2 // $Id: None.cc,v 1.3 2009/12/02 20:39:15 loizides Exp $
3 //
4 // ConversionRejectionExample
5 //
6 // Example analyzer for using the ConversionMatcher tool
7 //
8 // Authors: J.Bendavid
9 //--------------------------------------------------------------------------------------------------
10
11 #include <TMath.h>
12 #include "FWCore/Framework/interface/Frameworkfwd.h"
13 #include "FWCore/Framework/interface/EDAnalyzer.h"
14 #include "FWCore/Framework/interface/Event.h"
15 #include "FWCore/Framework/interface/MakerMacros.h"
16 #include "FWCore/ParameterSet/interface/ParameterSet.h"
17 #include "DataFormats/Common/interface/Handle.h"
18 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
19 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
20 #include "MitEdm/ConversionRejection/interface/ConversionMatcher.h"
21
22 namespace mitedm
23 {
24 class ConversionRejectionExample : public edm::EDAnalyzer {
25 public:
26 explicit ConversionRejectionExample(const edm::ParameterSet&) {}
27 ~ConversionRejectionExample() {};
28
29 private:
30 virtual void beginJob() {}
31 virtual void analyze(const edm::Event &e, const edm::EventSetup &es);
32 virtual void endJob() {}
33 };
34 }
35
36 using namespace mitedm;
37 using namespace edm;
38
39 //define this as a plug-in
40 DEFINE_FWK_MODULE(ConversionRejectionExample);
41
42 void ConversionRejectionExample::analyze(const edm::Event &e, const edm::EventSetup &es)
43 {
44
45 //get electron collection
46 Handle<reco::GsfElectronCollection> hElectrons;
47 e.getByLabel("gsfElectrons", hElectrons);
48 const reco::GsfElectronCollection *electronCol = hElectrons.product();
49
50 //get collection of reconstructed conversions
51 edm::Handle<std::vector<mitedm::DecayPart> > hConversions;
52 e.getByLabel("mvfConversionRemoval", hConversions);
53
54 //initialize ConversionMatcher with default cuts
55 mitedm::ConversionMatcher convMatcher;
56
57 for (reco::GsfElectronCollection::const_iterator it = electronCol->begin(); it!=electronCol->end(); ++it) {
58
59 //check if electron matches a conversion passing all of the cuts (then the electron should be rejected)
60 bool matchesGoodConversion = convMatcher.matchesGoodConversion(*it,hConversions);
61 printf("matchesGoodConversion = %i\n", matchesGoodConversion);
62
63 //dump some info on the good conversions
64 std::vector<edm::Ptr<DecayPart> > goodConversions = convMatcher.goodMatchedConversions(*it,hConversions);
65 for (std::vector<edm::Ptr<DecayPart> >::const_iterator jt = goodConversions.begin(); jt!=goodConversions.end(); ++jt) {
66 const DecayPart *conv = jt->get();
67 printf("Good Matched Conversion %i:\n",jt-goodConversions.begin());
68 printf(" radius = %5f\n",conv->position().rho());
69 printf(" lxy = %5f\n",conv->lxy());
70 printf(" lz = %5f\n",conv->lz());
71 printf(" chi2 = %5f\n",conv->chi2());
72 printf(" ndof = %i\n",conv->ndof());
73 printf(" prob = %5f\n",TMath::Prob(conv->chi2(),conv->ndof()));
74
75 //loop through daughters
76 for (int i=0; i<conv->nStableChild(); ++i) {
77 printf(" Daughter %i\n",i);
78 const StableData &sd = conv->getStableData(i);
79 const StablePart *sp = dynamic_cast<const StablePart*>(sd.originalPtr().get());
80 const reco::Track *trk = sp->track();
81
82 printf(" Hits before vertex = %i\n", sd.nWrongHits());
83 printf(" Track Algo = %i\n", trk->algo());
84 printf(" Track chi2 = %5f\n",trk->chi2());
85 printf(" Track ndof = %5f\n", trk->ndof());
86 printf(" Track prob = %5f\n", TMath::Prob(trk->chi2(),trk->ndof()));
87 }
88 }
89
90 //dump some info on all matched conversions
91 std::vector<edm::Ptr<DecayPart> > allConversions = convMatcher.allMatchedConversions(*it,hConversions);
92 for (std::vector<edm::Ptr<DecayPart> >::const_iterator jt = allConversions.begin(); jt!=allConversions.end(); ++jt) {
93 const DecayPart *conv = jt->get();
94 printf("Matched Conversion %i:\n",jt-allConversions.begin());
95 printf(" radius = %5f\n",conv->position().rho());
96 printf(" prob = %5f\n",TMath::Prob(conv->chi2(),conv->ndof()));
97 }
98 }
99
100 }