Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:30:52

0001 #ifndef SimPPS_RPIXDigiAnalyzer_h
0002 #define SimPPS_RPIXDigiAnalyzer_h
0003 
0004 #include "FWCore/Framework/interface/ConsumesCollector.h"
0005 #include "FWCore/ServiceRegistry/interface/Service.h"
0006 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0007 #include <FWCore/Framework/interface/one/EDAnalyzer.h>
0008 #include <DataFormats/CTPPSDetId/interface/CTPPSPixelDetId.h>
0009 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
0010 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
0011 #include "DataFormats/CTPPSDigi/interface/CTPPSPixelDigi.h"
0012 #include "DataFormats/CTPPSDigi/interface/CTPPSPixelDigiCollection.h"
0013 
0014 #include "SimDataFormats/CrossingFrame/interface/CrossingFrame.h"
0015 #include "SimDataFormats/CrossingFrame/interface/MixCollection.h"
0016 #include "FWCore/Framework/interface/EventSetup.h"
0017 #include "FWCore/Framework/interface/ESHandle.h"
0018 #include "DataFormats/Common/interface/Handle.h"
0019 #include "DataFormats/Common/interface/DetSetVector.h"
0020 #include <FWCore/Framework/interface/Event.h>
0021 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0022 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
0023 #include "CondFormats/PPSObjects/interface/PPSPixelTopology.h"
0024 #include "CondFormats/DataRecord/interface/PPSPixelTopologyRcd.h"
0025 
0026 #include <iostream>
0027 #include <string>
0028 
0029 #include "TH2D.h"
0030 
0031 #define SELECTED_PIXEL_ROW 89
0032 #define SELECTED_PIXEL_COLUMN 23
0033 #define SELECTED_UNITID 2014314496
0034 #define TG184 0.332655724
0035 
0036 #define USE_MIDDLE_OF_PIXEL_2
0037 #define CENTERX 1.05
0038 #define CENTERY -8.475
0039 
0040 using namespace edm;
0041 using namespace std;
0042 
0043 class PSimHit;
0044 
0045 namespace edm {
0046   class ParameterSet;
0047   class Event;
0048   class EventSetup;
0049 }  // namespace edm
0050 
0051 class PPSPixelDigiAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0052 public:
0053   explicit PPSPixelDigiAnalyzer(const edm::ParameterSet &pset);
0054   ~PPSPixelDigiAnalyzer() override;
0055   void endJob() override;
0056   void beginJob() override;
0057   void analyze(const edm::Event &event, const edm::EventSetup &eventSetup) override;
0058 
0059 private:
0060   TH2D *hAllHits;
0061   TH2D *hOneHitperEvent;
0062   TH2D *hOneHitperEvent2;
0063   TH2D *hOneHitperEventCenter;
0064   TH2D *hOneHitperEvent2Center;
0065   //TFile *file;
0066   std::string label_;
0067 
0068   int verbosity_;
0069   edm::EDGetTokenT<edm::PSimHitContainer> psim_token;
0070   edm::EDGetTokenT<edm::DetSetVector<CTPPSPixelDigi>> pixel_token;
0071   edm::ESGetToken<PPSPixelTopology, PPSPixelTopologyRcd> pixelTopologyToken_;
0072 
0073   unsigned int found_corresponding_digi_count_;
0074   unsigned int cumulative_cluster_size_[3];
0075 };
0076 
0077 PPSPixelDigiAnalyzer::PPSPixelDigiAnalyzer(const ParameterSet &pset)
0078     : hAllHits(nullptr),
0079       hOneHitperEvent(nullptr),
0080       hOneHitperEvent2(nullptr),
0081       hOneHitperEventCenter(nullptr),
0082       hOneHitperEvent2Center(nullptr) {
0083   label_ = pset.getUntrackedParameter<string>("label");
0084   verbosity_ = pset.getParameter<int>("Verbosity");
0085   edm::Service<TFileService> file;
0086 #ifdef USE_MIDDLE_OF_PIXEL
0087   hOneHitperEvent = file->make<TH2D>("OneHitperEvent", "One Hit per Event", 30, -8.511, -8.361, 20, 1, 1.1);
0088   hOneHitperEvent2 = file->make<TH2D>("OneHitperEvent2", "One Hit per Event 2", 30, -8.511, -8.361, 20, 1, 1.1);
0089 #else
0090   hOneHitperEvent = file->make<TH2D>("OneHitperEvent", "One Hit per Event", 30, -8.55, -8.4, 20, 1, 1.1);
0091   hOneHitperEvent2 = file->make<TH2D>("OneHitperEvent2", "One Hit per Event 2", 30, -8.55, -8.4, 20, 1, 1.1);
0092   hOneHitperEventCenter =
0093       file->make<TH2D>("OneHitperEventCenter", "One Hit per Event Center", 30, -0.075, 0.075, 20, -0.05, 0.05);
0094   hOneHitperEvent2Center =
0095       file->make<TH2D>("OneHitperEvent2Center", "Cluster Size 2", 30, -0.075, 0.075, 20, -0.05, 0.05);
0096 #endif
0097   file->cd();
0098   hAllHits = file->make<TH2D>("AllHits", "All Hits", 10, 1, 1.1, 10, -8.55, -8.4);
0099 
0100   psim_token = consumes<PSimHitContainer>(edm::InputTag("g4SimHits", "CTPPSPixelHits"));
0101   pixel_token = consumes<edm::DetSetVector<CTPPSPixelDigi>>(edm::InputTag(label_, ""));  //label=RPixDetDigitizer???
0102   pixelTopologyToken_ = esConsumes<PPSPixelTopology, PPSPixelTopologyRcd>();
0103 }
0104 
0105 PPSPixelDigiAnalyzer::~PPSPixelDigiAnalyzer() {}
0106 
0107 void PPSPixelDigiAnalyzer::beginJob() {
0108   found_corresponding_digi_count_ = 0;
0109   for (int a = 0; a < 3; a++)
0110     cumulative_cluster_size_[a] = 0;
0111 }
0112 void PPSPixelDigiAnalyzer::endJob() {
0113   edm::LogInfo("PPSPixelDigiAnalyzer") << "found_corresponding_digi_count_: " << found_corresponding_digi_count_;
0114   edm::LogInfo("PPSPixelDigiAnalyzer") << "Cumulative cluster size (1,2,>2) = " << cumulative_cluster_size_[0] << ", "
0115                                        << cumulative_cluster_size_[1] << ", " << cumulative_cluster_size_[2];
0116 }
0117 
0118 void PPSPixelDigiAnalyzer::analyze(const Event &event, const EventSetup &eventSetup) {
0119   if (verbosity_ > 0)
0120     edm::LogInfo("PPSPixelDigiAnalyzer") << "--- Run: " << event.id().run() << " Event: " << event.id().event();
0121 
0122   edm::LogInfo("PPSPixelDigiAnalyzer")
0123       << "                                                            I do love Pixels     ";
0124   Handle<PSimHitContainer> simHits;
0125   event.getByToken(psim_token, simHits);
0126 
0127   edm::Handle<edm::DetSetVector<CTPPSPixelDigi>> CTPPSPixelDigis;
0128   event.getByToken(pixel_token, CTPPSPixelDigis);
0129 
0130   edm::ESHandle<PPSPixelTopology> thePixelTopology = eventSetup.getHandle(pixelTopologyToken_);
0131 
0132   if (verbosity_ > 0)
0133     edm::LogInfo("PPSPixelDigiAnalyzer") << "\n=================== RPDA Starting SimHit access"
0134                                          << "  ===================";
0135 
0136   if (verbosity_ > 1)
0137     edm::LogInfo("PPSPixelDigiAnalyzer") << simHits->size();
0138 
0139   double selected_pixel_lower_x;
0140   double selected_pixel_lower_y;
0141   double selected_pixel_upper_x;
0142   double selected_pixel_upper_y;
0143   double myX = 0;
0144   double myY = 0;
0145 
0146   thePixelTopology->pixelRange(SELECTED_PIXEL_ROW,
0147                                SELECTED_PIXEL_COLUMN,
0148                                selected_pixel_lower_x,
0149                                selected_pixel_upper_x,
0150                                selected_pixel_lower_y,
0151                                selected_pixel_upper_y);
0152 
0153   double hit_inside_selected_pixel[2];
0154   bool found_hit_inside_selected_pixel = false;
0155 
0156   for (vector<PSimHit>::const_iterator hit = simHits->begin(); hit != simHits->end(); hit++) {
0157     LocalPoint entryP = hit->entryPoint();
0158     LocalPoint exitP = hit->exitPoint();
0159     LocalPoint midP((entryP.x() + exitP.x()) / 2., (entryP.y() + exitP.y()) / 2.);
0160 
0161 #ifdef USE_MIDDLE_OF_PIXEL
0162     if (entryP.x() > selected_pixel_lower_x && entryP.x() < selected_pixel_upper_x &&
0163         entryP.y() > (selected_pixel_lower_y + 0.115 * TG184) && entryP.y() < (selected_pixel_upper_y + 0.115 * TG184)
0164 #else
0165 #ifdef USE_MIDDLE_OF_PIXEL_2
0166     if (midP.x() > selected_pixel_lower_x && midP.x() < selected_pixel_upper_x && midP.y() > selected_pixel_lower_y &&
0167         midP.y() < selected_pixel_upper_y
0168 #else
0169     if (entryP.x() > selected_pixel_lower_x && entryP.x() < selected_pixel_upper_x &&
0170         entryP.y() > selected_pixel_lower_y && entryP.y() < selected_pixel_upper_y
0171 #endif
0172 #endif
0173         && hit->detUnitId() == SELECTED_UNITID) {
0174       hit_inside_selected_pixel[0] = entryP.x();
0175       hit_inside_selected_pixel[1] = entryP.y();
0176       found_hit_inside_selected_pixel = true;
0177 #ifdef USE_MIDDLE_OF_PIXEL_2
0178       hAllHits->Fill(midP.x(), midP.y());
0179       myX = midP.x();
0180       myY = midP.y();
0181 #else
0182       hAllHits->Fill(entryP.x(), entryP.y());
0183       myX = entryP.x();
0184       myY = entryP.y();
0185 #endif
0186       if (verbosity_ > 2)
0187         edm::LogInfo("PPSPixelDigiAnalyzer") << hit_inside_selected_pixel[0] << " " << hit_inside_selected_pixel[1];
0188     }
0189 
0190     //--------------
0191 
0192     if (verbosity_ > 1)
0193       if (hit->timeOfFlight() > 0) {
0194         edm::LogInfo("PPSPixelDigiAnalyzer")
0195             << "DetId: " << hit->detUnitId() << "PID: " << hit->particleType() << " TOF: " << hit->timeOfFlight()
0196             << " Proc Type: " << hit->processType() << " p: " << hit->pabs() << " x = " << entryP.x()
0197             << "   y = " << entryP.y() << "  z = " << entryP.z();
0198       }
0199   }
0200 
0201   if (verbosity_ > 0)
0202     edm::LogInfo("PPSPixelDigiAnalyzer") << "\n=================== RPDA Starting Digi access"
0203                                          << "  ===================";
0204   int numberOfDetUnits = 0;
0205 
0206   // Iterate on detector units
0207   edm::DetSetVector<CTPPSPixelDigi>::const_iterator DSViter = CTPPSPixelDigis->begin();
0208 
0209   for (; DSViter != CTPPSPixelDigis->end(); DSViter++) {
0210     ++numberOfDetUnits;
0211 
0212     DetId detIdObject(DSViter->detId());
0213     if (verbosity_ > 1)
0214       edm::LogInfo("PPSPixelDigiAnalyzer") << "DetId: " << DSViter->detId();
0215 
0216     bool found_corresponding_digi = false;
0217     unsigned int corresponding_digi_cluster_size = 0;
0218 
0219     // looping over digis in a unit id
0220     edm::DetSet<CTPPSPixelDigi>::const_iterator begin = (*DSViter).begin();
0221     edm::DetSet<CTPPSPixelDigi>::const_iterator end = (*DSViter).end();
0222 
0223     if (verbosity_ > 2) {
0224       edm::LogInfo("PPSPixelDigiAnalyzer") << "FF  " << DSViter->detId();
0225       for (edm::DetSet<CTPPSPixelDigi>::const_iterator di = begin; di != end; di++) {
0226         edm::LogInfo("PPSPixelDigiAnalyzer") << "           Digi row  " << di->row() << ", col " << di->column();
0227 
0228         // reconvert the digi to local coordinates
0229         double lx;
0230         double ly;
0231         double ux;
0232         double uy;
0233         unsigned int rr = di->row();
0234         unsigned int cc = di->column();
0235         thePixelTopology->pixelRange(rr, cc, lx, ux, ly, uy);
0236 
0237         edm::LogInfo("PPSPixelDigiAnalyzer")
0238             << " pixel boundaries x low up, y low up " << lx << " " << ux << " " << ly << " " << uy;
0239       }
0240     }
0241     if (DSViter->detId() == SELECTED_UNITID && found_hit_inside_selected_pixel) {
0242       for (edm::DetSet<CTPPSPixelDigi>::const_iterator di = begin; di != end; di++) {
0243         if (verbosity_ > 1)
0244           edm::LogInfo("PPSPixelDigiAnalyzer") << "           Digi row  " << di->row() << ", col " << di->column();
0245 
0246         if (di->row() == SELECTED_PIXEL_ROW && di->column() == SELECTED_PIXEL_COLUMN) {
0247           found_corresponding_digi_count_++;
0248           found_corresponding_digi = true;
0249           corresponding_digi_cluster_size = 1;
0250         }
0251       }
0252       //if coresponding digi found, re-loop to look for adjacent pixels
0253       if (found_corresponding_digi) {
0254         for (edm::DetSet<CTPPSPixelDigi>::const_iterator di = begin; di != end; di++) {
0255           if (verbosity_ > 1)
0256             edm::LogInfo("PPSPixelDigiAnalyzer") << "           Digi row  " << di->row() << ", col " << di->column();
0257 
0258           if ((di->row() == SELECTED_PIXEL_ROW + 1 && di->column() == SELECTED_PIXEL_COLUMN) ||
0259               (di->row() == SELECTED_PIXEL_ROW - 1 && di->column() == SELECTED_PIXEL_COLUMN) ||
0260               (di->row() == SELECTED_PIXEL_ROW && di->column() == SELECTED_PIXEL_COLUMN + 1) ||
0261               (di->row() == SELECTED_PIXEL_ROW && di->column() == SELECTED_PIXEL_COLUMN - 1)) {
0262             corresponding_digi_cluster_size++;
0263             edm::LogInfo("PPSPixelDigiAnalyzer") << "           Digi row  " << di->row() << ", col " << di->column();
0264           }
0265         }
0266       }
0267     }
0268     if (corresponding_digi_cluster_size > 0) {
0269       edm::LogInfo("PPSPixelDigiAnalyzer")
0270           << "corresponding_digi_cluster_size in the event: " << corresponding_digi_cluster_size;
0271       hOneHitperEvent->Fill(myY, myX);
0272       hOneHitperEventCenter->Fill(myY - CENTERY, myX - CENTERX);
0273       if (corresponding_digi_cluster_size < 3) {
0274         cumulative_cluster_size_[corresponding_digi_cluster_size - 1]++;
0275         if (corresponding_digi_cluster_size > 1) {
0276           hOneHitperEvent2->Fill(myY, myX);
0277           hOneHitperEvent2Center->Fill(myY - CENTERY, myX - CENTERX);
0278         }
0279       } else {
0280         cumulative_cluster_size_[2]++;
0281         hOneHitperEvent2->Fill(myY, myX);
0282         hOneHitperEvent2Center->Fill(myY - CENTERY, myX - CENTERX);
0283       }
0284     }
0285   }
0286 
0287   if (verbosity_ > 1)
0288     edm::LogInfo("PPSPixelDigiAnalyzer") << "numberOfDetUnits in the event: " << numberOfDetUnits;
0289 }
0290 
0291 #include "FWCore/PluginManager/interface/ModuleDef.h"
0292 #include "FWCore/Framework/interface/MakerMacros.h"
0293 
0294 DEFINE_FWK_MODULE(PPSPixelDigiAnalyzer);
0295 
0296 #endif