Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:10:59

0001 #ifndef EventFilter_SiPixelRawToDigi_interface_PixelUnpackingRegions_h
0002 #define EventFilter_SiPixelRawToDigi_interface_PixelUnpackingRegions_h
0003 
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "FWCore/Framework/interface/ESWatcher.h"
0007 #include "FWCore/Utilities/interface/ESGetToken.h"
0008 #include "DataFormats/Math/interface/Point3D.h"
0009 #include "DataFormats/Math/interface/Vector3D.h"
0010 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0011 #include "DataFormats/Candidate/interface/LeafCandidate.h"
0012 
0013 #include "CondFormats/DataRecord/interface/SiPixelFedCablingMapRcd.h"
0014 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingMap.h"
0015 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingTree.h"
0016 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0017 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0018 #include "FWCore/Framework/interface/ConsumesCollector.h"
0019 
0020 #include <cmath>
0021 #include <vector>
0022 #include <set>
0023 
0024 /** \class PixelUnpackingRegions
0025  *
0026  * Input: One or several collections of Candidate-based seeds with their objects
0027  *        defining the directions of unpacking regions; separate deltaPhi and maxZ 
0028  *        tolerances could be given to each input collection.
0029  * Output: FED ids and module detIds that need to be unpacked
0030  *
0031  */
0032 class PixelUnpackingRegions {
0033 public:
0034   /// container to define regions for objects of interest in each event by:
0035   ///   object direction
0036   ///   dphi max distance from region direction to center of a pixel module
0037   ///   maxZ max projected z of a pixel module (when projecting along region direction onto beamline)
0038   struct Region {
0039     Region(const math::XYZVector& dir, float dphi = 0.5f, float maxz = 24.f) : v(dir), dPhi(dphi), maxZ(maxz) {
0040       cosphi = v.x() / v.rho();
0041       sinphi = v.y() / v.rho();
0042       atantheta = v.z() / v.rho();
0043     }
0044     math::XYZVector v;
0045     float dPhi, maxZ;
0046     float cosphi, sinphi, atantheta;
0047   };
0048 
0049   PixelUnpackingRegions(const edm::ParameterSet&, edm::ConsumesCollector&& iC);
0050 
0051   ~PixelUnpackingRegions() {}
0052 
0053   /// has to be run during each event
0054   void run(const edm::Event& e, const edm::EventSetup& es);
0055 
0056   /// check whether a FED has to be unpacked
0057   bool mayUnpackFED(unsigned int fed_n) const;
0058 
0059   /// check whether a module has to be unpacked
0060   bool mayUnpackModule(unsigned int id) const;
0061 
0062   /// full set of module ids to unpack
0063   const std::set<unsigned int>* modulesToUnpack() const { return &modules_; }
0064 
0065   /// various informational accessors:
0066   unsigned int nFEDs() const { return feds_.size(); }
0067   unsigned int nModules() const { return modules_.size(); }
0068   unsigned int nBarrelModules() const;
0069   unsigned int nForwardModules() const;
0070   unsigned int nRegions() const { return nreg_; }
0071 
0072   struct Module {
0073     float phi;
0074     float x, y, z;
0075     unsigned int id;
0076     unsigned int fed;
0077 
0078     Module() {}
0079     Module(float ph) : phi(ph), x(0.f), y(0.f), z(0.f), id(0), fed(0) {}
0080 
0081     bool operator<(const Module& m) const {
0082       if (phi < m.phi)
0083         return true;
0084       if (phi == m.phi && id < m.id)
0085         return true;
0086       return false;
0087     }
0088   };
0089 
0090 private:
0091   // input parameters
0092   std::vector<edm::InputTag> inputs_;
0093   std::vector<double> dPhi_;
0094   std::vector<double> maxZ_;
0095   edm::InputTag beamSpotTag_;
0096 
0097   edm::EDGetTokenT<reco::BeamSpot> tBeamSpot;
0098   std::vector<edm::EDGetTokenT<reco::CandidateView>> tCandidateView;
0099 
0100   std::set<unsigned int> feds_;
0101   std::set<unsigned int> modules_;
0102   unsigned int nreg_;
0103 
0104   /// run by the run method: (re)initialize the cabling data when it's necessary
0105   void initialize(const edm::EventSetup& es);
0106 
0107   // add a new direction of a region of interest
0108   void addRegion(Region& r);
0109 
0110   // gather info into feds_ and modules_ from a range of a Module vector
0111   void gatherFromRange(Region& r, std::vector<Module>::const_iterator, std::vector<Module>::const_iterator);
0112 
0113   // addRegion for a local (BPIX or +-FPIX) container
0114   void addRegionLocal(Region& r, std::vector<Module>& container, const Module& lo, const Module& hi);
0115 
0116   // local containers of barrel and endcaps Modules sorted by phi
0117   std::vector<Module> phiBPIX_;
0118   std::vector<Module> phiFPIXp_;
0119   std::vector<Module> phiFPIXm_;
0120 
0121   std::unique_ptr<SiPixelFedCablingTree> cabling_;
0122   math::XYZPoint beamSpot_;
0123 
0124   edm::ESWatcher<SiPixelFedCablingMapRcd> watcherSiPixelFedCablingMap_;
0125   edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> trackerGeomToken_;
0126   edm::ESGetToken<SiPixelFedCablingMap, SiPixelFedCablingMapRcd> cablingMapToken_;
0127 };
0128 
0129 #endif  // EventFilter_SiPixelRawToDigi_interface_PixelUnpackingRegions_h