ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/Betchart/TopRefTuple/interface/Tuple_Gen.h
Revision: 1.4
Committed: Fri Nov 16 18:42:12 2012 UTC (12 years, 5 months ago) by bbetchar
Content type: text/plain
Branch: MAIN
CVS Tags: V00-03-02, V00-03-01, V00-02-02, V00-02-01, V00-02-00, V00-01-05, V00-01-04, V00-01-03, V00-01-02, V00-01-01, V00-01-00, HEAD
Changes since 1.3: +11 -8 lines
Log Message:
optionally store only status 3 particles

File Contents

# Content
1 #ifndef TUPLE_GEN
2 #define TUPLE_GEN
3
4 #include "FWCore/Framework/interface/EDProducer.h"
5 #include "FWCore/Framework/interface/Frameworkfwd.h"
6 #include "FWCore/Utilities/interface/InputTag.h"
7 #include "FWCore/Framework/interface/Event.h"
8 #include "FWCore/Framework/interface/ESHandle.h"
9 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
10 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
11 #include "fTypes.h"
12
13 #include <map>
14
15 template< typename T >
16 class Tuple_Gen : public edm::EDProducer {
17 public:
18 explicit Tuple_Gen(const edm::ParameterSet&);
19 private:
20 void produce(edm::Event &, const edm::EventSetup & );
21 void produceGenJets(edm::Event &);
22 static int index(const reco::Candidate*, const std::vector<const T*>&);
23 typedef fTypes::dPolarLorentzV LorentzVector;
24 const edm::InputTag inputTag;
25 const std::vector<edm::InputTag> jetCollections;
26 const std::string Prefix,Suffix;
27 const double GenStatus1PtCut;
28 const double GenJetPtCut;
29 const bool onlyStatus3;
30 };
31
32 template< typename T > Tuple_Gen<T>::
33 Tuple_Gen(const edm::ParameterSet& conf)
34 : inputTag(conf.getParameter<edm::InputTag>("InputTag")),
35 jetCollections(conf.getParameter<std::vector<edm::InputTag> >("JetCollections")),
36 Prefix(conf.getParameter<std::string>("Prefix")),
37 Suffix(conf.getParameter<std::string>("Suffix")),
38 GenStatus1PtCut(conf.getParameter<double>("GenStatus1PtCut")),
39 GenJetPtCut(conf.getParameter<double>("GenJetPtCut")),
40 onlyStatus3(conf.getParameter<bool>("OnlyStatus3"))
41 {
42 produces <unsigned int> (Prefix + "signalProcessID" + Suffix);
43 produces <bool> (Prefix + "GenInfoHandleValid" + Suffix);
44 produces <bool > (Prefix + "HandleValid" + Suffix);
45 produces <float> (Prefix + "pthat" + Suffix);
46 produces <int> (Prefix + "id1" + Suffix);
47 produces <int> (Prefix + "id2" + Suffix);
48 produces <float> (Prefix + "x1" + Suffix);
49 produces <float> (Prefix + "x2" + Suffix);
50 produces <float> (Prefix + "pdf1" + Suffix);
51 produces <float> (Prefix + "pdf2" + Suffix);
52 produces <std::vector<float> > (Prefix + "BinningValues" + Suffix);
53 produces <float> (Prefix + "Q" + Suffix);
54 produces <std::vector<LorentzVector> > ( Prefix + "P4" + Suffix );
55 produces <std::vector<int> > (Prefix + "PdgId" + Suffix);
56 produces <std::vector<int> > (Prefix + "Status" + Suffix);
57 produces <std::vector<int> > (Prefix + "MotherIndex" + Suffix);
58 produces <std::vector<int> > (Prefix + "MotherPdgId" + Suffix);
59
60
61 for(unsigned i=0; i<jetCollections.size(); ++i)
62 produces<std::vector<LorentzVector> >(Prefix + jetCollections[i].label() +"P4" + Suffix);
63 }
64
65 template< typename T > int Tuple_Gen<T>::
66 index(const reco::Candidate* item, const typename std::vector<const T*>& collection) {
67 typename std::vector<const T*>::const_iterator it(collection.begin()), begin(collection.begin()), end(collection.end());
68 for(; it!=end; it++) if ((*it)==item) return it-begin; //Compare addresses
69 return -1;
70 }
71
72 template< typename T > void Tuple_Gen<T>::
73 produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
74 produceGenJets(iEvent);
75
76 edm::Handle<std::vector<T> > collection; iEvent.getByLabel(inputTag,collection);
77 edm::Handle<GenEventInfoProduct> geninfo; iEvent.getByLabel("generator",geninfo);
78
79 std::auto_ptr<unsigned int> signalProcessID(new unsigned int(geninfo->signalProcessID()));
80 std::auto_ptr<float> Q (new float(geninfo->pdf()->scalePDF));
81 std::auto_ptr<int> id1 (new int( geninfo->pdf()->id.first));
82 std::auto_ptr<int> id2 (new int( geninfo->pdf()->id.second));
83 std::auto_ptr<float> x1 (new float( geninfo->pdf()->x.first));
84 std::auto_ptr<float> x2 (new float( geninfo->pdf()->x.second));
85 std::auto_ptr<float> pdf1 (new float( geninfo->pdf()->xPDF.first));
86 std::auto_ptr<float> pdf2 (new float( geninfo->pdf()->xPDF.second));
87
88 std::auto_ptr<bool> handleValid ( new bool(collection.isValid()) );
89 std::auto_ptr<bool> genInfoValid ( new bool( geninfo.isValid() && !geninfo->binningValues().empty()));
90 std::auto_ptr<float> pthat (new float(*genInfoValid ? geninfo->binningValues()[0] : -1.));
91 std::auto_ptr<std::vector<float> > binningValues (*genInfoValid ? new std::vector<float>(geninfo->binningValues().begin(),
92 geninfo->binningValues().end()) : new std::vector<float>());
93 std::auto_ptr<std::vector<LorentzVector> > p4 ( new std::vector<LorentzVector>() ) ;
94 std::auto_ptr<std::vector<int> > status ( new std::vector<int>() ) ;
95 std::auto_ptr<std::vector<int> > pdgId ( new std::vector<int>() ) ;
96 std::auto_ptr<std::vector<int> > motherIndex ( new std::vector<int>() ) ;
97 std::auto_ptr<std::vector<int> > motherPdgId ( new std::vector<int>() ) ;
98
99 std::vector<const T*> self;
100 std::vector<const reco::Candidate*> mom;
101
102 if(collection.isValid()){
103 for(typename std::vector<T>::const_iterator it = collection->begin(); it != collection->end(); ++it) {
104 if ( it->status() == 3 || // any status 3 genParticle
105 ( !onlyStatus3 &&
106 ( abs(it->pdgId()) == 11 // any electron
107 || abs(it->pdgId()) == 13 // any muon
108 || abs(it->pdgId()) == 15 // any tau
109 || ( it->status() == 1 // status 1 particles
110 && it->pt() > GenStatus1PtCut) // above threshold
111 ))) {
112 p4->push_back(LorentzVector(it->pt(),it->eta(),it->phi(),it->mass()));
113 status->push_back(it->status());
114 pdgId->push_back(it->pdgId());
115 motherPdgId->push_back( it->numberOfMothers() ? it->mother()->pdgId() : 0 );
116 self.push_back(&*it);
117 mom.push_back( it->numberOfMothers() ? it->mother(): 0);
118 }
119 }
120 } //collection
121
122
123 for(typename std::vector<const reco::Candidate*>::const_iterator it = mom.begin(); it!=mom.end(); ++it)
124 motherIndex->push_back( index(*it,self) );
125
126 iEvent.put( handleValid, Prefix + "HandleValid" + Suffix);
127 iEvent.put( genInfoValid, Prefix + "GenInfoHandleValid" + Suffix);
128 iEvent.put( pthat, Prefix + "pthat" + Suffix);
129 iEvent.put( binningValues,Prefix + "BinningValues" + Suffix);
130 iEvent.put( p4, Prefix + "P4" + Suffix );
131 iEvent.put( status, Prefix + "Status" + Suffix );
132 iEvent.put( pdgId, Prefix + "PdgId" + Suffix );
133 iEvent.put( motherIndex, Prefix + "MotherIndex" + Suffix );
134 iEvent.put( motherPdgId, Prefix + "MotherPdgId" + Suffix );
135 iEvent.put( signalProcessID, Prefix + "signalProcessID" + Suffix );
136 iEvent.put( Q, Prefix + "Q" + Suffix );
137 iEvent.put( x1, Prefix + "x1" + Suffix );
138 iEvent.put( x2, Prefix + "x2" + Suffix );
139 iEvent.put( pdf1, Prefix + "pdf1" + Suffix );
140 iEvent.put( pdf2, Prefix + "pdf2" + Suffix );
141 iEvent.put( id1, Prefix + "id1" + Suffix );
142 iEvent.put( id2, Prefix + "id2" + Suffix );
143 }
144
145 template< typename T > void Tuple_Gen<T>::
146 produceGenJets(edm::Event& iEvent) {
147 for(unsigned i=0; i<jetCollections.size(); ++i) {
148 std::auto_ptr<std::vector<LorentzVector> > p4(new std::vector<LorentzVector>());
149 edm::Handle<edm::View<reco::GenJet> > genjets;
150 iEvent.getByLabel(jetCollections[i], genjets);
151 if(genjets.isValid())
152 for(edm::View<reco::GenJet>::const_iterator it(genjets->begin()), end(genjets->end()); it!=end; ++it) {
153 if (it->pt() >= GenJetPtCut) p4->push_back(LorentzVector(it->pt(),it->eta(),it->phi(),it->mass()));
154 }
155 iEvent.put(p4, Prefix + jetCollections[i].label() + "P4" + Suffix);
156 }
157 }
158
159 #endif