CMS 3D CMS Logo

Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes

LHEProducer Class Reference

Inheritance diagram for LHEProducer:
edm::EDFilter edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 LHEProducer (const edm::ParameterSet &params)
virtual ~LHEProducer ()

Protected Member Functions

virtual void beginJob ()
virtual bool beginRun (edm::Run &run, const edm::EventSetup &es)
virtual void endJob ()
virtual bool endRun (edm::Run &run, const edm::EventSetup &es)
virtual bool filter (edm::Event &event, const edm::EventSetup &es)

Private Member Functions

double matching (const HepMC::GenEvent *event, bool shower) const
void onBeforeHadronisation ()
void onInit ()
bool showeredEvent (const boost::shared_ptr< HepMC::GenEvent > &event)

Private Attributes

BranchingRatios branchingRatios
unsigned int eventsToPrint
double extCrossSect
double extFilterEff
std::auto_ptr< Hadronisationhadronisation
unsigned int index
std::auto_ptr< JetMatchingjetMatching
bool matchingDone
bool matchSummary
boost::shared_ptr< LHEEventpartonLevel
std::vector< int > removeResonances
boost::shared_ptr< LHERunInforunInfo
double weight

Detailed Description

Definition at line 33 of file LHEProducer.cc.


Constructor & Destructor Documentation

LHEProducer::LHEProducer ( const edm::ParameterSet params) [explicit]

Definition at line 69 of file LHEProducer.cc.

References branchingRatios, SurfaceDeformationFactory::create(), edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), hadronisation, i, jetMatching, matchSummary, onBeforeHadronisation(), onInit(), removeResonances, lhef::BranchingRatios::set(), and showeredEvent().

                                                      :
        eventsToPrint(params.getUntrackedParameter<unsigned int>("eventsToPrint", 0)),
        extCrossSect(params.getUntrackedParameter<double>("crossSection", -1.0)),
        extFilterEff(params.getUntrackedParameter<double>("filterEfficiency", -1.0)),
        matchSummary(false)
{
        hadronisation = Hadronisation::create(
                params.getParameter<edm::ParameterSet>("hadronisation"));

        if (params.exists("removeResonances"))
                removeResonances =
                        params.getParameter<std::vector<int> >(
                                                        "removeResonances");

        std::set<std::string> matchingCapabilities;
        if (params.exists("jetMatching")) {
                edm::ParameterSet jetParams =
                        params.getUntrackedParameter<edm::ParameterSet>(
                                                                "jetMatching");
                jetMatching = JetMatching::create(jetParams);

                matchingCapabilities = jetMatching->capabilities();
                hadronisation->matchingCapabilities(matchingCapabilities);
        }

        produces<edm::HepMCProduct>();
        produces<GenEventInfoProduct>();
        produces<GenRunInfoProduct, edm::InRun>();

        if (jetMatching.get()) {
                if (params.getUntrackedParameter<bool>(
                                        "preferShowerVetoCallback", true))
                        hadronisation->onShoweredEvent().connect(
                                sigc::mem_fun(*this,
                                              &LHEProducer::showeredEvent));
                hadronisation->onInit().connect(
                                sigc::mem_fun(*this, &LHEProducer::onInit));
                hadronisation->onBeforeHadronisation().connect(
                        sigc::mem_fun(*this,
                                      &LHEProducer::onBeforeHadronisation));

                matchSummary = matchingCapabilities.count("matchSummary");
                if (matchSummary) {
                        produces< std::vector<double> >("matchDeltaR");
                        produces< std::vector<double> >("matchDeltaPRel");
                }
        }

        // force total branching ratio for QCD/QED to 1
        for(int i = 0; i < 6; i++)
                branchingRatios.set(i);
        for(int i = 9; i < 23; i++)
                branchingRatios.set(i);
}
LHEProducer::~LHEProducer ( ) [virtual]

Definition at line 124 of file LHEProducer.cc.

{
}

Member Function Documentation

void LHEProducer::beginJob ( void  ) [protected, virtual]

Reimplemented from edm::EDFilter.

Definition at line 128 of file LHEProducer.cc.

References hadronisation.

{
        hadronisation->init();
}
bool LHEProducer::beginRun ( edm::Run run,
const edm::EventSetup es 
) [protected, virtual]

Reimplemented from edm::EDFilter.

Definition at line 139 of file LHEProducer.cc.

References edm::Run::getByLabel(), index, and runInfo.

{
        edm::Handle<LHERunInfoProduct> product;
        run.getByLabel("source", product);

        runInfo.reset(new LHERunInfo(*product));
        index = 0;

        return true;
}
void LHEProducer::endJob ( void  ) [protected, virtual]

Reimplemented from edm::EDFilter.

Definition at line 133 of file LHEProducer.cc.

References hadronisation, and jetMatching.

{
        hadronisation.reset();
        jetMatching.reset();
}
bool LHEProducer::endRun ( edm::Run run,
const edm::EventSetup es 
) [protected, virtual]

Reimplemented from edm::EDFilter.

Definition at line 150 of file LHEProducer.cc.

References PYTHIA6_MinBias_2360GeV_cff_py_GEN_FASTSIM::crossSection, lhef::LHERunInfo::XSec::error, extCrossSect, extFilterEff, hadronisation, edm::Run::put(), runInfo, and lhef::LHERunInfo::XSec::value.

{
        hadronisation->statistics();

        LHERunInfo::XSec crossSection;
        if (runInfo) {
                crossSection = runInfo->xsec();
                runInfo->statistics();
        }

        std::auto_ptr<GenRunInfoProduct> runInfo(new GenRunInfoProduct);

        runInfo->setInternalXSec(
                        GenRunInfoProduct::XSec(crossSection.value,
                                                crossSection.error));
        runInfo->setExternalXSecLO(extCrossSect);
        runInfo->setFilterEfficiency(extFilterEff);

        run.put(runInfo);

        runInfo.reset();

        return true;
}
bool LHEProducer::filter ( edm::Event event,
const edm::EventSetup es 
) [protected, virtual]

Implements edm::EDFilter.

Definition at line 175 of file LHEProducer.cc.

References beamvalidation::br, branchingRatios, eventsToPrint, lhef::BranchingRatios::getFactor(), hadronisation, index, info, jetMatching, matching(), matchingDone, matchSummary, partonLevel, removeResonances, query::result, runInfo, and weight.

{
        std::auto_ptr<edm::HepMCProduct> result(new edm::HepMCProduct);

        edm::Handle<LHEEventProduct> product;
        event.getByLabel("source", product);

        partonLevel.reset(new LHEEvent(runInfo, *product));
        if (!removeResonances.empty())
                partonLevel->removeResonances(removeResonances);

        if (product->pdf())
                partonLevel->setPDF(
                        std::auto_ptr<LHEEvent::PDF>(
                                new LHEEvent::PDF(*product->pdf())));

        hadronisation->setEvent(partonLevel);

        double br = branchingRatios.getFactor(hadronisation.get(),
                                              partonLevel);

        matchingDone = false;
        weight = 1.0;
        std::auto_ptr<HepMC::GenEvent> hadronLevel(hadronisation->hadronize());

        if (!hadronLevel.get()) {
                if (matchingDone) {
                        if (weight == 0.0)
                                partonLevel->count(LHERunInfo::kSelected, br);
                        else
                                partonLevel->count(LHERunInfo::kKilled,
                                                   br, weight);
                } else
                        partonLevel->count(LHERunInfo::kTried, br);
        }

        if (!matchingDone && jetMatching.get() && hadronLevel.get())
                weight = matching(hadronLevel.get(), false);

        if (weight == 0.0) {
                edm::LogInfo("Generator|LHEInterface")
                        << "Event got rejected by the "
                           "jet matching." << std::endl;

                if (hadronLevel.get()) {
                        partonLevel->count(LHERunInfo::kSelected);
                        hadronLevel.reset();
                }
        }

        if (!hadronLevel.get()) {
                event.put(result);
                std::auto_ptr<GenEventInfoProduct> info(
                                                new GenEventInfoProduct);
                event.put(info);
                return false;
        }

        partonLevel->count(LHERunInfo::kAccepted, br, weight);

        hadronLevel->set_event_number(++index);

        if (eventsToPrint) {
                eventsToPrint--;
                hadronLevel->print();
        }

        std::auto_ptr<GenEventInfoProduct> info(
                                new GenEventInfoProduct(hadronLevel.get()));
        result->addHepMCData(hadronLevel.release());
        event.put(result);
        event.put(info);

        if (jetMatching.get() && matchSummary) {
                std::auto_ptr< std::vector<double> > matchDeltaR(
                                                new std::vector<double>);
                std::auto_ptr< std::vector<double> > matchDeltaPRel(
                                                new std::vector<double>);

                typedef std::vector<JetMatching::JetPartonMatch> Matches;
                Matches matches = jetMatching->getMatchSummary();

                for(Matches::const_iterator iter = matches.begin();
                    iter != matches.end(); ++iter) {
                        if (!iter->isMatch())
                                continue;

                        matchDeltaR->push_back(iter->delta);
                        matchDeltaPRel->push_back(iter->jet.rho() /
                                                  iter->parton.rho() - 1.0);
                }

                event.put(matchDeltaR, "matchDeltaR");
                event.put(matchDeltaPRel, "matchDeltaPRel");
        }

        return true;
}
double LHEProducer::matching ( const HepMC::GenEvent *  event,
bool  shower 
) const [private]

Definition at line 274 of file LHEProducer.cc.

References event(), jetMatching, and partonLevel.

Referenced by filter(), and showeredEvent().

{
        if (!jetMatching.get())
                return 1.0;

        return jetMatching->match(partonLevel->asHepMCEvent().get(),
                                  event, shower);
}
void LHEProducer::onBeforeHadronisation ( ) [private]

Definition at line 295 of file LHEProducer.cc.

References jetMatching, and partonLevel.

Referenced by LHEProducer().

{
        jetMatching->beforeHadronisation(partonLevel);
}
void LHEProducer::onInit ( ) [private]

Definition at line 290 of file LHEProducer.cc.

References jetMatching, and runInfo.

Referenced by LHEProducer().

{
        jetMatching->init(runInfo);
}
bool LHEProducer::showeredEvent ( const boost::shared_ptr< HepMC::GenEvent > &  event) [private]

Definition at line 283 of file LHEProducer.cc.

References matching(), matchingDone, and weight.

Referenced by LHEProducer().

{
        weight = matching(event.get(), true);
        matchingDone = true;
        return weight == 0.0;
}

Member Data Documentation

Definition at line 66 of file LHEProducer.cc.

Referenced by filter(), and LHEProducer().

unsigned int LHEProducer::eventsToPrint [private]

Definition at line 52 of file LHEProducer.cc.

Referenced by filter().

double LHEProducer::extCrossSect [private]

Definition at line 57 of file LHEProducer.cc.

Referenced by endRun().

double LHEProducer::extFilterEff [private]

Definition at line 58 of file LHEProducer.cc.

Referenced by endRun().

std::auto_ptr<Hadronisation> LHEProducer::hadronisation [private]

Definition at line 54 of file LHEProducer.cc.

Referenced by beginJob(), endJob(), endRun(), filter(), and LHEProducer().

unsigned int LHEProducer::index [private]

Definition at line 63 of file LHEProducer.cc.

Referenced by beginRun(), and filter().

std::auto_ptr<JetMatching> LHEProducer::jetMatching [private]

Definition at line 55 of file LHEProducer.cc.

Referenced by endJob(), filter(), LHEProducer(), matching(), onBeforeHadronisation(), and onInit().

bool LHEProducer::matchingDone [private]

Definition at line 64 of file LHEProducer.cc.

Referenced by filter(), and showeredEvent().

bool LHEProducer::matchSummary [private]

Definition at line 59 of file LHEProducer.cc.

Referenced by filter(), and LHEProducer().

boost::shared_ptr<LHEEvent> LHEProducer::partonLevel [private]

Definition at line 61 of file LHEProducer.cc.

Referenced by filter(), matching(), and onBeforeHadronisation().

std::vector<int> LHEProducer::removeResonances [private]

Definition at line 53 of file LHEProducer.cc.

Referenced by filter(), and LHEProducer().

boost::shared_ptr<LHERunInfo> LHEProducer::runInfo [private]

Definition at line 62 of file LHEProducer.cc.

Referenced by beginRun(), endRun(), filter(), and onInit().

double LHEProducer::weight [private]

Definition at line 65 of file LHEProducer.cc.

Referenced by filter(), and showeredEvent().