Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Run3-gex133A Try to avoid compilation issues in 2 test codes of GeneratorInterface #37780

Merged
merged 2 commits into from May 5, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
65 changes: 32 additions & 33 deletions GeneratorInterface/Pythia6Interface/test/HZZ4muAnalyzer.cc
Expand Up @@ -13,47 +13,49 @@
#include "CommonTools/UtilAlgos/interface/TFileService.h"
#include "TH1.h"

#include "FWCore/Framework/interface/EDAnalyzer.h"
#include "FWCore/Framework/interface/one/EDAnalyzer.h"

class HZZ4muAnalyzer : public edm::EDAnalyzer {
class HZZ4muAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
public:
//
explicit HZZ4muAnalyzer(const edm::ParameterSet&);
virtual ~HZZ4muAnalyzer() {} // no need to delete ROOT stuff
~HZZ4muAnalyzer() = default; // no need to delete ROOT stuff
// as it'll be deleted upon closing TFile

virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
virtual void beginJob() override;
void analyze(const edm::Event&, const edm::EventSetup&) override;
void beginJob() override;

private:
const edm::EDGetTokenT<GenEventInfoProduct> tokenGenInfo_;
const edm::EDGetTokenT<edm::HepMCProduct> tokenHepMC_;
TH1D* fHist2muMass;
TH1D* fHist4muMass;
TH1D* fHistZZMass;
};

using namespace edm;
using namespace std;

HZZ4muAnalyzer::HZZ4muAnalyzer(const ParameterSet& pset) : fHist2muMass(0), fHist4muMass(0), fHistZZMass(0) {
HZZ4muAnalyzer::HZZ4muAnalyzer(const edm::ParameterSet& pset)
: tokenGenInfo_(consumes<GenEventInfoProduct>(edm::InputTag("generator"))),
tokenHepMC_(consumes<edm::HepMCProduct>(edm::InputTag("VtxSmeared"))),
fHist2muMass(nullptr),
fHist4muMass(nullptr),
fHistZZMass(nullptr) {
// actually, pset is NOT in use - we keep it here just for illustratory putposes
usesResource(TFileService::kSharedResource);
}

void HZZ4muAnalyzer::beginJob() {
Service<TFileService> fs;
edm::Service<TFileService> fs;
fHist2muMass = fs->make<TH1D>("Hist2muMass", "2-mu inv. mass", 100, 60., 120.);
fHist4muMass = fs->make<TH1D>("Hist4muMass", "4-mu inv. mass", 100, 170., 210.);
fHistZZMass = fs->make<TH1D>("HistZZMass", "ZZ inv. mass", 100, 170., 210.);

return;
}

void HZZ4muAnalyzer::analyze(const Event& e, const EventSetup&) {
void HZZ4muAnalyzer::analyze(const edm::Event& e, const edm::EventSetup&) {
// here's an example of accessing GenEventInfoProduct
Handle<GenEventInfoProduct> GenInfoHandle;
e.getByLabel("generator", GenInfoHandle);
const edm::Handle<GenEventInfoProduct>& GenInfoHandle = e.getHandle(tokenGenInfo_);
double qScale = GenInfoHandle->qScale();
double pthat = (GenInfoHandle->hasBinningValues() ? (GenInfoHandle->binningValues())[0] : 0.0);
cout << " qScale = " << qScale << " pthat = " << pthat << endl;
std::cout << " qScale = " << qScale << " pthat = " << pthat << std::endl;
//
// this (commented out) code below just exemplifies how to access certain info
//
Expand All @@ -69,11 +71,8 @@ void HZZ4muAnalyzer::analyze(const Event& e, const EventSetup&) {

// here's an example of accessing particles in the event record (HepMCProduct)
//
Handle<HepMCProduct> EvtHandle;

// find initial (unsmeared, unfiltered,...) HepMCProduct
//
e.getByLabel("VtxSmeared", EvtHandle);
const edm::Handle<edm::HepMCProduct>& EvtHandle = e.getHandle(tokenHepMC_);

const HepMC::GenEvent* Evt = EvtHandle->GetEvent();

Expand Down Expand Up @@ -106,34 +105,34 @@ void HZZ4muAnalyzer::analyze(const Event& e, const EventSetup&) {
}

if (HiggsDecVtx == 0) {
cout << " There is NO Higgs in this event ! " << endl;
std::cout << " There is NO Higgs in this event ! " << std::endl;
return;
}

if (e.id().event() == 1) {
cout << " " << endl;
cout << " We do some example printouts in the event 1 " << endl;
cout << " Higgs decay found at the vertex " << HiggsDecVtx->barcode() << " (barcode)" << endl;
std::cout << " " << std::endl;
std::cout << " We do some example printouts in the event 1 " << std::endl;
std::cout << " Higgs decay found at the vertex " << HiggsDecVtx->barcode() << " (barcode)" << std::endl;

vector<HepMC::GenParticle*> HiggsChildren;
std::vector<HepMC::GenParticle*> HiggsChildren;

for (HepMC::GenVertex::particles_out_const_iterator H0in = HiggsDecVtx->particles_out_const_begin();
H0in != HiggsDecVtx->particles_out_const_end();
H0in++) {
HiggsChildren.push_back(*H0in);
}
cout << " Number of Higgs (immediate) children = " << HiggsChildren.size() << endl;
std::cout << " Number of Higgs (immediate) children = " << HiggsChildren.size() << std::endl;
for (unsigned int ic = 0; ic < HiggsChildren.size(); ic++) {
HiggsChildren[ic]->print();
}
}

// select and store stable descendants of the Higgs
//
vector<HepMC::GenParticle*> StableHiggsDesc;
std::vector<HepMC::GenParticle*> StableHiggsDesc;

if (e.id().event() == 1)
cout << " Now let us list all descendents of the Higgs" << endl;
std::cout << " Now let us list all descendents of the Higgs" << std::endl;
for (HepMC::GenVertex::particle_iterator des = HiggsDecVtx->particles_begin(HepMC::descendants);
des != HiggsDecVtx->particles_end(HepMC::descendants);
des++) {
Expand All @@ -147,21 +146,21 @@ void HZZ4muAnalyzer::analyze(const Event& e, const EventSetup&) {
double XMass2part = 0.;
double XMass4part = 0.;
double XMass2pairs = 0.;
vector<HepMC::FourVector> Mom2partCont;
std::vector<HepMC::FourVector> Mom2partCont;

// browse the array of stable descendants
// and do 2-mu inv.mass
//
for (unsigned int i = 0; i < StableHiggsDesc.size(); i++) {
// skip other than mu
//
if (abs(StableHiggsDesc[i]->pdg_id()) != 13)
if (std::abs(StableHiggsDesc[i]->pdg_id()) != 13)
continue;

for (unsigned int j = i + 1; j < StableHiggsDesc.size(); j++) {
// skip other than mu
//
if (abs(StableHiggsDesc[j]->pdg_id()) != 13)
if (std::abs(StableHiggsDesc[j]->pdg_id()) != 13)
continue;
//
// skip same charge combo's
Expand Down Expand Up @@ -205,11 +204,11 @@ void HZZ4muAnalyzer::analyze(const Event& e, const EventSetup&) {
XMass4part = HepMC::FourVector(px4, py4, pz4, e4).m();
fHist4muMass->Fill(XMass4part);
}
//cout << " 4-part inv. mass = " << XMass4part << endl ;
//std::cout << " 4-part inv. mass = " << XMass4part << std::endl ;

// make 2-pairs (ZZ) inv.mass
//
//cout << " selected Z-candidates in this event : " << Mom2partCont.size() << endl ;
//std::cout << " selected Z-candidates in this event : " << Mom2partCont.size() << std::endl ;
for (unsigned int i = 0; i < Mom2partCont.size(); i++) {
for (unsigned int j = i + 1; j < Mom2partCont.size(); j++) {
// Mom2pairs = Mom2partCont[i] + Mom2partCont[j] ;
Expand Down
37 changes: 15 additions & 22 deletions GeneratorInterface/Pythia8Interface/test/ZJetsAnalyzer.cc
Expand Up @@ -14,7 +14,7 @@
#include "FWCore/ServiceRegistry/interface/Service.h"
#include "CommonTools/UtilAlgos/interface/TFileService.h"

#include "FWCore/Framework/interface/EDAnalyzer.h"
#include "FWCore/Framework/interface/one/EDAnalyzer.h"

#include "TH1.h"

Expand All @@ -31,21 +31,22 @@ struct ParticlePtGreater {
}
};

class ZJetsAnalyzer : public edm::EDAnalyzer {
class ZJetsAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
public:
//
explicit ZJetsAnalyzer(const edm::ParameterSet&);
virtual ~ZJetsAnalyzer(); // no need to delete ROOT stuff
// as it'll be deleted upon closing TFile
virtual ~ZJetsAnalyzer() = default; // no need to delete ROOT stuff
// as it'll be deleted upon closing TFile

virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
virtual void beginJob() override;
virtual void endRun(const edm::Run&, const edm::EventSetup&) override;
void analyze(const edm::Event&, const edm::EventSetup&) override;
void beginJob() override;
void beginRun(const edm::Run&, const edm::EventSetup&) override {}
void endRun(const edm::Run&, const edm::EventSetup&) override;

private:
edm::EDGetTokenT<GenEventInfoProduct> tokenGenEvent_;
edm::EDGetTokenT<edm::HepMCProduct> tokenHepMC_;
edm::EDGetTokenT<GenRunInfoProduct> tokenGenRun_;
const edm::EDGetTokenT<GenEventInfoProduct> tokenGenEvent_;
const edm::EDGetTokenT<edm::HepMCProduct> tokenHepMC_;
const edm::EDGetTokenT<GenRunInfoProduct> tokenGenRun_;

LeptonAnalyserHepMC LA;
JetInputHepMC JetInput;
Expand All @@ -66,11 +67,10 @@ ZJetsAnalyzer::ZJetsAnalyzer(const edm::ParameterSet& pset)
tokenGenRun_(consumes<GenRunInfoProduct, edm::InRun>(
edm::InputTag(pset.getUntrackedParameter("moduleLabel", std::string("generator")), ""))),
fHist2muMass(0) {
usesResource(TFileService::kSharedResource);
// actually, pset is NOT in use - we keep it here just for illustratory putposes
}

ZJetsAnalyzer::~ZJetsAnalyzer() { ; }

void ZJetsAnalyzer::beginJob() {
edm::Service<TFileService> fs;
fHist2muMass = fs->make<TH1D>("Hist2muMass", "2-mu inv. mass", 100, 60., 120.);
Expand All @@ -91,9 +91,7 @@ void ZJetsAnalyzer::endRun(const edm::Run& r, const edm::EventSetup&) {
std::ofstream testi("testi.dat");
double val, errval;

edm::Handle<GenRunInfoProduct> genRunInfoProduct;
//r.getByLabel("generator", genRunInfoProduct );
r.getByToken(tokenGenRun_, genRunInfoProduct);
const edm::Handle<GenRunInfoProduct>& genRunInfoProduct = r.getHandle(tokenGenRun_);

val = (double)genRunInfoProduct->crossSection();
std::cout << std::endl;
Expand Down Expand Up @@ -127,9 +125,7 @@ void ZJetsAnalyzer::analyze(const edm::Event& e, const edm::EventSetup&) {
icategories[0]++;

// here's an example of accessing GenEventInfoProduct
edm::Handle<GenEventInfoProduct> GenInfoHandle;
//e.getByLabel( "generator", GenInfoHandle );
e.getByToken(tokenGenEvent_, GenInfoHandle);
const edm::Handle<GenEventInfoProduct>& GenInfoHandle = e.getHandle(tokenGenEvent_);

double qScale = GenInfoHandle->qScale();
double pthat = (GenInfoHandle->hasBinningValues() ? (GenInfoHandle->binningValues())[0] : 0.0);
Expand All @@ -149,10 +145,7 @@ void ZJetsAnalyzer::analyze(const edm::Event& e, const edm::EventSetup&) {

// here's an example of accessing particles in the event record (HepMCProduct)
//
edm::Handle<edm::HepMCProduct> EvtHandle;
// find initial (unsmeared, unfiltered,...) HepMCProduct
//e.getByLabel("VtxSmeared", EvtHandle);
e.getByToken(tokenHepMC_, EvtHandle);
const edm::Handle<edm::HepMCProduct>& EvtHandle = e.getHandle(tokenHepMC_);

const HepMC::GenEvent* Evt = EvtHandle->GetEvent();

Expand Down