File indexing completed on 2023-03-17 11:18:38
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021 #include <vector>
0022 #include <iostream>
0023 #include <TMath.h>
0024 #include <TRandom3.h>
0025
0026
0027 #include "FWCore/Framework/interface/Frameworkfwd.h"
0028 #include "FWCore/Framework/interface/global/EDProducer.h"
0029 #include "FWCore/Framework/interface/Event.h"
0030 #include "FWCore/Framework/interface/MakerMacros.h"
0031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0032 #include "DataFormats/Math/interface/Point3D.h"
0033
0034
0035 #include "DataFormats/CastorReco/interface/CastorCell.h"
0036 #include "DataFormats/HcalRecHit/interface/CastorRecHit.h"
0037 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
0038
0039
0040
0041
0042
0043 class CastorCellProducer : public edm::global::EDProducer<> {
0044 public:
0045 explicit CastorCellProducer(const edm::ParameterSet&);
0046
0047 private:
0048 void beginJob() override;
0049 void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0050 void endJob() override;
0051
0052
0053 typedef math::XYZPointD Point;
0054 typedef ROOT::Math::RhoZPhiPoint CellPoint;
0055 typedef std::vector<reco::CastorCell> CastorCellCollection;
0056 edm::EDGetTokenT<CastorRecHitCollection> tok_input_;
0057 };
0058
0059
0060
0061
0062 namespace {
0063 constexpr double MYR2D = 180 / M_PI;
0064 }
0065
0066
0067
0068
0069
0070
0071
0072
0073 CastorCellProducer::CastorCellProducer(const edm::ParameterSet& iConfig) {
0074 tok_input_ = consumes<CastorRecHitCollection>(
0075 edm::InputTag(iConfig.getUntrackedParameter<std::string>("inputprocess", "castorreco")));
0076
0077 produces<CastorCellCollection>();
0078
0079 }
0080
0081
0082
0083
0084
0085
0086
0087 void CastorCellProducer::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0088 using namespace edm;
0089 using namespace reco;
0090 using namespace TMath;
0091
0092
0093
0094 edm::Handle<CastorRecHitCollection> InputRecHits;
0095 iEvent.getByToken(tok_input_, InputRecHits);
0096
0097 auto OutputCells = std::make_unique<CastorCellCollection>();
0098
0099
0100
0101 LogDebug("CastorCellProducer") << "1. entering CastorCellProducer ";
0102
0103 for (size_t i = 0; i < InputRecHits->size(); ++i) {
0104 const CastorRecHit& rh = (*InputRecHits)[i];
0105 int sector = rh.id().sector();
0106 int module = rh.id().module();
0107 double energy = rh.energy();
0108 int zside = rh.id().zside();
0109
0110
0111 double zCell = 0.;
0112 double phiCell;
0113 double rhoCell;
0114
0115
0116 if (module < 3) {
0117
0118 if (module == 1)
0119 zCell = 14415;
0120 if (module == 2)
0121 zCell = 14464;
0122 } else {
0123
0124 zCell = 14534 + (module - 3) * 92;
0125 }
0126
0127
0128 double castorphi[16];
0129 for (int j = 0; j < 16; j++) {
0130 castorphi[j] = -2.94524 + j * 0.3927;
0131 }
0132 if (sector > 8) {
0133 phiCell = castorphi[sector - 9];
0134 } else {
0135 phiCell = castorphi[sector + 7];
0136 }
0137
0138
0139 if (zside <= 0)
0140 zCell = -1 * zCell;
0141
0142
0143 rhoCell = 88.5;
0144
0145
0146 CellPoint tempcellposition(rhoCell, zCell, phiCell);
0147 Point cellposition(tempcellposition);
0148
0149 LogDebug("CastorCellProducer") << "cell number: " << i + 1 << std::endl
0150 << "rho: " << cellposition.rho() << " phi: " << cellposition.phi() * MYR2D
0151 << " eta: " << cellposition.eta() << std::endl
0152 << "x: " << cellposition.x() << " y: " << cellposition.y()
0153 << " z: " << cellposition.z();
0154
0155 if (energy > 0.) {
0156 CastorCell newCell(energy, cellposition);
0157 OutputCells->push_back(newCell);
0158 }
0159
0160 }
0161
0162 LogDebug("CastorCellProducer") << "total number of cells in the event: " << InputRecHits->size();
0163
0164 iEvent.put(std::move(OutputCells));
0165 }
0166
0167
0168 void CastorCellProducer::beginJob() { LogDebug("CastorCellProducer") << "Starting CastorCellProducer"; }
0169
0170
0171 void CastorCellProducer::endJob() { LogDebug("CastorCellProducer") << "Ending CastorCellProducer"; }
0172
0173
0174 DEFINE_FWK_MODULE(CastorCellProducer);