Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:06:49

0001 // -*- C++ -*-
0002 //
0003 // Package:   BeamSplash
0004 // Class:     BeamSPlash
0005 //
0006 //
0007 // Original Author:  Luca Malgeri
0008 
0009 #include <memory>
0010 #include <vector>
0011 #include <map>
0012 #include <set>
0013 
0014 // user include files
0015 
0016 #include "FWCore/Utilities/interface/InputTag.h"
0017 #include "FWCore/Framework/interface/Frameworkfwd.h"
0018 #include "FWCore/Framework/interface/Event.h"
0019 #include "FWCore/Framework/interface/MakerMacros.h"
0020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0022 
0023 #include "DPGAnalysis/Skims/interface/SelectHFMinBias.h"
0024 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
0025 #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
0026 
0027 using namespace edm;
0028 using namespace std;
0029 
0030 SelectHFMinBias::SelectHFMinBias(const edm::ParameterSet& iConfig) {}
0031 
0032 SelectHFMinBias::~SelectHFMinBias() {}
0033 
0034 bool SelectHFMinBias::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0035   edm::Handle<CaloTowerCollection> towers;
0036   iEvent.getByLabel("towerMaker", towers);
0037 
0038   int negTowers = 0;
0039   int posTowers = 0;
0040   for (CaloTowerCollection::const_iterator cal = towers->begin(); cal != towers->end(); ++cal) {
0041     for (unsigned int i = 0; i < cal->constituentsSize(); i++) {
0042       const DetId id = cal->constituent(i);
0043       if (id.det() == DetId::Hcal) {
0044         HcalSubdetector subdet = (HcalSubdetector(id.subdetId()));
0045         if (subdet == HcalForward) {
0046           if (cal->energy() > 3. && cal->eta() < -3.)
0047             negTowers++;
0048           if (cal->energy() > 3. && cal->eta() > 3.)
0049             posTowers++;
0050         }
0051       }
0052     }
0053   }
0054   if (negTowers > 0 && posTowers > 0)
0055     return true;
0056 
0057   return false;
0058 }
0059 
0060 //define this as a plug-in
0061 DEFINE_FWK_MODULE(SelectHFMinBias);