Skip to content

Commit

Permalink
Merge pull request #1672 from Martin-Grunewald/HLTPrescaler
Browse files Browse the repository at this point in the history
Make HLTPrescaler a stream filter with a global cache accounting
  • Loading branch information
nclopezo committed Dec 5, 2013
2 parents 5e5e80d + 93294ea commit 48bd9b5
Show file tree
Hide file tree
Showing 3 changed files with 37 additions and 11 deletions.
22 changes: 18 additions & 4 deletions HLTrigger/HLTcore/interface/HLTPrescaler.h
Expand Up @@ -12,26 +12,39 @@
* \author Philipp Schieferdecker
*/

#include <atomic>
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EDFilter.h"
#include "FWCore/Framework/interface/stream/EDFilter.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/PrescaleService/interface/PrescaleService.h"
#include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
namespace edm {
class ConfigurationDescriptions;
}

namespace trigger {
struct Efficiency {
Efficiency(): eventCount_(0),acceptCount_(0) { }
mutable std::atomic<unsigned int> eventCount_;
mutable std::atomic<unsigned int> acceptCount_;
};
}

class HLTPrescaler : public edm::EDFilter
class HLTPrescaler : public edm::stream::EDFilter<edm::GlobalCache<trigger::Efficiency> >
{
public:
//
// construction/destruction
//
explicit HLTPrescaler(edm::ParameterSet const& iConfig);
explicit HLTPrescaler(edm::ParameterSet const& iConfig, const trigger::Efficiency* efficiency);
virtual ~HLTPrescaler();

static std::unique_ptr<trigger::Efficiency> initializeGlobalCache(edm::ParameterSet const&) {
return std::unique_ptr<trigger::Efficiency>(new trigger::Efficiency());
}



//
// member functions
Expand All @@ -40,7 +53,8 @@ class HLTPrescaler : public edm::EDFilter
virtual void beginLuminosityBlock(edm::LuminosityBlock const&lb,
edm::EventSetup const& iSetup) override;
virtual bool filter(edm::Event& iEvent,edm::EventSetup const& iSetup) override;
virtual void endJob() override;
virtual void endStream() override;
static void globalEndJob(const trigger::Efficiency* efficiency);


private:
Expand Down
24 changes: 19 additions & 5 deletions HLTrigger/HLTcore/plugins/HLTPrescaler.cc
Expand Up @@ -29,7 +29,7 @@ const unsigned int HLTPrescaler::prescaleSeed_ = 65537;
///////////////////////////////////////////////////////////////////////////////

//_____________________________________________________________________________
HLTPrescaler::HLTPrescaler(edm::ParameterSet const& iConfig) :
HLTPrescaler::HLTPrescaler(edm::ParameterSet const& iConfig, const trigger::Efficiency* efficiency) :
prescaleSet_(0)
, prescaleFactor_(1)
, eventCount_(0)
Expand All @@ -53,7 +53,6 @@ HLTPrescaler::~HLTPrescaler()

}


////////////////////////////////////////////////////////////////////////////////
// implementation of member functions
////////////////////////////////////////////////////////////////////////////////
Expand Down Expand Up @@ -127,12 +126,27 @@ bool HLTPrescaler::filter(edm::Event& iEvent, const edm::EventSetup&)


//_____________________________________________________________________________
void HLTPrescaler::endJob()
void HLTPrescaler::endStream()
{
//since these are std::atomic, it is safe to increment them
// even if multiple endStreams are being called.
globalCache()->eventCount_ += eventCount_;
globalCache()->acceptCount_ += acceptCount_;
return;
}

//_____________________________________________________________________________
void HLTPrescaler::globalEndJob(const trigger::Efficiency* efficiency)
{
unsigned int accept(efficiency->acceptCount_);
unsigned int event (efficiency->eventCount_);
edm::LogInfo("PrescaleSummary")
<< acceptCount_<< "/" <<eventCount_
<< accept << "/" << event
<< " ("
<< 100.*acceptCount_/static_cast<double>(std::max(1u,eventCount_))
<< 100.*accept/static_cast<double>(std::max(1u,event))
<< "% of events accepted).";
return;
}

#include "FWCore/Framework/interface/MakerMacros.h"
DEFINE_FWK_MODULE(HLTPrescaler);
2 changes: 0 additions & 2 deletions HLTrigger/HLTcore/plugins/SealModule.cc
Expand Up @@ -3,7 +3,6 @@

#include "FWCore/Framework/interface/MakerMacros.h"

#include "HLTrigger/HLTcore/interface/HLTPrescaler.h"
#include "HLTrigger/HLTcore/interface/HLTEventAnalyzerAOD.h"
#include "HLTrigger/HLTcore/interface/HLTEventAnalyzerRAW.h"
#include "HLTrigger/HLTcore/interface/TriggerSummaryAnalyzerAOD.h"
Expand All @@ -12,7 +11,6 @@
#include "HLTrigger/HLTcore/interface/TriggerSummaryProducerRAW.h"
#include "HLTrigger/HLTcore/interface/HLTPrescaleRecorder.h"

DEFINE_FWK_MODULE(HLTPrescaler);
DEFINE_FWK_MODULE(HLTEventAnalyzerAOD);
DEFINE_FWK_MODULE(HLTEventAnalyzerRAW);
DEFINE_FWK_MODULE(TriggerSummaryAnalyzerAOD);
Expand Down

0 comments on commit 48bd9b5

Please sign in to comment.