Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-02-14 14:24:18

0001 // -*- C++ -*-
0002 //
0003 // Package:    Castor
0004 // Class:      CastorCellProducer
0005 //
0006 
0007 /**\class CastorCellProducer CastorCellProducer.cc RecoLocalCalo/Castor/src/CastorCellProducer.cc
0008 
0009  Description: CastorCell Reconstruction Producer. Produce CastorCells from CastorRecHits.
0010  Implementation:
0011 */
0012 
0013 //
0014 // Original Author:  Hans Van Haevermaet, Benoit Roland
0015 //         Created:  Wed Jul  9 14:00:40 CEST 2008
0016 //
0017 //
0018 
0019 // system include
0020 #include <memory>
0021 #include <vector>
0022 #include <iostream>
0023 #include <TMath.h>
0024 #include <TRandom3.h>
0025 
0026 // user include
0027 #include "FWCore/Framework/interface/Frameworkfwd.h"
0028 #include "FWCore/Framework/interface/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 // Castor Object include
0035 #include "DataFormats/CastorReco/interface/CastorCell.h"
0036 #include "DataFormats/HcalRecHit/interface/CastorRecHit.h"
0037 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
0038 
0039 //
0040 // class declaration
0041 //
0042 
0043 class CastorCellProducer : public edm::EDProducer {
0044 public:
0045   explicit CastorCellProducer(const edm::ParameterSet&);
0046   ~CastorCellProducer() override;
0047 
0048 private:
0049   void beginJob() override;
0050   void produce(edm::Event&, const edm::EventSetup&) override;
0051   void endJob() override;
0052 
0053   // member data
0054   typedef math::XYZPointD Point;
0055   typedef ROOT::Math::RhoZPhiPoint CellPoint;
0056   typedef std::vector<reco::CastorCell> CastorCellCollection;
0057   edm::EDGetTokenT<CastorRecHitCollection> tok_input_;
0058 };
0059 
0060 //
0061 // constants, enums and typedefs
0062 //
0063 
0064 const double MYR2D = 180 / M_PI;
0065 
0066 //
0067 // static data member definitions
0068 //
0069 
0070 //
0071 // constructor and destructor
0072 //
0073 
0074 CastorCellProducer::CastorCellProducer(const edm::ParameterSet& iConfig) {
0075   tok_input_ = consumes<CastorRecHitCollection>(
0076       edm::InputTag(iConfig.getUntrackedParameter<std::string>("inputprocess", "castorreco")));
0077   // register your products
0078   produces<CastorCellCollection>();
0079   // now do what ever other initialization is needed
0080 }
0081 
0082 CastorCellProducer::~CastorCellProducer() {
0083   // do anything here that needs to be done at desctruction time
0084   // (e.g. close files, deallocate resources etc.)
0085 }
0086 
0087 //
0088 // member functions
0089 //
0090 
0091 // ------------ method called to produce the data  ------------
0092 
0093 void CastorCellProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0094   using namespace edm;
0095   using namespace reco;
0096   using namespace TMath;
0097 
0098   // Produce CastorCells from CastorRecHits
0099 
0100   edm::Handle<CastorRecHitCollection> InputRecHits;
0101   iEvent.getByToken(tok_input_, InputRecHits);
0102 
0103   auto OutputCells = std::make_unique<CastorCellCollection>();
0104 
0105   // looping over all CastorRecHits
0106 
0107   LogDebug("CastorCellProducer") << "1. entering CastorCellProducer ";
0108 
0109   for (size_t i = 0; i < InputRecHits->size(); ++i) {
0110     const CastorRecHit& rh = (*InputRecHits)[i];
0111     int sector = rh.id().sector();
0112     int module = rh.id().module();
0113     double energy = rh.energy();
0114     int zside = rh.id().zside();
0115 
0116     // define CastorCell properties
0117     double zCell = 0.;
0118     double phiCell;
0119     double rhoCell;
0120 
0121     // set z position of the cell
0122     if (module < 3) {
0123       // starting in EM section
0124       if (module == 1)
0125         zCell = 14415;
0126       if (module == 2)
0127         zCell = 14464;
0128     } else {
0129       // starting in HAD section
0130       zCell = 14534 + (module - 3) * 92;
0131     }
0132 
0133     // set phi position of the cell
0134     double castorphi[16];
0135     for (int j = 0; j < 16; j++) {
0136       castorphi[j] = -2.94524 + j * 0.3927;
0137     }
0138     if (sector > 8) {
0139       phiCell = castorphi[sector - 9];
0140     } else {
0141       phiCell = castorphi[sector + 7];
0142     }
0143 
0144     // add condition to select in eta sides
0145     if (zside <= 0)
0146       zCell = -1 * zCell;
0147 
0148     // set rho position of the cell (inner radius 3.7cm, outer radius 14cm)
0149     rhoCell = 88.5;
0150 
0151     // store cell position
0152     CellPoint tempcellposition(rhoCell, zCell, phiCell);
0153     Point cellposition(tempcellposition);
0154 
0155     LogDebug("CastorCellProducer") << "cell number: " << i + 1 << std::endl
0156                                    << "rho: " << cellposition.rho() << " phi: " << cellposition.phi() * MYR2D
0157                                    << " eta: " << cellposition.eta() << std::endl
0158                                    << "x: " << cellposition.x() << " y: " << cellposition.y()
0159                                    << " z: " << cellposition.z();
0160 
0161     if (energy > 0.) {
0162       CastorCell newCell(energy, cellposition);
0163       OutputCells->push_back(newCell);
0164     }
0165 
0166   }  // end loop over CastorRecHits
0167 
0168   LogDebug("CastorCellProducer") << "total number of cells in the event: " << InputRecHits->size();
0169 
0170   iEvent.put(std::move(OutputCells));
0171 }
0172 
0173 // ------------ method called once each job just before starting event loop  ------------
0174 void CastorCellProducer::beginJob() { LogDebug("CastorCellProducer") << "Starting CastorCellProducer"; }
0175 
0176 // ------------ method called once each job just after ending the event loop  ------------
0177 void CastorCellProducer::endJob() { LogDebug("CastorCellProducer") << "Ending CastorCellProducer"; }
0178 
0179 //define this as a plug-in
0180 DEFINE_FWK_MODULE(CastorCellProducer);