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
#include "EventFilter/RPCRawToDigi/plugins/RPCDigiMerger.h"

#include <memory>

#include "DataFormats/Common/interface/Handle.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
#include "FWCore/Utilities/interface/CRC16.h"
#include "FWCore/Utilities/interface/InputTag.h"

#include "DataFormats/RPCDigi/interface/RPCDigiCollection.h"

using namespace edm;
using namespace std;

RPCDigiMerger::RPCDigiMerger(edm::ParameterSet const& config)
    : bx_minTwinMux_(config.getParameter<int>("bxMinTwinMux")),
      bx_maxTwinMux_(config.getParameter<int>("bxMaxTwinMux")),
      bx_minOMTF_(config.getParameter<int>("bxMinOMTF")),
      bx_maxOMTF_(config.getParameter<int>("bxMaxOMTF")),
      bx_minCPPF_(config.getParameter<int>("bxMinCPPF")),
      bx_maxCPPF_(config.getParameter<int>("bxMaxCPPF")) {
  produces<RPCDigiCollection>();
  simRPC_token_ = consumes<RPCDigiCollection>(config.getParameter<edm::InputTag>("inputTagSimRPCDigis"));
  // protection against empty InputTag to allow for Data/MC compatibility
  if (not config.getParameter<edm::InputTag>("inputTagTwinMuxDigis").label().empty()) {
    twinMux_token_ = consumes<RPCDigiCollection>(config.getParameter<edm::InputTag>("inputTagTwinMuxDigis"));
  }
  if (not config.getParameter<edm::InputTag>("inputTagOMTFDigis").label().empty()) {
    omtf_token_ = consumes<RPCDigiCollection>(config.getParameter<edm::InputTag>("inputTagOMTFDigis"));
  }
  if (not config.getParameter<edm::InputTag>("inputTagCPPFDigis").label().empty()) {
    cppf_token_ = consumes<RPCDigiCollection>(config.getParameter<edm::InputTag>("inputTagCPPFDigis"));
  }
}

RPCDigiMerger::~RPCDigiMerger() {}

void RPCDigiMerger::fillDescriptions(edm::ConfigurationDescriptions& descs) {
  edm::ParameterSetDescription desc;
  desc.add<edm::InputTag>("inputTagSimRPCDigis", edm::InputTag("simMuonRPCDigis", ""));
  desc.add<edm::InputTag>("inputTagTwinMuxDigis", edm::InputTag("", ""));
  desc.add<edm::InputTag>("inputTagOMTFDigis", edm::InputTag("", ""));
  desc.add<edm::InputTag>("inputTagCPPFDigis", edm::InputTag("", ""));
  desc.add<edm::InputTag>("InputLabel", edm::InputTag(" "));
  desc.add<int>("bxMinTwinMux", -2);
  desc.add<int>("bxMaxTwinMux", 2);
  desc.add<int>("bxMinOMTF", -3);
  desc.add<int>("bxMaxOMTF", 4);
  desc.add<int>("bxMinCPPF", -2);
  desc.add<int>("bxMaxCPPF", 2);

  descs.add("rpcDigiMerger", desc);
}

void RPCDigiMerger::produce(edm::Event& event, edm::EventSetup const& setup) {
  // Get the digis
  // new RPCDigiCollection
  std::unique_ptr<RPCDigiCollection> rpc_digi_collection(new RPCDigiCollection());

  //Check if its Data
  if (not(cppf_token_.isUninitialized() && omtf_token_.isUninitialized() && twinMux_token_.isUninitialized())) {
    // loop over TwinMux digis
    // protection against empty InputTag to allow for Data/MC compatibility
    if (not twinMux_token_.isUninitialized()) {
      Handle<RPCDigiCollection> TwinMux_digis;
      event.getByToken(twinMux_token_, TwinMux_digis);
      for (const auto&& rpcdgIt : (*TwinMux_digis)) {
        // The layerId
        const RPCDetId& rpcId = rpcdgIt.first;
        // Get the iterators over the digis associated with this LayerId
        const RPCDigiCollection::Range& range = rpcdgIt.second;
        rpc_digi_collection->put(range, rpcId);
      }
    }
    // loop over CPPF digis
    // protection against empty InputTag to allow for Data/MC compatibility
    if (not cppf_token_.isUninitialized()) {
      Handle<RPCDigiCollection> CPPF_digis;
      event.getByToken(cppf_token_, CPPF_digis);
      for (const auto&& rpcdgIt : (*CPPF_digis)) {
        // The layerId
        const RPCDetId& rpcId = rpcdgIt.first;
        // Get the iterators over the digis associated with this LayerId
        const RPCDigiCollection::Range& range = rpcdgIt.second;
        rpc_digi_collection->put(range, rpcId);
      }
    }
    // loop over OMTF digis
    // protection against empty InputTag to allow for Data/MC compatibility
    if (not omtf_token_.isUninitialized()) {
      Handle<RPCDigiCollection> OMTF_digis;
      event.getByToken(omtf_token_, OMTF_digis);
      for (const auto& rpcdgIt : (*OMTF_digis)) {
        // The layerId
        const RPCDetId& rpcId = rpcdgIt.first;
        // Get the iterators over the digis associated with this LayerId
        const RPCDigiCollection::Range& range = rpcdgIt.second;
        // accepts only rings: RE-2_R3 ; RE-1_R3 ; RE+1_R3 ; RE+2_R3 ;
        if (((rpcId.region() == -1 || rpcId.region() == 1) && (rpcId.ring() == 3) &&
             (rpcId.station() == 1 || rpcId.station() == 2))) {
          rpc_digi_collection->put(range, rpcId);
        }
      }
    }
  } else {  //its MC
    // SimRPCDigis collection
    Handle<RPCDigiCollection> SimRPC_digis;
    event.getByToken(simRPC_token_, SimRPC_digis);

    RPCDetId rpc_det_id;
    std::vector<RPCDigi> local_rpc_digis;

    // loop over SimRPC digis
    for (const auto& rpc_digi : (*SimRPC_digis)) {
      // The layerId
      const RPCDetId& rpcId = rpc_digi.first;
      // Get the iterators over the digis associated with this LayerId
      const RPCDigiCollection::Range& range = rpc_digi.second;

      if (rpcId != rpc_det_id) {
        if (!local_rpc_digis.empty()) {
          rpc_digi_collection->put(RPCDigiCollection::Range(local_rpc_digis.begin(), local_rpc_digis.end()),
                                   rpc_det_id);
          local_rpc_digis.clear();
        }
        rpc_det_id = rpcId;
      }
      for (std::vector<RPCDigi>::const_iterator id = range.first; id != range.second; id++) {
        const RPCDigi& dit = (*id);
        //Barrel
        if (rpcId.region() == 0) {
          //TwinMux
          if (dit.bx() >= bx_minTwinMux_ && dit.bx() <= bx_maxTwinMux_) {
            local_rpc_digis.push_back(dit);
          }
        }
        //EndCap
        if (rpcId.region() == -1 || rpcId.region() == 1) {
          //OMTF
          if (rpcId.ring() == 3 && (rpcId.station() == 1 || rpcId.station() == 2) && dit.bx() >= bx_minOMTF_ &&
              dit.bx() <= bx_maxOMTF_) {
            local_rpc_digis.push_back(dit);
          }
          //CPPF
          if (((rpcId.ring() == 2) || (rpcId.ring() == 3 && (rpcId.station() == 3 || rpcId.station() == 4))) &&
              (dit.bx() >= bx_minCPPF_ && dit.bx() <= bx_maxCPPF_)) {
            local_rpc_digis.push_back(dit);
          }
        }
      }
    }
    if (!local_rpc_digis.empty()) {
      rpc_digi_collection->put(RPCDigiCollection::Range(local_rpc_digis.begin(), local_rpc_digis.end()), rpc_det_id);
    }
  }
  // "put" into the event
  event.put(std::move(rpc_digi_collection));
}

#include "FWCore/Framework/interface/MakerMacros.h"
DEFINE_FWK_MODULE(RPCDigiMerger);