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

Initial framework setup #1

Merged
merged 4 commits into from
May 29, 2014
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.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
8 changes: 8 additions & 0 deletions MicroAODFormats/BuildFile.xml
Original file line number Diff line number Diff line change
@@ -0,0 +1,8 @@
<use name="DataFormats/Candidate"/>
<use name="DataFormats/PatCandidates"/>
<use name="CommonTools/CandUtils"/>
<use name="DataFormats/VertexReco"/>
<use name="rootrflx"/>
<export>
<lib name="1"/>
</export>
33 changes: 33 additions & 0 deletions MicroAODFormats/interface/DiPhotonCandidate.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,33 @@
#ifndef FLASHgg_DiPhotonCandidate_h
#define FLASHgg_DiPhotonCandidate_h

#include "DataFormats/Candidate/interface/CompositeCandidate.h"
#include "flashgg/MicroAODFormats/interface/Photon.h"
#include "DataFormats/VertexReco/interface/Vertex.h"


namespace flashgg {

class DiPhotonCandidate : public reco::CompositeCandidate {
public:
DiPhotonCandidate();
DiPhotonCandidate(edm::Ptr<flashgg::Photon>,edm::Ptr<flashgg::Photon>,edm::Ptr<reco::Vertex>);
~DiPhotonCandidate();

const reco::Vertex & getVertex() const {return *vertex_;};
const flashgg::Photon & getPhoton(size_type) const;

// These are an alternative to getPhoton above, but if we put them in we have to enforce ordering
// const flashggPhoton& leadingPhoton() const { return static_cast<const flashggPhoton&>(*daughter(0)); }
// const flashggPhoton& subLeadingPhoton() const { return static_cast<const flashggPhoton&>(*daughter(1)); }

private:
edm::Ptr<reco::Vertex> vertex_;

};


}


#endif
21 changes: 21 additions & 0 deletions MicroAODFormats/interface/Photon.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,21 @@
#ifndef FLASHgg_Photon_h
#define FLASHgg_Photon_h

#include "DataFormats/PatCandidates/interface/Photon.h"

namespace flashgg {

class Photon : public pat::Photon {

public:
Photon();
Photon(const pat::Photon &);
~Photon();
void setTestVariable(unsigned int i) { testVariable_ = i; }
unsigned int const getTestVariable() const {return testVariable_;};
private:
unsigned int testVariable_;
};
}

#endif
29 changes: 29 additions & 0 deletions MicroAODFormats/src/DiPhotonCandidate.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,29 @@
#include "flashgg/MicroAODFormats/interface/DiPhotonCandidate.h"
#include "flashgg/MicroAODFormats/interface/Photon.h"
#include "CommonTools/CandUtils/interface/AddFourMomenta.h"

//using namespace edm;
//using namespace std;
//using namespace reco;
using namespace flashgg;

DiPhotonCandidate::DiPhotonCandidate(){}

DiPhotonCandidate::~DiPhotonCandidate(){}

DiPhotonCandidate::DiPhotonCandidate(edm::Ptr<flashgg::Photon> photon1,edm::Ptr<flashgg::Photon> photon2, edm::Ptr<reco::Vertex> vertex) {
addDaughter(*photon1);
addDaughter(*photon2);
vertex_ = vertex;

// Adding momenta
// Needs its own object - but why?
// Copied from example
AddFourMomenta addP4;
addP4.set(*this);
}

const flashgg::Photon & DiPhotonCandidate::getPhoton(size_type i) const {
return static_cast<const flashgg::Photon&>(*daughter(i));
}

13 changes: 13 additions & 0 deletions MicroAODFormats/src/Photon.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,13 @@
#include "flashgg/MicroAODFormats/interface/Photon.h"

using namespace flashgg;

Photon::Photon() : pat::Photon() {
testVariable_ = 0;
}

Photon::Photon(const pat::Photon& aPhoton ) : pat::Photon(aPhoton) {
testVariable_ = 0;
}

Photon::~Photon() {}
12 changes: 12 additions & 0 deletions MicroAODFormats/src/classes.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,12 @@
#include "flashgg/MicroAODFormats/interface/Photon.h"
#include "flashgg/MicroAODFormats/interface/DiPhotonCandidate.h"
#include "DataFormats/Common/interface/Wrapper.h"
#include "DataFormats/VertexReco/interface/Vertex.h"

namespace { struct dictionary {
std::vector<flashgg::Photon> dummy0;
edm::Wrapper<std::vector<flashgg::Photon> > dummy1;
std::vector<flashgg::DiPhotonCandidate> dummy2;
edm::Wrapper<std::vector<flashgg::DiPhotonCandidate> > dummy3;
edm::Ptr<reco::Vertex> dummy4;
};}
11 changes: 11 additions & 0 deletions MicroAODFormats/src/classes_def.xml
Original file line number Diff line number Diff line change
@@ -0,0 +1,11 @@
<lcgdict>
<class name="edm::Ptr<reco::Vertex>"/>
<class name="flashgg::Photon"/>
<class name="std::vector<flashgg::Photon>"/>
<class name="edm::Wrapper<std::vector<flashgg::Photon> >"/>
<class name="flashgg::DiPhotonCandidate"/>
<class name="std::vector<flashgg::DiPhotonCandidate>"/>
<class name="edm::Wrapper<std::vector<flashgg::DiPhotonCandidate> >"/>
</lcgdict>


8 changes: 8 additions & 0 deletions MicroAODProducers/plugins/BuildFile.xml
Original file line number Diff line number Diff line change
@@ -0,0 +1,8 @@
<library name="FlashggProducerPlugins" file="*.cc">
<use name="DataFormats/Candidate"/>
<use name="DataFormats/PatCandidates"/>
<use name="CommonTools/CandUtils"/>
<use name="DataFormats/VertexReco"/>
<use name="flashgg/MicroAODFormats"/>
<flags EDM_PLUGIN="1"/>
</library>
63 changes: 63 additions & 0 deletions MicroAODProducers/plugins/DiPhotonProducer.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,63 @@
#include "FWCore/Framework/interface/EDProducer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "DataFormats/Common/interface/Handle.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Utilities/interface/EDMException.h"
#include "CommonTools/UtilAlgos/interface/TFileService.h"
#include "flashgg/MicroAODFormats/interface/DiPhotonCandidate.h"
#include "DataFormats/VertexReco/interface/Vertex.h"

using namespace edm;
using namespace std;

namespace flashgg {

class DiPhotonProducer : public EDProducer {

public:
DiPhotonProducer( const ParameterSet & );
private:
void produce( Event &, const EventSetup & ) override;
EDGetTokenT<View<reco::Vertex> > vertexToken_;
EDGetTokenT<View<flashgg::Photon> > photonToken_;
};

DiPhotonProducer::DiPhotonProducer(const ParameterSet & iConfig) :
vertexToken_(consumes<View<reco::Vertex> >(iConfig.getUntrackedParameter<InputTag> ("VertexTag", InputTag("offlineSlimmedPrimaryVertices")))),
photonToken_(consumes<View<flashgg::Photon> >(iConfig.getUntrackedParameter<InputTag> ("PhotonTag", InputTag("flashggPhotons"))))
{
produces<vector<flashgg::DiPhotonCandidate> >();
}

void DiPhotonProducer::produce( Event & evt, const EventSetup & ) {

Handle<View<reco::Vertex> > primaryVertices;
evt.getByToken(vertexToken_,primaryVertices);
const PtrVector<reco::Vertex>& pvPointers = primaryVertices->ptrVector();

Handle<View<flashgg::Photon> > photons;
evt.getByToken(photonToken_,photons);
const PtrVector<flashgg::Photon>& photonPointers = photons->ptrVector();

auto_ptr<vector<DiPhotonCandidate> > diPhotonColl(new vector<DiPhotonCandidate>);

for (unsigned int i = 0 ; i < photonPointers.size() ; i++) {
Ptr<flashgg::Photon> pp1 = photonPointers[i];
for (unsigned int j = i+1 ; j < photonPointers.size() ; j++) {
Ptr<flashgg::Photon> pp2 = photonPointers[j];
for (unsigned int k = 0 ; k < pvPointers.size() ; k++) {
Ptr<reco::Vertex> pvx = pvPointers[k];
diPhotonColl->push_back(DiPhotonCandidate(pp1,pp2,pvx));
}
}
}

evt.put(diPhotonColl);
}
}

typedef flashgg::DiPhotonProducer FlashggDiPhotonProducer;
DEFINE_FWK_MODULE(FlashggDiPhotonProducer);
55 changes: 55 additions & 0 deletions MicroAODProducers/plugins/PhotonProducer.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,55 @@
#include "FWCore/Framework/interface/EDProducer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "DataFormats/Common/interface/Handle.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "DataFormats/VertexReco/interface/Vertex.h"
#include "DataFormats/PatCandidates/interface/Photon.h"
#include "flashgg/MicroAODFormats/interface/Photon.h"



using namespace std;
using namespace edm;

namespace flashgg {

class PhotonProducer : public EDProducer {

public:
PhotonProducer( const ParameterSet & );
private:
void produce( Event &, const EventSetup & ) override;
EDGetTokenT<View<pat::Photon> > photonToken_;
};


PhotonProducer::PhotonProducer(const ParameterSet & iConfig) :
photonToken_(consumes<View<pat::Photon> >(iConfig.getUntrackedParameter<InputTag> ("PhotonTag", InputTag("slimmedPhotons"))))
{
produces<vector<flashgg::Photon> >();
}

void PhotonProducer::produce( Event & evt, const EventSetup & ) {

Handle<View<pat::Photon> > photons;
evt.getByToken(photonToken_,photons);
const PtrVector<pat::Photon>& photonPointers = photons->ptrVector();

auto_ptr<vector<flashgg::Photon> > photonColl(new vector<flashgg::Photon>);

for (unsigned int i = 0 ; i < photonPointers.size() ; i++) {
Ptr<pat::Photon> pp = photonPointers[i];
flashgg::Photon fg = flashgg::Photon(*pp);
fg.setTestVariable(i); // The index of the photon is as good an example of distinctive test data as any
photonColl->push_back(fg);
}

evt.put(photonColl);
}
}

typedef flashgg::PhotonProducer FlashggPhotonProducer;
DEFINE_FWK_MODULE(FlashggPhotonProducer);
47 changes: 47 additions & 0 deletions MicroAODProducers/test/scan_test_output.C
Original file line number Diff line number Diff line change
@@ -0,0 +1,47 @@
{
// Set up FW Lite for automatic loading of CMS libraries
// and data formats. As you may have other user-defined setup
// in your rootlogon.C, the CMS setup is executed only if the CMS
// environment is set up.
//
TString cmsswbase = getenv("CMSSW_BASE");
if (cmsswbase.Length() > 0) {
//
// The CMSSW environment is defined (this is true even for FW Lite)
// so set up the rest.
//
cout << "Loading FW Lite setup." << endl;
gSystem->Load("libFWCoreFWLite.so");
AutoLibraryLoader::enable();
gSystem->Load("libDataFormatsFWLite.so");
gSystem->Load("libDataFormatsPatCandidates.so");
}

TFile f("myOutputFile.root");
TTree *Events = f.Get("Events");
Events->Print();
Events->SetScanField(0);
Events->Scan("flashggDiPhotonCandidates_flashggDiPhotons__FLASHggTEST.obj.pt_:flashggPhotons_flashggPhotons__FLASHggTEST.obj.pt_:flashggPhotons_flashggPhotons__FLASHggTEST.obj.testVariable_");

#include "DataFormats/FWLite/interface/Handle.h"

fwlite::Event ev(&f);

int n = 0;

for( ev.toBegin(); ! ev.atEnd(); ++ev) {
fwlite::Handle<std::vector<flashgg::DiPhotonCandidate> > objs;
objs.getByLabel(ev,"flashggDiPhotons");
// now can access data
cout << "We have " << objs.ptr().size() << " diPhotons in event " << n++ << endl;

// Actually we can't access very much data
// FWLite doesn't really know about flashgg objects, so the code below doesn't really work.
// We'll leave figuring this out for another day

// for(std::vector<DiPhotonCandidate>::const_iterator dipho=objs->begin(); dipho != objs->end(); ++dipho) {
// cout << " DiPhoton with pt " << dipho->pt() << endl;
// }
}

}
20 changes: 20 additions & 0 deletions MicroAODProducers/test/simple_Producer_test.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,20 @@
import FWCore.ParameterSet.Config as cms

process = cms.Process("FLASHggTEST")

process.load("FWCore.MessageService.MessageLogger_cfi")

process.maxEvents = cms.untracked.PSet( input = cms.untracked.int32(10) )

process.source = cms.Source("PoolSource",fileNames=cms.untracked.vstring("file:/afs/cern.ch/work/g/gpetrucc/micro/70x/CMSSW_7_0_4/src/miniProd/WH_ZH_HToGG_M-125_13TeV_pythia6_PAT.root"))

process.flashggPhotons = cms.EDProducer('FlashggPhotonProducer',PhotonTag=cms.untracked.InputTag('slimmedPhotons'))
process.flashggDiPhotons = cms.EDProducer('FlashggDiPhotonProducer',PhotonTag=cms.untracked.InputTag('flashggPhotons'),VertexTag=cms.untracked.InputTag('offlineSlimmedPrimaryVertices'))

process.out = cms.OutputModule("PoolOutputModule", fileName = cms.untracked.string('myOutputFile.root'),
outputCommands = cms.untracked.vstring("drop *","keep *_flashgg*_*_*","keep *_offlineSlimmedPrimaryVertices_*_*")
)

process.p = cms.Path(process.flashggPhotons*process.flashggDiPhotons)

process.e = cms.EndPath(process.out)
20 changes: 20 additions & 0 deletions README.md
Original file line number Diff line number Diff line change
@@ -1,2 +1,22 @@
flashgg
=======
1. Set up MiniAOD with CMSSW_7_0_4 as described here: https://twiki.cern.ch/twiki/bin/viewauth/CMS/MiniAOD#Recipe (it takes a while)

2. Fork flashgg repository on the web here: https://github.com/cms-analysis/flashgg

3. In CMSSW_7_0_4/src, do commands something like:

```
git clone [email protected]:sethzenz/flashgg.git flashgg
git remote add upstream https://github.com/cms-analysis/flashgg
# see https://help.github.com/articles/fork-a-repo for more about this
```

Now build, a very basic workflow test, and an extremely primitive FWLite script:

```
scram b
cd flashgg/MicroAODProducers/test
cmsRun simple_Producer_test.py
root -b -q scan_test_output.C
```
1 change: 0 additions & 1 deletion analysis/README.md

This file was deleted.

1 change: 0 additions & 1 deletion microAOD/README.md

This file was deleted.