ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitProd/TreeFiller/src/FillerMCEventInfo.cc
Revision: 1.5
Committed: Sun Mar 22 10:00:29 2009 UTC (16 years, 1 month ago) by loizides
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_009a, Mit_009, Mit_008
Changes since 1.4: +31 -27 lines
Log Message:
Fix printout.

File Contents

# User Rev Content
1 loizides 1.5 // $Id: FillerMCEventInfo.cc,v 1.4 2009/03/19 17:28:50 loizides Exp $
2 loizides 1.1
3     #include "MitProd/TreeFiller/interface/FillerMCEventInfo.h"
4     #include "FWCore/MessageLogger/interface/MessageLogger.h"
5     #include "DataFormats/Common/interface/Handle.h"
6 loizides 1.3 #include "DataFormats/HepMCCandidate/interface/PdfInfo.h"
7 loizides 1.4 #include "SimDataFormats/HepMCProduct/interface/HepMCProduct.h"
8 loizides 1.1 #include "MitAna/DataTree/interface/Names.h"
9     #include "MitAna/DataTree/interface/MCEventInfo.h"
10    
11     using namespace std;
12     using namespace edm;
13     using namespace mithep;
14    
15     //--------------------------------------------------------------------------------------------------
16 loizides 1.3 FillerMCEventInfo::FillerMCEventInfo(const ParameterSet &cfg, const char *name, bool active) :
17 loizides 1.1 BaseFiller(cfg,"MCEventInfo",active),
18     evtName_(Conf().getUntrackedParameter<string>("evtName",Names::gkMCEvtInfoBrn)),
19 loizides 1.4 genHepMCEvName_(Conf().getUntrackedParameter<string>("genHepMCEventEdmName","source")),
20 loizides 1.3 genEvWeightName_(Conf().getUntrackedParameter<string>("genEventWeightEdmName","genEventWeight")),
21     genEvScaleName_(Conf().getUntrackedParameter<string>("genEventScaleEdmName","genEventScale")),
22     genEvProcIdName_(Conf().getUntrackedParameter<string>("genEventProcIdEdmName","genEventProcID")),
23     genPdfInfoName_(Conf().getUntrackedParameter<string>("genPdfInfoEdmName","genEventPdfInfo")),
24 loizides 1.1 eventInfo_(new MCEventInfo())
25     {
26     // Constructor.
27     }
28    
29     //--------------------------------------------------------------------------------------------------
30     FillerMCEventInfo::~FillerMCEventInfo()
31     {
32     // Destructor.
33    
34     delete eventInfo_;
35     }
36    
37     //--------------------------------------------------------------------------------------------------
38     void FillerMCEventInfo::BookDataBlock(TreeWriter &tws)
39     {
40     // Create run info tre and book our branches.
41    
42 loizides 1.2 tws.AddBranch(evtName_,&eventInfo_);
43     OS()->add<mithep::MCEventInfo>(eventInfo_,evtName_);
44 loizides 1.1 }
45    
46     //--------------------------------------------------------------------------------------------------
47     void FillerMCEventInfo::FillDataBlock(const edm::Event &event,
48     const edm::EventSetup &setup)
49     {
50     // Fill our data structures.
51    
52 loizides 1.3 if (event.isRealData()) {
53     PrintErrorAndExit("Expected monte-carlo record, but did not get it. Aborting.");
54     }
55    
56 loizides 1.4 Handle<edm::HepMCProduct> hHepMCProduct;
57 loizides 1.5 if (genHepMCEvName_.empty() || !GetProductSafe(genHepMCEvName_, hHepMCProduct, event)) {
58 loizides 1.4
59 loizides 1.5 if (!genEvWeightName_.empty()) {
60     Handle<double> genEventWeight;
61     GetProduct(genEvWeightName_, genEventWeight, event);
62     eventInfo_->SetWeight(*genEventWeight);
63     }
64    
65     if (!genEvScaleName_.empty()) {
66     Handle<double> genEventScale;
67     GetProduct(genEvScaleName_, genEventScale, event);
68     eventInfo_->SetScale(*genEventScale);
69     }
70    
71     if (!genEvProcIdName_.empty()) {
72     Handle<int> genEventProcId;
73     GetProduct(genEvProcIdName_, genEventProcId, event);
74     eventInfo_->SetProcessId(*genEventProcId);
75     }
76    
77     if (!genPdfInfoName_.empty()) {
78     Handle<reco::PdfInfo> genPdfInfo;
79     GetProduct(genPdfInfoName_, genPdfInfo, event);
80     eventInfo_->SetId1(genPdfInfo->id1);
81     eventInfo_->SetId2(genPdfInfo->id2);
82     eventInfo_->SetPdf1(genPdfInfo->pdf1);
83     eventInfo_->SetPdf2(genPdfInfo->pdf2);
84     eventInfo_->SetScalePdf(genPdfInfo->scalePDF);
85     eventInfo_->SetX1(genPdfInfo->pdf1);
86     eventInfo_->SetX2(genPdfInfo->pdf2);
87     }
88 loizides 1.4 } else {
89     const HepMC::GenEvent *genEvt = hHepMCProduct->GetEvent();
90     eventInfo_->SetScale(genEvt->event_scale());
91     eventInfo_->SetProcessId(genEvt->signal_process_id());
92     HepMC::WeightContainer wc = genEvt->weights();
93     if (wc.size() > 0)
94     eventInfo_->SetWeight(wc[0]);
95     const HepMC::PdfInfo *genPdfInfo = genEvt->pdf_info();
96     eventInfo_->SetId1(genPdfInfo->id1());
97     eventInfo_->SetId2(genPdfInfo->id2());
98     eventInfo_->SetPdf1(genPdfInfo->pdf1());
99     eventInfo_->SetPdf2(genPdfInfo->pdf2());
100     eventInfo_->SetScalePdf(genPdfInfo->scalePDF());
101     eventInfo_->SetX1(genPdfInfo->pdf1());
102     eventInfo_->SetX2(genPdfInfo->pdf2());
103     }
104 loizides 1.1 }