Fork me on GitHub

Changeset 3241a0e in git for readers/DelphesCMSFWLite.cpp


Ignore:
Timestamp:
May 21, 2015, 5:43:25 PM (9 years ago)
Author:
Pavel Demin <pavel.demin@…>
Branches:
ImprovedOutputFile, Timing, dual_readout, llp, master
Children:
e8070b6
Parents:
95aa610
Message:

add LHE weights DelphesCMSFWLite

File:
1 edited

Legend:

Unmodified
Added
Removed
  • readers/DelphesCMSFWLite.cpp

    r95aa610 r3241a0e  
    5656#include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
    5757#include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
     58#include "SimDataFormats/GeneratorProducts/interface/LHEEventProduct.h"
     59#include "SimDataFormats/GeneratorProducts/interface/WeightsInfo.h"
    5860
    5961using namespace std;
     
    6264
    6365void ConvertInput(fwlite::Event &event, Long64_t eventCounter,
    64   ExRootTreeBranch *branchEvent, DelphesFactory *factory,
    65   TObjArray *allParticleOutputArray, TObjArray *stableParticleOutputArray,
    66   TObjArray *partonOutputArray)
     66  ExRootTreeBranch *branchEvent, ExRootTreeBranch *branchRwgt,
     67  DelphesFactory *factory, TObjArray *allParticleOutputArray,
     68  TObjArray *stableParticleOutputArray, TObjArray *partonOutputArray)
    6769{
    6870  fwlite::Handle< GenEventInfoProduct > handleGenEventInfo;
    6971
     72  fwlite::Handle< LHEEventProduct > handleLHEEvent;
     73
    7074  fwlite::Handle< vector< reco::GenParticle > > handleParticle;
    7175  vector< reco::GenParticle >::const_iterator itParticle;
     
    7579
    7680  handleGenEventInfo.getByLabel(event, "generator");
     81  handleLHEEvent.getByLabel(event, "source");
    7782  handleParticle.getByLabel(event, "genParticles");
    7883
    7984  HepMCEvent *element;
     85  Weight *weight;
    8086  Candidate *candidate;
    8187  TDatabasePDG *pdg;
     
    8692  Double_t px, py, pz, e, mass;
    8793  Double_t x, y, z;
     94
     95  const vector< gen::WeightsInfo > &vectorWeightsInfo = handleLHEEvent->weights();
     96  vector< gen::WeightsInfo >::const_iterator itWeightsInfo;
    8897
    8998  element = static_cast<HepMCEvent *>(branchEvent->NewEntry());
     
    108117  element->ReadTime = 0.0;
    109118  element->ProcTime = 0.0;
     119
     120  for(itWeightsInfo = vectorWeightsInfo.begin(); itWeightsInfo != vectorWeightsInfo.end(); ++itWeightsInfo)
     121  {
     122    weight = static_cast<Weight *>(branchRwgt->NewEntry());
     123    weight->Weight = itWeightsInfo->wgt;
     124  }
    110125
    111126  pdg = TDatabasePDG::Instance();
     
    186201  TStopwatch eventStopWatch;
    187202  ExRootTreeWriter *treeWriter = 0;
    188   ExRootTreeBranch *branchEvent = 0;
     203  ExRootTreeBranch *branchEvent = 0, *branchRwgt = 0;
    189204  ExRootConfReader *confReader = 0;
    190205  Delphes *modularDelphes = 0;
     
    226241
    227242    branchEvent = treeWriter->NewBranch("Event", HepMCEvent::Class());
     243    branchRwgt = treeWriter->NewBranch("Rwgt", Weight::Class());
    228244
    229245    confReader = new ExRootConfReader;
     
    268284      for(event.toBegin(); !event.atEnd() && !interrupted; ++event)
    269285      {
    270         ConvertInput(event, eventCounter, branchEvent, factory,
     286        ConvertInput(event, eventCounter, branchEvent, branchRwgt, factory,
    271287          allParticleOutputArray, stableParticleOutputArray, partonOutputArray);
    272288        modularDelphes->ProcessTask();
Note: See TracChangeset for help on using the changeset viewer.