File indexing completed on 2024-04-06 12:25:51
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019 #include "FWCore/Framework/interface/Frameworkfwd.h"
0020 #include "FWCore/Framework/interface/stream/EDProducer.h"
0021 #include "FWCore/Utilities/interface/ESGetToken.h"
0022 #include "FWCore/Framework/interface/EventSetup.h"
0023 #include "FWCore/Framework/interface/Event.h"
0024 #include "FWCore/Framework/interface/MakerMacros.h"
0025
0026 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0027 #include "FWCore/Utilities/interface/StreamID.h"
0028 #include "FWCore/Utilities/interface/Exception.h"
0029
0030 #include "DataFormats/HcalRecHit/interface/HBHERecHitAuxSetter.h"
0031 #include "DataFormats/HcalRecHit/interface/CaloRecHitAuxSetter.h"
0032
0033 #include "RecoLocalCalo/HcalRecAlgos/interface/parsePlan1RechitCombiner.h"
0034
0035 #include "Geometry/Records/interface/HcalRecNumberingRecord.h"
0036 #include "Geometry/CaloTopology/interface/HcalTopology.h"
0037
0038
0039
0040
0041 class HBHEPlan1Combiner : public edm::stream::EDProducer<> {
0042 public:
0043 explicit HBHEPlan1Combiner(const edm::ParameterSet&);
0044 ~HBHEPlan1Combiner() override;
0045
0046 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0047
0048 private:
0049 void produce(edm::Event&, const edm::EventSetup&) override;
0050
0051
0052
0053
0054 edm::EDGetTokenT<HBHERecHitCollection> tok_rechits_;
0055 bool ignorePlan1Topology_;
0056 bool usePlan1Mode_;
0057
0058
0059 std::unique_ptr<AbsPlan1RechitCombiner> combiner_;
0060
0061
0062 edm::ESGetToken<HcalTopology, HcalRecNumberingRecord> htopoToken_;
0063 };
0064
0065
0066
0067
0068 HBHEPlan1Combiner::HBHEPlan1Combiner(const edm::ParameterSet& conf)
0069 : ignorePlan1Topology_(conf.getParameter<bool>("ignorePlan1Topology")),
0070 usePlan1Mode_(conf.getParameter<bool>("usePlan1Mode")),
0071 combiner_(parsePlan1RechitCombiner(conf.getParameter<edm::ParameterSet>("algorithm"))) {
0072
0073 if (!combiner_.get())
0074 throw cms::Exception("HBHEPlan1BadConfig") << "Invalid Plan1RechitCombiner algorithm configuration" << std::endl;
0075
0076
0077 tok_rechits_ = consumes<HBHERecHitCollection>(conf.getParameter<edm::InputTag>("hbheInput"));
0078 produces<HBHERecHitCollection>();
0079
0080
0081 htopoToken_ = esConsumes<HcalTopology, HcalRecNumberingRecord>();
0082 }
0083
0084 HBHEPlan1Combiner::~HBHEPlan1Combiner() {
0085
0086
0087 }
0088
0089
0090
0091
0092
0093
0094 void HBHEPlan1Combiner::produce(edm::Event& iEvent, const edm::EventSetup& eventSetup) {
0095 using namespace edm;
0096
0097
0098 const HcalTopology& htopo = eventSetup.getData(htopoToken_);
0099 combiner_->setTopo(&htopo);
0100
0101
0102 const bool plan1Mode = ignorePlan1Topology_ ? usePlan1Mode_ : htopo.getMergePositionFlag();
0103
0104
0105 Handle<HBHERecHitCollection> inputRechits;
0106 iEvent.getByToken(tok_rechits_, inputRechits);
0107
0108
0109 std::unique_ptr<HBHERecHitCollection> outputRechits = std::make_unique<HBHERecHitCollection>();
0110 outputRechits->reserve(inputRechits->size());
0111
0112
0113
0114 combiner_->clear();
0115 for (typename HBHERecHitCollection::const_iterator it = inputRechits->begin(); it != inputRechits->end(); ++it) {
0116
0117 if (plan1Mode && CaloRecHitAuxSetter::getBit(it->auxPhase1(), HBHERecHitAuxSetter::OFF_TDC_TIME))
0118 combiner_->add(*it);
0119 else
0120 outputRechits->push_back(*it);
0121 }
0122
0123
0124 combiner_->combine(&*outputRechits);
0125
0126
0127 iEvent.put(std::move(outputRechits));
0128 }
0129
0130 #define add_param_set(name) \
0131 edm::ParameterSetDescription name; \
0132 name.setAllowAnything(); \
0133 desc.add<edm::ParameterSetDescription>(#name, name)
0134
0135
0136 void HBHEPlan1Combiner::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0137 edm::ParameterSetDescription desc;
0138
0139 desc.add<edm::InputTag>("hbheInput");
0140 desc.add<bool>("ignorePlan1Topology");
0141 desc.add<bool>("usePlan1Mode");
0142
0143 add_param_set(algorithm);
0144
0145 descriptions.addDefault(desc);
0146 }
0147
0148
0149 DEFINE_FWK_MODULE(HBHEPlan1Combiner);