Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:32:12

0001 #include "Validation/EventGenerator/interface/WeightManager.h"
0002 #include "FWCore/Framework/interface/Event.h"
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004 
0005 #include "DataFormats/Common/interface/Handle.h"
0006 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0007 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
0008 
0009 using namespace edm;
0010 
0011 WeightManager::WeightManager(const ParameterSet& iConfig, edm::ConsumesCollector iC)
0012     : _useHepMC(iConfig.getParameter<bool>("UseWeightFromHepMC")) {
0013   if (_useHepMC) {
0014     _hepmcCollection = iConfig.getParameter<InputTag>("hepmcCollection");
0015     hepmcCollectionToken_ = iC.consumes<HepMCProduct>(_hepmcCollection);
0016   } else {
0017     _genEventInfos = iConfig.getParameter<std::vector<InputTag>>("genEventInfos");
0018     for (unsigned int i = 0; i < _genEventInfos.size(); i++)
0019       genEventInfosTokens_.push_back(iC.consumes<GenEventInfoProduct>(_genEventInfos[i]));
0020   }
0021 }
0022 
0023 double WeightManager::weight(const Event& iEvent) {
0024   if (_useHepMC) {
0025     edm::Handle<HepMCProduct> evt;
0026     iEvent.getByToken(hepmcCollectionToken_, evt);
0027     const HepMC::GenEvent* myGenEvent = evt->GetEvent();
0028 
0029     double weight = 1.;
0030     if (!myGenEvent->weights().empty())
0031       weight = myGenEvent->weights()[0];
0032     return weight;
0033   } else {
0034     double weight = 1.;
0035     for (unsigned int i = 0; i < genEventInfosTokens_.size(); ++i) {
0036       edm::Handle<GenEventInfoProduct> info;
0037       iEvent.getByToken(genEventInfosTokens_[i], info);
0038       weight *= info->weight();
0039     }
0040     return weight;
0041   }
0042 }
0043 
0044 std::vector<std::vector<double>> WeightManager::weightsCollection(const Event& iEvent) {
0045   std::vector<std::vector<double>> weightsCollection;
0046   for (unsigned int i = 0; i < genEventInfosTokens_.size(); ++i) {
0047     edm::Handle<GenEventInfoProduct> info;
0048     iEvent.getByToken(genEventInfosTokens_[i], info);
0049     weightsCollection.push_back(info->weights());
0050   }
0051 
0052   return weightsCollection;
0053 }