HcalLaserHFFilter2012

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167
// -*- C++ -*-
//
// Package:    HcalLaserHFFilter2012
// Class:      HcalLaserHFFilter2012
//
/**\class HcalLaserHFFilter2012 HcalLaserHFFilter2012.cc UserCode/HcalLaserHFFilter2012/src/HcalLaserHFFilter2012.cc

 Description: [one line class summary]

 Implementation:
     [Notes on implementation]
*/
//
// Original Author:
//         Created:  Fri Oct 19 13:15:44 EDT 2012
//
//

// system include files
#include <memory>

// user include files
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/one/EDFilter.h"

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

#include "FWCore/ParameterSet/interface/ParameterSet.h"

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

#include <cmath>
#include <iostream>
#include <sstream>
#include <fstream>

#include "FWCore/Framework/interface/Run.h"
#include "FWCore/Framework/interface/LuminosityBlock.h"

#include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"

#include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
#include "CalibFormats/HcalObjects/interface/HcalCoderDb.h"

//
// class declaration
//

class HcalLaserHFFilter2012 : public edm::one::EDFilter<> {
public:
  explicit HcalLaserHFFilter2012(const edm::ParameterSet&);
  ~HcalLaserHFFilter2012() override;

  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);

private:
  bool filter(edm::Event&, const edm::EventSetup&) override;
  void endJob() override;

  // ----------member data ---------------------------
  bool verbose_;  // if set to true, then the run:LS:event for any event failing the cut will be printed out
  std::string
      prefix_;  // prefix will be printed before any event if verbose mode is true, in order to make searching for events easier
  int minCalibChannelsHFLaser_;  // set minimum number of HF Calib events that causes an event to be considered a bad (i.e., HF laser) event
  edm::InputTag digiLabel_;
  edm::EDGetTokenT<HcalCalibDigiCollection> tok_calib_;

  bool WriteBadToFile_;
  bool forceFilterTrue_;
  std::ofstream outfile_;
};

//
// constants, enums and typedefs
//

//
// static data member definitions
//

//
// constructors and destructor
//
HcalLaserHFFilter2012::HcalLaserHFFilter2012(const edm::ParameterSet& ps) {
  //now do what ever initialization is needed
  verbose_ = ps.getUntrackedParameter<bool>("verbose", false);
  prefix_ = ps.getUntrackedParameter<std::string>("prefix", "");
  minCalibChannelsHFLaser_ = ps.getUntrackedParameter<int>("minCalibChannelsHFLaser", 10);
  edm::InputTag digi_default("hcalDigis");
  digiLabel_ = ps.getUntrackedParameter<edm::InputTag>("digiLabel", digi_default);
  tok_calib_ = consumes<HcalCalibDigiCollection>(digiLabel_);
  WriteBadToFile_ = ps.getUntrackedParameter<bool>("WriteBadToFile", false);
  if (WriteBadToFile_)
    outfile_.open("badHcalLaserList_hffilter.txt");
  forceFilterTrue_ = ps.getUntrackedParameter<bool>("forceFilterTrue", false);

}  // HcalLaserHFFilter2012::HcalLaserHFFilter2012  constructor

HcalLaserHFFilter2012::~HcalLaserHFFilter2012() {
  // do anything here that needs to be done at destruction time
  // (e.g. close files, deallocate resources etc.)
}

//
// member functions
//

// ------------ method called on each new Event  ------------
bool HcalLaserHFFilter2012::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
  // Step 1:: try to get calib digi collection.
  // Return true if collection not found?  Or false?  What should default behavior be?
  edm::Handle<HcalCalibDigiCollection> calib_digi;
  if (!(iEvent.getByToken(tok_calib_, calib_digi))) {
    edm::LogWarning("HcalLaserHFFilter2012") << digiLabel_ << " calib_digi not available";
    return true;
  }

  if (!(calib_digi.isValid())) {
    edm::LogWarning("HcalLaserHFFilter2012") << digiLabel_ << " calib_digi is not valid";
    return true;
  }

  // Step 2:  Count HF digi calib channels
  int ncalibHF = 0;  // this will track number of HF digi channels

  for (HcalCalibDigiCollection::const_iterator Calibiter = calib_digi->begin(); Calibiter != calib_digi->end();
       ++Calibiter) {
    const HcalCalibDataFrame digi = (const HcalCalibDataFrame)(*Calibiter);
    if (digi.zsMarkAndPass())
      continue;  // skip digis labeled as "mark and pass" in NZS events
    HcalCalibDetId myid = (HcalCalibDetId)digi.id();
    if (myid.hcalSubdet() != HcalForward)
      continue;
    ++ncalibHF;
    if (ncalibHF >= minCalibChannelsHFLaser_) {
      if (verbose_)
        std::cout << prefix_ << iEvent.id().run() << ":" << iEvent.luminosityBlock() << ":" << iEvent.id().event()
                  << std::endl;
      if (WriteBadToFile_)
        outfile_ << iEvent.id().run() << ":" << iEvent.luminosityBlock() << ":" << iEvent.id().event() << std::endl;
      if (forceFilterTrue_)
        return true;  // if special input boolean set, always return true, regardless of filter decision
      else
        return false;
    }
  }

  return true;
}  // HcalLaserHFFilter2012::filter

// ------------ method called once each job just after ending the event loop  ------------
void HcalLaserHFFilter2012::endJob() {
  if (WriteBadToFile_)
    outfile_.close();
}

// ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
void HcalLaserHFFilter2012::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
  //The following says we do not know what parameters are allowed so do no validation
  // Please change this to state exactly what you do use, even if it is no parameters
  edm::ParameterSetDescription desc;
  desc.setUnknown();
  descriptions.addDefault(desc);
}
//define this as a plug-in
DEFINE_FWK_MODULE(HcalLaserHFFilter2012);