ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/Betchart/TopRefTuple/interface/Tuple_Gen.h
Revision: 1.2
Committed: Thu Nov 8 01:47:44 2012 UTC (12 years, 5 months ago) by bbetchar
Content type: text/plain
Branch: MAIN
Changes since 1.1: +8 -71 lines
Log Message:
add lumiTree; add gen

File Contents

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