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

Run2-hcx139 Add a HLT filter for Laser events in HCAL #20049

Merged
merged 5 commits into from Aug 9, 2017
Merged
Show file tree
Hide file tree
Changes from 1 commit
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
53 changes: 53 additions & 0 deletions HLTrigger/special/interface/HLTHcalLaserMisfireFilter.h
@@ -0,0 +1,53 @@
#ifndef HLTriggerspecialHLTHcalLaserMisfireFilter_h
#define HLTriggerspecialHLTHcalLaserMisfireFilter_h
// -*- C++ -*-
//
// Package: HLTHcalLaserMisfireFilter
// Class: HLTHcalLaserMisfireFilter
//
/**\class HLTHcalLaserMisfireFilter HLTHcalLaserMisfireFilter.cc filter/HLTHcalCalibTypeFilter/src/HLTHcalCalibTypeFilter.cc

Description: Filter to select HCAL Laser fired events

Implementation:
<Notes on implementation>
*/

// include files
#include "HLTrigger/HLTcore/interface/HLTFilter.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"

#include "DataFormats/FEDRawData/interface/FEDRawData.h"
#include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
#include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
#include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"

namespace edm {
class ConfigurationDescriptions;
}

//
// class declaration
//

class HLTHcalLaserMisfireFilter : public HLTFilter {
public:
explicit HLTHcalLaserMisfireFilter(const edm::ParameterSet&);
~HLTHcalLaserMisfireFilter();
virtual bool hltFilter(edm::Event&, const edm::EventSetup&, trigger::TriggerFilterObjectWithRefs & filterproduct) const override;
static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);

private:

// ----------member data ---------------------------
edm::InputTag inputHBHE_, inputHF_;
edm::EDGetTokenT<HBHEDigiCollection> inputTokenHBHE_;
edm::EDGetTokenT<QIE10DigiCollection> inputTokenHF_;
double minFracDiffHBHELaser_, minFracHFLaser_;
int minADCHBHE_, minADCHF_;
bool testMode_;
};

#endif
136 changes: 136 additions & 0 deletions HLTrigger/special/src/HLTHcalLaserMisfireFilter.cc
@@ -0,0 +1,136 @@
// -*- C++ -*-
//
// Package: HLTHcalLaserMisfireFilter
// Class: HLTHcalLaserMisfireFilter
//
/**\class HLTHcalLaserMisfireFilter HLTHcalLaserMisfireFilter.cc

Description: Filter to select misfires of the HCAL laser.

Implementation:
Three HBHE RBX's ("bad") do not see the laser signal. Laser events are
selected by finding events with a large fraction of digis in HBHE above
a given threshold that do not have signals in the "bad" RBXes.

For HF events are selected if a large fraction of digis in HF are above
a given threshold
*/

// system include files
#include <memory>

// user include files
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"

#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/Common/interface/RefToBase.h"
#include "DataFormats/FEDRawData/interface/FEDNumbering.h"
#include "DataFormats/HcalDigi/interface/HcalCalibrationEventTypes.h"
#include "EventFilter/HcalRawToDigi/interface/HcalDCCHeader.h"
#include "HLTrigger/special/interface/HLTHcalLaserMisfireFilter.h"

//
// constructors and destructor
//
HLTHcalLaserMisfireFilter::HLTHcalLaserMisfireFilter(const edm::ParameterSet& config) : HLTFilter(config) {
inputHBHE_ = config.getParameter<edm::InputTag>("InputHBHE");
inputHF_ = config.getParameter<edm::InputTag>("InputHF");
minFracDiffHBHELaser_ = config.getParameter<double>("minFracDiffHBHELaser");
minFracHFLaser_ = config.getParameter<double>("minFracHFLaser");
minADCHBHE_ = config.getParameter<int>("minADCHBHE");
minADCHF_ = config.getParameter<int>("minADCHF"),
testMode_ = config.getUntrackedParameter<bool>("testMode",false);

inputTokenHBHE_ = consumes<HBHEDigiCollection>(inputHBHE_);
inputTokenHF_ = consumes<QIE10DigiCollection>(inputHF_);

}

HLTHcalLaserMisfireFilter::~HLTHcalLaserMisfireFilter() { }

void HLTHcalLaserMisfireFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<edm::InputTag>("InputHBHE",edm::InputTag("source"));
desc.add<edm::InputTag>("InputHF",edm::InputTag("source"));
desc.add<int>("minADCHBHE",10);
desc.add<int>("minADCHF",10);
desc.add<double>("minFracDiffHBHELaser",0.3);
desc.add<double>("minFracHFLaser",0.3);
desc.addUntracked<bool>("testMode",false);
descriptions.add("hltHcalLaserMisfireFilter",desc);
}

bool HLTHcalLaserMisfireFilter::hltFilter(edm::Event& iEvent, const edm::EventSetup&, trigger::TriggerFilterObjectWithRefs & filterproduct) const {

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

It dioes not make sense that this is an HLTFilter, rather than a normal EDFilter, as saveTags() does not make sense for DIGI collections.

// the filter object
if (saveTags()) {
filterproduct.addCollectionTag(inputHBHE_);
filterproduct.addCollectionTag(inputHF_);
}

edm::Handle<HBHEDigiCollection> hbhe_digi;
iEvent.getByToken(inputTokenHBHE_, hbhe_digi);

edm::Handle<QIE10DigiCollection> hf_digi;
iEvent.getByToken(inputTokenHF_, hf_digi);

// Count digis in good, bad RBXes. ('bad' RBXes see no laser signal)
double badrbxfracHBHE(0), goodrbxfracHBHE(0), goodrbxfracHF(0);
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

declare these where you first use them.

int NbadHBHE = 72*3; // 3 bad RBXes, 72 channels each
int NgoodHBHE= 2592*2-NbadHBHE; // remaining HBHE channels are 'good'
int NgoodHF = 864*4;

for (auto hbhe = hbhe_digi->begin(); hbhe != hbhe_digi->end(); ++ hbhe){
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

hi @bsunanda - please use a range based loop. Thanks

const HBHEDataFrame digi = (const HBHEDataFrame)(*hbhe);
HcalDetId myid=(HcalDetId)digi.id();
bool isbad(false); // assume channel is not bad

bool passCut(false);
int maxdigiHB(0);
for (int i=0; i<digi.size(); i++)
if(digi.sample(i).adc() > maxdigiHB) maxdigiHB = digi.sample(i).adc();
if (maxdigiHB > minADCHBHE_) passCut = true;

if ( myid.subdet()==HcalBarrel && myid.ieta()<0 ) {
if (myid.iphi()>=15 && myid.iphi()<=18) isbad=true;
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

please document these numbers (maybe they should be constants next to the "72" above? I presume if one changes the other has to change too?

else if (myid.iphi()>=27 && myid.iphi()<=34) isbad=true;
}

if (passCut) {
if (isbad) {
badrbxfracHBHE += 1.;
} else goodrbxfracHBHE += 1.;
}
}
goodrbxfracHBHE /= NgoodHBHE;
badrbxfracHBHE /= NbadHBHE;

for (auto hf = hf_digi->begin(); hf != hf_digi->end(); ++hf) {
const QIE10DataFrame digi = (const QIE10DataFrame)(*hf);
bool passCut(false);
int maxdigiHF(0);
for (int i=0; i<digi.samples(); i++)
if(digi[i].adc() > maxdigiHF) maxdigiHF = digi[i].adc();
if (maxdigiHF > minADCHF_) passCut = true;

if (passCut) {
goodrbxfracHF += 1.;
}
}
goodrbxfracHF /= NgoodHF;
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

for now, it seems this should be called NHF. as there is no concept of bad HF in the code?


if (testMode_)
edm::LogVerbatim("Report")
<< "******************************************************************\n"
<< "goodrbxfracHBHE: " << goodrbxfracHBHE << " badrbxfracHBHE: "
<< badrbxfracHBHE << " Size " << hbhe_digi->size() << "\n"
<< "goodrbxfracHF: " << goodrbxfracHF << " Size " << hf_digi->size()
<< "\n******************************************************************";

if (((goodrbxfracHBHE-badrbxfracHBHE) < minFracDiffHBHELaser_) ||
(goodrbxfracHF < minFracHFLaser_)) return false;

return true;
}
2 changes: 2 additions & 0 deletions HLTrigger/special/src/SealModule.cc
Expand Up @@ -33,6 +33,7 @@
#include "HLTrigger/special/interface/HLTPixelActivityHFSumEnergyFilter.h"
#include "HLTrigger/special/interface/HLTHFAsymmetryFilter.h"
#include "HLTrigger/special/interface/HLTTrackerHaloFilter.h"
#include "HLTrigger/special/interface/HLTHcalLaserMisfireFilter.h"

DEFINE_FWK_MODULE(HLTPixlMBFilt);
DEFINE_FWK_MODULE(HLTPixlMBForAlignmentFilter);
Expand All @@ -51,6 +52,7 @@ DEFINE_FWK_MODULE(HLTHcalNZSFilter);

DEFINE_FWK_MODULE(HLTCSCOverlapFilter);
DEFINE_FWK_MODULE(HLTCSCRing2or3Filter);
DEFINE_FWK_MODULE(HLTHcalLaserMisfireFilter);

typedef HLTCountNumberOfObject<SiStripRecHit2DCollection> HLTCountNumberOfSingleRecHit;
DEFINE_FWK_MODULE(HLTCountNumberOfSingleRecHit);
Expand Down