Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-08-16 22:25:32

0001 // -*- C++ -*-
0002 //
0003 // Package:    CalibTracker/SiPixelLorentzAnglePCLWorker
0004 // Class:      SiPixelLorentzAnglePCLWorker
0005 //
0006 /**\class SiPixelLorentzAnglePCLWorker SiPixelLorentzAnglePCLWorker.cc CalibTracker/SiPixelLorentzAnglePCLWorker/src/SiPixelLorentzAnglePCLWorker.cc
0007  Description: generates the intermediate ALCAPROMPT dataset for the measurement of the SiPixel Lorentz Angle in the Prompt Calibration Loop
0008  Implementation:
0009      Books and fills 2D histograms of the drift vs depth in bins of pixel module rings to be fed into the SiPixelLorentzAnglePCLHarvester
0010 */
0011 //
0012 // Original Author:  mmusich
0013 //         Created:  Sat, 29 May 2021 14:46:19 GMT
0014 //
0015 //
0016 
0017 // system includes
0018 #include <string>
0019 #include <fmt/printf.h>
0020 
0021 // user include files
0022 #include "CalibTracker/Records/interface/SiPixelTemplateDBObjectESProducerRcd.h"
0023 #include "CalibTracker/SiPixelLorentzAngle/interface/SiPixelLorentzAngleCalibrationStruct.h"
0024 #include "CondFormats/SiPixelObjects/interface/SiPixelTemplateDBObject.h"
0025 #include "CondFormats/SiPixelTransient/interface/SiPixelTemplate.h"
0026 #include "CondFormats/SiPixelTransient/interface/SiPixelTemplateDefs.h"
0027 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0028 #include "DQMServices/Core/interface/DQMStore.h"
0029 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
0030 #include "DataFormats/GeometryVector/interface/LocalVector.h"
0031 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
0032 #include "DataFormats/TrackReco/interface/TrackExtra.h"
0033 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0034 #include "DataFormats/TrackerCommon/interface/PixelBarrelName.h"
0035 #include "DataFormats/TrackerCommon/interface/PixelEndcapName.h"
0036 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0037 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
0038 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
0039 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
0040 #include "FWCore/Framework/interface/Event.h"
0041 #include "FWCore/Framework/interface/ESWatcher.h"
0042 #include "FWCore/Framework/interface/Frameworkfwd.h"
0043 #include "FWCore/Framework/interface/MakerMacros.h"
0044 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0045 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
0046 #include "Geometry/CommonTopologies/interface/StripTopology.h"
0047 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0048 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0049 #include "Geometry/TrackerGeometryBuilder/interface/PixelTopologyMap.h"
0050 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0051 #include "Geometry/CommonTopologies/interface/SimplePixelTopology.h"
0052 #include "RecoTracker/TransientTrackingRecHit/interface/TkTransientTrackingRecHitBuilder.h"
0053 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
0054 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
0055 #include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
0056 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0057 
0058 // ROOT includes
0059 #include <TTree.h>
0060 #include <TFile.h>
0061 #include <fstream>
0062 
0063 //
0064 // class declaration
0065 //
0066 
0067 static const int maxpix = 1000;
0068 struct Pixinfo {
0069   int npix;
0070   float row[maxpix];
0071   float col[maxpix];
0072   float adc[maxpix];
0073   float x[maxpix];
0074   float y[maxpix];
0075 };
0076 
0077 struct Hit {
0078   float x;
0079   float y;
0080   double alpha;
0081   double beta;
0082   double gamma;
0083 };
0084 struct Clust {
0085   float x;
0086   float y;
0087   float charge;
0088   int size_x;
0089   int size_y;
0090   int maxPixelCol;
0091   int maxPixelRow;
0092   int minPixelCol;
0093   int minPixelRow;
0094 };
0095 struct Rechit {
0096   float x;
0097   float y;
0098 };
0099 
0100 class SiPixelLorentzAnglePCLWorker : public DQMEDAnalyzer {
0101 public:
0102   explicit SiPixelLorentzAnglePCLWorker(const edm::ParameterSet&);
0103   ~SiPixelLorentzAnglePCLWorker() override = default;
0104 
0105   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0106 
0107 private:
0108   void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0109 
0110   void analyze(edm::Event const&, edm::EventSetup const&) override;
0111 
0112   void dqmBeginRun(edm::Run const&, edm::EventSetup const&) override;
0113 
0114   void dqmEndRun(edm::Run const&, edm::EventSetup const&);
0115 
0116   const Pixinfo fillPix(const SiPixelCluster& LocPix, const PixelTopology* topol) const;
0117   const std::pair<LocalPoint, LocalPoint> surface_deformation(const PixelTopology* topol,
0118                                                               TrajectoryStateOnSurface& tsos,
0119                                                               const SiPixelRecHit* recHitPix) const;
0120   // ------------ member data ------------
0121   SiPixelLorentzAngleCalibrationHistograms iHists;
0122 
0123   // template stuff
0124   edm::ESWatcher<SiPixelTemplateDBObjectESProducerRcd> watchSiPixelTemplateRcd_;
0125   const SiPixelTemplateDBObject* templateDBobject_;
0126   std::vector<SiPixelTemplateStore> thePixelTemp_;
0127 
0128   std::string folder_;
0129   bool notInPCL_;
0130   std::string filename_;
0131   std::vector<std::string> newmodulelist_;
0132 
0133   // tree branches barrel
0134   int run_;
0135   long int event_;
0136   int lumiblock_;
0137   int bx_;
0138   int orbit_;
0139   int module_;
0140   int ladder_;
0141   int layer_;
0142   int isflipped_;
0143   float pt_;
0144   float eta_;
0145   float phi_;
0146   double chi2_;
0147   double ndof_;
0148   Pixinfo pixinfo_;
0149   Hit simhit_, trackhit_;
0150   Clust clust_;
0151   Rechit rechit_;
0152   Rechit rechitCorr_;
0153   float trackhitCorrX_;
0154   float trackhitCorrY_;
0155   float qScale_;
0156   float rQmQt_;
0157 
0158   // tree branches forward
0159   int sideF_;
0160   int diskF_;
0161   int bladeF_;
0162   int panelF_;
0163   int moduleF_;
0164   Pixinfo pixinfoF_;
0165   Hit simhitF_, trackhitF_;
0166   Clust clustF_;
0167   Rechit rechitF_;
0168   Rechit rechitCorrF_;
0169   float trackhitCorrXF_;
0170   float trackhitCorrYF_;
0171   float qScaleF_;
0172   float rQmQtF_;
0173 
0174   // parameters from config file
0175   double ptmin_;
0176   double normChi2Max_;
0177   std::vector<int> clustSizeYMin_;
0178   int clustSizeXMax_;
0179   double residualMax_;
0180   double clustChargeMaxPerLength_;
0181   int hist_depth_;
0182   int hist_drift_;
0183 
0184   std::unique_ptr<TFile> hFile_;
0185   std::unique_ptr<TTree> SiPixelLorentzAngleTreeBarrel_;
0186   std::unique_ptr<TTree> SiPixelLorentzAngleTreeForward_;
0187 
0188   // es consumes
0189   edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomEsToken_;
0190   edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> topoEsToken_;
0191   edm::ESGetToken<SiPixelTemplateDBObject, SiPixelTemplateDBObjectESProducerRcd> siPixelTemplateEsToken_;
0192   edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> topoPerEventEsToken_;
0193   edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomPerEventEsToken_;
0194 
0195   // event consumes
0196   edm::EDGetTokenT<TrajTrackAssociationCollection> t_trajTrack;
0197 };
0198 
0199 //
0200 // constructors and destructor
0201 //
0202 SiPixelLorentzAnglePCLWorker::SiPixelLorentzAnglePCLWorker(const edm::ParameterSet& iConfig)
0203     : folder_(iConfig.getParameter<std::string>("folder")),
0204       notInPCL_(iConfig.getParameter<bool>("notInPCL")),
0205       filename_(iConfig.getParameter<std::string>("fileName")),
0206       newmodulelist_(iConfig.getParameter<std::vector<std::string>>("newmodulelist")),
0207       ptmin_(iConfig.getParameter<double>("ptMin")),
0208       normChi2Max_(iConfig.getParameter<double>("normChi2Max")),
0209       clustSizeYMin_(iConfig.getParameter<std::vector<int>>("clustSizeYMin")),
0210       clustSizeXMax_(iConfig.getParameter<int>("clustSizeXMax")),
0211       residualMax_(iConfig.getParameter<double>("residualMax")),
0212       clustChargeMaxPerLength_(iConfig.getParameter<double>("clustChargeMaxPerLength")),
0213       hist_depth_(iConfig.getParameter<int>("binsDepth")),
0214       hist_drift_(iConfig.getParameter<int>("binsDrift")),
0215       geomEsToken_(esConsumes<edm::Transition::BeginRun>()),
0216       topoEsToken_(esConsumes<edm::Transition::BeginRun>()),
0217       siPixelTemplateEsToken_(esConsumes<edm::Transition::BeginRun>()),
0218       topoPerEventEsToken_(esConsumes()),
0219       geomPerEventEsToken_(esConsumes()) {
0220   t_trajTrack = consumes<TrajTrackAssociationCollection>(iConfig.getParameter<edm::InputTag>("src"));
0221 
0222   // now do what ever initialization is needed
0223   int bufsize = 64000;
0224 
0225   //    create tree structure
0226   //    Barrel pixel
0227   if (notInPCL_) {
0228     hFile_ = std::make_unique<TFile>(filename_.c_str(), "RECREATE");
0229     SiPixelLorentzAngleTreeBarrel_ =
0230         std::make_unique<TTree>("SiPixelLorentzAngleTreeBarrel_", "SiPixel LorentzAngle tree barrel", bufsize);
0231     SiPixelLorentzAngleTreeBarrel_->Branch("run", &run_, "run/I", bufsize);
0232     SiPixelLorentzAngleTreeBarrel_->Branch("event", &event_, "event/l", bufsize);
0233     SiPixelLorentzAngleTreeBarrel_->Branch("lumiblock", &lumiblock_, "lumiblock/I", bufsize);
0234     SiPixelLorentzAngleTreeBarrel_->Branch("bx", &bx_, "bx/I", bufsize);
0235     SiPixelLorentzAngleTreeBarrel_->Branch("orbit", &orbit_, "orbit/I", bufsize);
0236     SiPixelLorentzAngleTreeBarrel_->Branch("module", &module_, "module/I", bufsize);
0237     SiPixelLorentzAngleTreeBarrel_->Branch("ladder", &ladder_, "ladder/I", bufsize);
0238     SiPixelLorentzAngleTreeBarrel_->Branch("layer", &layer_, "layer/I", bufsize);
0239     SiPixelLorentzAngleTreeBarrel_->Branch("isflipped", &isflipped_, "isflipped/I", bufsize);
0240     SiPixelLorentzAngleTreeBarrel_->Branch("pt", &pt_, "pt/F", bufsize);
0241     SiPixelLorentzAngleTreeBarrel_->Branch("eta", &eta_, "eta/F", bufsize);
0242     SiPixelLorentzAngleTreeBarrel_->Branch("phi", &phi_, "phi/F", bufsize);
0243     SiPixelLorentzAngleTreeBarrel_->Branch("chi2", &chi2_, "chi2/D", bufsize);
0244     SiPixelLorentzAngleTreeBarrel_->Branch("ndof", &ndof_, "ndof/D", bufsize);
0245     SiPixelLorentzAngleTreeBarrel_->Branch("trackhit", &trackhit_, "x/F:y/F:alpha/D:beta/D:gamma_/D", bufsize);
0246     SiPixelLorentzAngleTreeBarrel_->Branch("npix", &pixinfo_.npix, "npix/I", bufsize);
0247     SiPixelLorentzAngleTreeBarrel_->Branch("rowpix", pixinfo_.row, "row[npix]/F", bufsize);
0248     SiPixelLorentzAngleTreeBarrel_->Branch("colpix", pixinfo_.col, "col[npix]/F", bufsize);
0249     SiPixelLorentzAngleTreeBarrel_->Branch("adc", pixinfo_.adc, "adc[npix]/F", bufsize);
0250     SiPixelLorentzAngleTreeBarrel_->Branch("xpix", pixinfo_.x, "x[npix]/F", bufsize);
0251     SiPixelLorentzAngleTreeBarrel_->Branch("ypix", pixinfo_.y, "y[npix]/F", bufsize);
0252 
0253     SiPixelLorentzAngleTreeBarrel_->Branch(
0254         "clust",
0255         &clust_,
0256         "x/F:y/F:charge/F:size_x/I:size_y/I:maxPixelCol/I:maxPixelRow:minPixelCol/I:minPixelRow/I",
0257         bufsize);
0258     SiPixelLorentzAngleTreeBarrel_->Branch("rechit", &rechit_, "x/F:y/F", bufsize);
0259     SiPixelLorentzAngleTreeBarrel_->Branch("rechit_corr", &rechitCorr_, "x/F:y/F", bufsize);
0260     SiPixelLorentzAngleTreeBarrel_->Branch("trackhitcorr_x", &trackhitCorrX_, "trackhitcorr_x/F", bufsize);
0261     SiPixelLorentzAngleTreeBarrel_->Branch("trackhitcorr_y", &trackhitCorrY_, "trackhitcorr_y/F", bufsize);
0262     SiPixelLorentzAngleTreeBarrel_->Branch("qScale", &qScale_, "qScale/F", bufsize);
0263     SiPixelLorentzAngleTreeBarrel_->Branch("rQmQt", &rQmQt_, "rQmQt/F", bufsize);
0264     //    Forward pixel
0265 
0266     SiPixelLorentzAngleTreeForward_ =
0267         std::make_unique<TTree>("SiPixelLorentzAngleTreeForward_", "SiPixel LorentzAngle tree forward", bufsize);
0268     SiPixelLorentzAngleTreeForward_->Branch("run", &run_, "run/I", bufsize);
0269     SiPixelLorentzAngleTreeForward_->Branch("event", &event_, "event/l", bufsize);
0270     SiPixelLorentzAngleTreeForward_->Branch("lumiblock", &lumiblock_, "lumiblock/I", bufsize);
0271     SiPixelLorentzAngleTreeForward_->Branch("bx", &bx_, "bx/I", bufsize);
0272     SiPixelLorentzAngleTreeForward_->Branch("orbit", &orbit_, "orbit/I", bufsize);
0273     SiPixelLorentzAngleTreeForward_->Branch("side", &sideF_, "side/I", bufsize);
0274     SiPixelLorentzAngleTreeForward_->Branch("disk", &diskF_, "disk/I", bufsize);
0275     SiPixelLorentzAngleTreeForward_->Branch("blade", &bladeF_, "blade/I", bufsize);
0276     SiPixelLorentzAngleTreeForward_->Branch("panel", &panelF_, "panel/I", bufsize);
0277     SiPixelLorentzAngleTreeForward_->Branch("module", &moduleF_, "module/I", bufsize);
0278     SiPixelLorentzAngleTreeForward_->Branch("pt", &pt_, "pt/F", bufsize);
0279     SiPixelLorentzAngleTreeForward_->Branch("eta", &eta_, "eta/F", bufsize);
0280     SiPixelLorentzAngleTreeForward_->Branch("phi", &phi_, "phi/F", bufsize);
0281     SiPixelLorentzAngleTreeForward_->Branch("chi2", &chi2_, "chi2/D", bufsize);
0282     SiPixelLorentzAngleTreeForward_->Branch("ndof", &ndof_, "ndof/D", bufsize);
0283     SiPixelLorentzAngleTreeForward_->Branch("trackhit", &trackhitF_, "x/F:y/F:alpha/D:beta/D:gamma_/D", bufsize);
0284     SiPixelLorentzAngleTreeForward_->Branch("npix", &pixinfoF_.npix, "npix/I", bufsize);
0285     SiPixelLorentzAngleTreeForward_->Branch("rowpix", pixinfoF_.row, "row[npix]/F", bufsize);
0286     SiPixelLorentzAngleTreeForward_->Branch("colpix", pixinfoF_.col, "col[npix]/F", bufsize);
0287     SiPixelLorentzAngleTreeForward_->Branch("adc", pixinfoF_.adc, "adc[npix]/F", bufsize);
0288     SiPixelLorentzAngleTreeForward_->Branch("xpix", pixinfoF_.x, "x[npix]/F", bufsize);
0289     SiPixelLorentzAngleTreeForward_->Branch("ypix", pixinfoF_.y, "y[npix]/F", bufsize);
0290 
0291     SiPixelLorentzAngleTreeForward_->Branch(
0292         "clust",
0293         &clustF_,
0294         "x/F:y/F:charge/F:size_x/I:size_y/I:maxPixelCol/I:maxPixelRow:minPixelCol/I:minPixelRow/I",
0295         bufsize);
0296     SiPixelLorentzAngleTreeForward_->Branch("rechit", &rechitF_, "x/F:y/F", bufsize);
0297     SiPixelLorentzAngleTreeForward_->Branch("rechit_corr", &rechitCorrF_, "x/F:y/F", bufsize);
0298     SiPixelLorentzAngleTreeForward_->Branch("trackhitcorr_x", &trackhitCorrXF_, "trackhitcorr_x/F", bufsize);
0299     SiPixelLorentzAngleTreeForward_->Branch("trackhitcorr_y", &trackhitCorrYF_, "trackhitcorr_y/F", bufsize);
0300     SiPixelLorentzAngleTreeForward_->Branch("qScale", &qScaleF_, "qScale/F", bufsize);
0301     SiPixelLorentzAngleTreeForward_->Branch("rQmQt", &rQmQtF_, "rQmQt/F", bufsize);
0302   }
0303 }
0304 
0305 //
0306 // member functions
0307 //
0308 
0309 // ------------ method called for each event  ------------
0310 
0311 void SiPixelLorentzAnglePCLWorker::analyze(edm::Event const& iEvent, edm::EventSetup const& iSetup) {
0312   // Retrieve tracker topology from geometry
0313   const TrackerTopology* const tTopo = &iSetup.getData(topoPerEventEsToken_);
0314 
0315   // Retrieve track geometry
0316   const TrackerGeometry* tracker = &iSetup.getData(geomPerEventEsToken_);
0317 
0318   // get the association map between tracks and trajectories
0319   edm::Handle<TrajTrackAssociationCollection> trajTrackCollectionHandle;
0320   iEvent.getByToken(t_trajTrack, trajTrackCollectionHandle);
0321 
0322   module_ = -1;
0323   layer_ = -1;
0324   ladder_ = -1;
0325   isflipped_ = -1;
0326   pt_ = -999;
0327   eta_ = 999;
0328   phi_ = 999;
0329   pixinfo_.npix = 0;
0330 
0331   run_ = iEvent.id().run();
0332   event_ = iEvent.id().event();
0333   lumiblock_ = iEvent.luminosityBlock();
0334   bx_ = iEvent.bunchCrossing();
0335   orbit_ = iEvent.orbitNumber();
0336 
0337   if (!trajTrackCollectionHandle->empty()) {
0338     for (TrajTrackAssociationCollection::const_iterator it = trajTrackCollectionHandle->begin();
0339          it != trajTrackCollectionHandle->end();
0340          ++it) {
0341       const reco::Track& track = *it->val;
0342       const Trajectory& traj = *it->key;
0343 
0344       // get the trajectory measurements
0345       std::vector<TrajectoryMeasurement> tmColl = traj.measurements();
0346       pt_ = track.pt();
0347       eta_ = track.eta();
0348       phi_ = track.phi();
0349       chi2_ = traj.chiSquared();
0350       ndof_ = traj.ndof();
0351 
0352       if (pt_ < ptmin_)
0353         continue;
0354 
0355       iHists.h_trackEta_->Fill(eta_);
0356       iHists.h_trackPhi_->Fill(phi_);
0357       iHists.h_trackPt_->Fill(pt_);
0358       iHists.h_trackChi2_->Fill(chi2_ / ndof_);
0359       iHists.h_tracks_->Fill(0);
0360       bool pixeltrack = false;
0361 
0362       // iterate over trajectory measurements
0363       for (const auto& itTraj : tmColl) {
0364         if (!itTraj.updatedState().isValid())
0365           continue;
0366         const TransientTrackingRecHit::ConstRecHitPointer& recHit = itTraj.recHit();
0367         if (!recHit->isValid() || recHit->geographicalId().det() != DetId::Tracker)
0368           continue;
0369         unsigned int subDetID = (recHit->geographicalId().subdetId());
0370         if (subDetID == PixelSubdetector::PixelBarrel || subDetID == PixelSubdetector::PixelEndcap) {
0371           if (!pixeltrack) {
0372             iHists.h_tracks_->Fill(1);
0373           }
0374           pixeltrack = true;
0375         }
0376 
0377         if (subDetID == PixelSubdetector::PixelBarrel) {
0378           DetId detIdObj = recHit->geographicalId();
0379           const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(tracker->idToDet(detIdObj));
0380           if (!theGeomDet)
0381             continue;
0382 
0383           const PixelTopology* topol = &(theGeomDet->specificTopology());
0384 
0385           float ypitch_ = topol->pitch().second;
0386           float width_ = theGeomDet->surface().bounds().thickness();
0387 
0388           if (!topol)
0389             continue;
0390 
0391           layer_ = tTopo->pxbLayer(detIdObj);
0392           ladder_ = tTopo->pxbLadder(detIdObj);
0393           module_ = tTopo->pxbModule(detIdObj);
0394 
0395           float tmp1 = theGeomDet->surface().toGlobal(Local3DPoint(0., 0., 0.)).perp();
0396           float tmp2 = theGeomDet->surface().toGlobal(Local3DPoint(0., 0., 1.)).perp();
0397 
0398           isflipped_ = (tmp2 < tmp1) ? 1 : 0;
0399 
0400           const SiPixelRecHit* recHitPix = dynamic_cast<const SiPixelRecHit*>((*recHit).hit());
0401           if (!recHitPix)
0402             continue;
0403           rechit_.x = recHitPix->localPosition().x();
0404           rechit_.y = recHitPix->localPosition().y();
0405           SiPixelRecHit::ClusterRef const& cluster = recHitPix->cluster();
0406 
0407           pixinfo_ = fillPix(*cluster, topol);
0408 
0409           // fill entries in clust_
0410 
0411           clust_.x = (cluster)->x();
0412           clust_.y = (cluster)->y();
0413           clust_.charge = (cluster->charge()) / 1000.;  // clust_.charge: in the unit of 1000e
0414           clust_.size_x = cluster->sizeX();
0415           clust_.size_y = cluster->sizeY();
0416           clust_.maxPixelCol = cluster->maxPixelCol();
0417           clust_.maxPixelRow = cluster->maxPixelRow();
0418           clust_.minPixelCol = cluster->minPixelCol();
0419           clust_.minPixelRow = cluster->minPixelRow();
0420 
0421           // fill the trackhit info
0422           TrajectoryStateOnSurface tsos = itTraj.updatedState();
0423           if (!tsos.isValid()) {
0424             edm::LogWarning("SiPixelLorentzAnglePCLWorker") << "tsos not valid";
0425             continue;
0426           }
0427           LocalVector trackdirection = tsos.localDirection();
0428           LocalPoint trackposition = tsos.localPosition();
0429 
0430           if (trackdirection.z() == 0)
0431             continue;
0432           // the local position and direction
0433           trackhit_.alpha = atan2(trackdirection.z(), trackdirection.x());
0434           trackhit_.beta = atan2(trackdirection.z(), trackdirection.y());
0435           trackhit_.gamma = atan2(trackdirection.x(), trackdirection.y());
0436           trackhit_.x = trackposition.x();
0437           trackhit_.y = trackposition.y();
0438 
0439           // get qScale_ = templ.qscale() and  templ.r_qMeas_qTrue();
0440           float cotalpha = trackdirection.x() / trackdirection.z();
0441           float cotbeta = trackdirection.y() / trackdirection.z();
0442           float cotbeta_min = clustSizeYMin_[layer_ - 1] * ypitch_ / width_;
0443           if (fabs(cotbeta) <= cotbeta_min)
0444             continue;
0445           double drdz = sqrt(1. + cotalpha * cotalpha + cotbeta * cotbeta);
0446           double clusterCharge_cut = clustChargeMaxPerLength_ * drdz;
0447 
0448           auto detId = detIdObj.rawId();
0449           int DetId_index = -1;
0450 
0451           const auto& newModIt = (std::find(iHists.BPixnewDetIds_.begin(), iHists.BPixnewDetIds_.end(), detId));
0452           bool isNewMod = (newModIt != iHists.BPixnewDetIds_.end());
0453           if (isNewMod) {
0454             DetId_index = std::distance(iHists.BPixnewDetIds_.begin(), newModIt);
0455           }
0456 
0457           if (notInPCL_) {
0458             // fill the template from the store (from dqmBeginRun)
0459             SiPixelTemplate theTemplate(thePixelTemp_);
0460 
0461             float locBx = (cotbeta < 0.) ? -1 : 1.;
0462             float locBz = (cotalpha < 0.) ? -locBx : locBx;
0463 
0464             int TemplID = templateDBobject_->getTemplateID(detId);
0465             theTemplate.interpolate(TemplID, cotalpha, cotbeta, locBz, locBx);
0466             qScale_ = theTemplate.qscale();
0467             rQmQt_ = theTemplate.r_qMeas_qTrue();
0468           }
0469 
0470           // Surface deformation
0471           const auto& lp_pair = surface_deformation(topol, tsos, recHitPix);
0472 
0473           LocalPoint lp_track = lp_pair.first;
0474           LocalPoint lp_rechit = lp_pair.second;
0475 
0476           rechitCorr_.x = lp_rechit.x();
0477           rechitCorr_.y = lp_rechit.y();
0478           trackhitCorrX_ = lp_track.x();
0479           trackhitCorrY_ = lp_track.y();
0480 
0481           if (notInPCL_) {
0482             SiPixelLorentzAngleTreeBarrel_->Fill();
0483           }
0484 
0485           // is one pixel in cluster a large pixel ? (hit will be excluded)
0486           bool large_pix = false;
0487           for (int j = 0; j < pixinfo_.npix; j++) {
0488             int colpos = static_cast<int>(pixinfo_.col[j]);
0489             if (pixinfo_.row[j] == 0 || pixinfo_.row[j] == 79 || pixinfo_.row[j] == 80 || pixinfo_.row[j] == 159 ||
0490                 colpos % 52 == 0 || colpos % 52 == 51) {
0491               large_pix = true;
0492             }
0493           }
0494 
0495           double residualsq = (trackhitCorrX_ - rechitCorr_.x) * (trackhitCorrX_ - rechitCorr_.x) +
0496                               (trackhitCorrY_ - rechitCorr_.y) * (trackhitCorrY_ - rechitCorr_.y);
0497 
0498           double xlim1 = trackhitCorrX_ - width_ * cotalpha / 2.;
0499           double hypitch_ = ypitch_ / 2.;
0500           double ylim1 = trackhitCorrY_ - width_ * cotbeta / 2.;
0501           double ylim2 = trackhitCorrY_ + width_ * cotbeta / 2.;
0502 
0503           int clustSizeY_cut = clustSizeYMin_[layer_ - 1];
0504 
0505           if (!large_pix && (chi2_ / ndof_) < normChi2Max_ && cluster->sizeY() >= clustSizeY_cut &&
0506               residualsq < residualMax_ * residualMax_ && cluster->charge() < clusterCharge_cut &&
0507               cluster->sizeX() < clustSizeXMax_) {
0508             // iterate over pixels in hit
0509             for (int j = 0; j < pixinfo_.npix; j++) {
0510               // use trackhits and include bowing correction
0511               float ypixlow = pixinfo_.y[j] - hypitch_;
0512               float ypixhigh = pixinfo_.y[j] + hypitch_;
0513               if (cotbeta > 0.) {
0514                 if (ylim1 > ypixlow)
0515                   ypixlow = ylim1;
0516                 if (ylim2 < ypixhigh)
0517                   ypixhigh = ylim2;
0518               } else {
0519                 if (ylim2 > ypixlow)
0520                   ypixlow = ylim2;
0521                 if (ylim1 < ypixhigh)
0522                   ypixhigh = ylim1;
0523               }
0524               float ypixavg = 0.5f * (ypixlow + ypixhigh);
0525 
0526               float dx = (pixinfo_.x[j] - xlim1) * siPixelLACalibration::cmToum;  // dx: in the unit of micrometer
0527               float dy = (ypixavg - ylim1) * siPixelLACalibration::cmToum;        // dy: in the unit of micrometer
0528               float depth = dy * tan(trackhit_.beta);
0529               float drift = dx - dy * tan(trackhit_.gamma);
0530 
0531               if (isNewMod == false) {
0532                 int i_index = module_ + (layer_ - 1) * iHists.nModules_[layer_ - 1];
0533                 iHists.h_drift_depth_adc_[i_index]->Fill(drift, depth, pixinfo_.adc[j]);
0534                 iHists.h_drift_depth_adc2_[i_index]->Fill(drift, depth, pixinfo_.adc[j] * pixinfo_.adc[j]);
0535                 iHists.h_drift_depth_noadc_[i_index]->Fill(drift, depth, 1.);
0536                 iHists.h_bySectOccupancy_->Fill(i_index - 1);  // histogram starts at 0
0537 
0538                 if (tracker->getDetectorType(subDetID) == TrackerGeometry::ModuleType::Ph1PXB) {
0539                   if ((module_ == 3 || module_ == 5) && (layer_ == 3 || layer_ == 4)) {
0540                     int i_index_merge = i_index + 1;
0541                     iHists.h_drift_depth_adc_[i_index_merge]->Fill(drift, depth, pixinfo_.adc[j]);
0542                     iHists.h_drift_depth_adc2_[i_index_merge]->Fill(drift, depth, pixinfo_.adc[j] * pixinfo_.adc[j]);
0543                     iHists.h_drift_depth_noadc_[i_index_merge]->Fill(drift, depth, 1.);
0544                     iHists.h_bySectOccupancy_->Fill(i_index_merge - 1);
0545                   }
0546                   if ((module_ == 4 || module_ == 6) && (layer_ == 3 || layer_ == 4)) {
0547                     int i_index_merge = i_index - 1;
0548                     iHists.h_drift_depth_adc_[i_index_merge]->Fill(drift, depth, pixinfo_.adc[j]);
0549                     iHists.h_drift_depth_adc2_[i_index_merge]->Fill(drift, depth, pixinfo_.adc[j] * pixinfo_.adc[j]);
0550                     iHists.h_drift_depth_noadc_[i_index_merge]->Fill(drift, depth, 1.);
0551                     iHists.h_bySectOccupancy_->Fill(i_index_merge - 1);
0552                   }
0553                 }
0554 
0555               } else {
0556                 int new_index = iHists.nModules_[iHists.nlay - 1] +
0557                                 (iHists.nlay - 1) * iHists.nModules_[iHists.nlay - 1] + 1 + DetId_index;
0558 
0559                 iHists.h_drift_depth_adc_[new_index]->Fill(drift, depth, pixinfo_.adc[j]);
0560                 iHists.h_drift_depth_adc2_[new_index]->Fill(drift, depth, pixinfo_.adc[j] * pixinfo_.adc[j]);
0561                 iHists.h_drift_depth_noadc_[new_index]->Fill(drift, depth, 1.);
0562                 iHists.h_bySectOccupancy_->Fill(new_index - 1);  // histogram starts at 0
0563               }
0564             }
0565           }
0566         } else if (subDetID == PixelSubdetector::PixelEndcap) {
0567           DetId detIdObj = recHit->geographicalId();
0568           const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(tracker->idToDet(detIdObj));
0569           if (!theGeomDet)
0570             continue;
0571 
0572           const PixelTopology* topol = &(theGeomDet->specificTopology());
0573 
0574           if (!topol)
0575             continue;
0576 
0577           sideF_ = tTopo->pxfSide(detIdObj);
0578           diskF_ = tTopo->pxfDisk(detIdObj);
0579           bladeF_ = tTopo->pxfBlade(detIdObj);
0580           panelF_ = tTopo->pxfPanel(detIdObj);
0581           moduleF_ = tTopo->pxfModule(detIdObj);
0582 
0583           const SiPixelRecHit* recHitPix = dynamic_cast<const SiPixelRecHit*>((*recHit).hit());
0584           if (!recHitPix)
0585             continue;
0586           rechitF_.x = recHitPix->localPosition().x();
0587           rechitF_.y = recHitPix->localPosition().y();
0588           SiPixelRecHit::ClusterRef const& cluster = recHitPix->cluster();
0589 
0590           pixinfoF_ = fillPix(*cluster, topol);
0591 
0592           // fill entries in clust_
0593 
0594           clustF_.x = (cluster)->x();
0595           clustF_.y = (cluster)->y();
0596           clustF_.charge = (cluster->charge()) / 1000.;  // clustF_.charge: in the unit of 1000e
0597           clustF_.size_x = cluster->sizeX();
0598           clustF_.size_y = cluster->sizeY();
0599           clustF_.maxPixelCol = cluster->maxPixelCol();
0600           clustF_.maxPixelRow = cluster->maxPixelRow();
0601           clustF_.minPixelCol = cluster->minPixelCol();
0602           clustF_.minPixelRow = cluster->minPixelRow();
0603 
0604           // fill the trackhit info
0605           TrajectoryStateOnSurface tsos = itTraj.updatedState();
0606           if (!tsos.isValid()) {
0607             edm::LogWarning("SiPixelLorentzAnglePCLWorker") << "tsos not valid";
0608             continue;
0609           }
0610           LocalVector trackdirection = tsos.localDirection();
0611           LocalPoint trackposition = tsos.localPosition();
0612 
0613           if (trackdirection.z() == 0)
0614             continue;
0615           // the local position and direction
0616           trackhitF_.alpha = atan2(trackdirection.z(), trackdirection.x());
0617           trackhitF_.beta = atan2(trackdirection.z(), trackdirection.y());
0618           trackhitF_.gamma = atan2(trackdirection.x(), trackdirection.y());
0619           trackhitF_.x = trackposition.x();
0620           trackhitF_.y = trackposition.y();
0621 
0622           float cotalpha = trackdirection.x() / trackdirection.z();
0623           float cotbeta = trackdirection.y() / trackdirection.z();
0624 
0625           auto detId = detIdObj.rawId();
0626 
0627           if (notInPCL_) {
0628             // fill the template from the store (from dqmBeginRun)
0629             SiPixelTemplate theTemplate(thePixelTemp_);
0630 
0631             float locBx = (cotbeta < 0.) ? -1 : 1.;
0632             float locBz = (cotalpha < 0.) ? -locBx : locBx;
0633 
0634             int TemplID = templateDBobject_->getTemplateID(detId);
0635             theTemplate.interpolate(TemplID, cotalpha, cotbeta, locBz, locBx);
0636             qScaleF_ = theTemplate.qscale();
0637             rQmQtF_ = theTemplate.r_qMeas_qTrue();
0638           }
0639 
0640           // Surface deformation
0641           const auto& lp_pair = surface_deformation(topol, tsos, recHitPix);
0642 
0643           LocalPoint lp_track = lp_pair.first;
0644           LocalPoint lp_rechit = lp_pair.second;
0645 
0646           rechitCorrF_.x = lp_rechit.x();
0647           rechitCorrF_.y = lp_rechit.y();
0648           trackhitCorrXF_ = lp_track.x();
0649           trackhitCorrYF_ = lp_track.y();
0650           if (notInPCL_) {
0651             SiPixelLorentzAngleTreeForward_->Fill();
0652           }
0653         }
0654       }  //end iteration over trajectory measurements
0655     }    //end iteration over trajectories
0656   }
0657 }
0658 
0659 void SiPixelLorentzAnglePCLWorker::dqmBeginRun(edm::Run const& run, edm::EventSetup const& iSetup) {
0660   // geometry
0661   const TrackerGeometry* geom = &iSetup.getData(geomEsToken_);
0662   const TrackerTopology* tTopo = &iSetup.getData(topoEsToken_);
0663 
0664   if (notInPCL_) {
0665     // Initialize 1D templates
0666     if (watchSiPixelTemplateRcd_.check(iSetup)) {
0667       templateDBobject_ = &iSetup.getData(siPixelTemplateEsToken_);
0668       if (!SiPixelTemplate::pushfile(*templateDBobject_, thePixelTemp_)) {
0669         edm::LogError("SiPixelLorentzAnglePCLWorker")
0670             << "Templates not filled correctly. Check the sqlite file. Using SiPixelTemplateDBObject version "
0671             << (*templateDBobject_).version();
0672       }
0673     }
0674   }
0675 
0676   PixelTopologyMap map = PixelTopologyMap(geom, tTopo);
0677   iHists.nlay = geom->numberOfLayers(PixelSubdetector::PixelBarrel);
0678   iHists.nModules_.resize(iHists.nlay);
0679   for (int i = 0; i < iHists.nlay; i++) {
0680     iHists.nModules_[i] = map.getPXBModules(i + 1);
0681   }
0682 
0683   // list of modules already filled, then return (we already entered here)
0684   if (!iHists.BPixnewDetIds_.empty() || !iHists.FPixnewDetIds_.empty())
0685     return;
0686 
0687   if (!newmodulelist_.empty()) {
0688     for (auto const& modulename : newmodulelist_) {
0689       if (modulename.find("BPix_") != std::string::npos) {
0690         PixelBarrelName bn(modulename, true);
0691         const auto& detId = bn.getDetId(tTopo);
0692         iHists.BPixnewmodulename_.push_back(modulename);
0693         iHists.BPixnewDetIds_.push_back(detId.rawId());
0694         iHists.BPixnewModule_.push_back(bn.moduleName());
0695         iHists.BPixnewLayer_.push_back(bn.layerName());
0696       } else if (modulename.find("FPix_") != std::string::npos) {
0697         PixelEndcapName en(modulename, true);
0698         const auto& detId = en.getDetId(tTopo);
0699         iHists.FPixnewmodulename_.push_back(modulename);
0700         iHists.FPixnewDetIds_.push_back(detId.rawId());
0701         iHists.FPixnewDisk_.push_back(en.diskName());
0702         iHists.FPixnewBlade_.push_back(en.bladeName());
0703       }
0704     }
0705   }
0706 }
0707 
0708 void SiPixelLorentzAnglePCLWorker::bookHistograms(DQMStore::IBooker& iBooker,
0709                                                   edm::Run const& run,
0710                                                   edm::EventSetup const& iSetup) {
0711   // book the by partition monitoring
0712   const auto maxSect = iHists.nlay * iHists.nModules_[iHists.nlay - 1] + (int)iHists.BPixnewDetIds_.size();
0713 
0714   iBooker.setCurrentFolder(fmt::sprintf("%s/SectorMonitoring", folder_.data()));
0715   iHists.h_bySectOccupancy_ = iBooker.book1D(
0716       "h_bySectorOccupancy", "hit occupancy by sector;pixel sector;hits on track", maxSect, -0.5, maxSect - 0.5);
0717 
0718   iBooker.setCurrentFolder(folder_);
0719   static constexpr double min_depth_ = -100.;
0720   static constexpr double max_depth_ = 400.;
0721   static constexpr double min_drift_ = -500.;
0722   static constexpr double max_drift_ = 500.;
0723 
0724   // book the mean values projections and set the bin names of the by sector monitoring
0725   char name[128];
0726   char title[256];
0727   for (int i_layer = 1; i_layer <= iHists.nlay; i_layer++) {
0728     for (int i_module = 1; i_module <= iHists.nModules_[i_layer - 1]; i_module++) {
0729       unsigned int i_index = i_module + (i_layer - 1) * iHists.nModules_[i_layer - 1];
0730       std::string binName = fmt::sprintf("BPix Lay%i Mod%i", i_layer, i_module);
0731       LogDebug("SiPixelLorentzAnglePCLWorker") << " i_index: " << i_index << " bin name: " << binName
0732                                                << " (i_layer: " << i_layer << " i_module:" << i_module << ")";
0733 
0734       iHists.h_bySectOccupancy_->setBinLabel(i_index, binName);
0735 
0736       sprintf(name, "h_mean_layer%i_module%i", i_layer, i_module);
0737       sprintf(title,
0738               "average drift vs depth layer%i module%i; production depth [#mum]; #LTdrift#GT [#mum]",
0739               i_layer,
0740               i_module);
0741       iHists.h_mean_[i_index] = iBooker.book1D(name, title, hist_depth_, min_depth_, max_depth_);
0742     }
0743   }
0744   for (int i = 0; i < (int)iHists.BPixnewDetIds_.size(); i++) {
0745     sprintf(name, "h_BPixnew_mean_%s", iHists.BPixnewmodulename_[i].c_str());
0746     sprintf(title,
0747             "average drift vs depth %s; production depth [#mum]; #LTdrift#GT [#mum]",
0748             iHists.BPixnewmodulename_[i].c_str());
0749     int new_index = iHists.nModules_[iHists.nlay - 1] + (iHists.nlay - 1) * iHists.nModules_[iHists.nlay - 1] + 1 + i;
0750     iHists.h_mean_[new_index] = iBooker.book1D(name, title, hist_depth_, min_depth_, max_depth_);
0751 
0752     LogDebug("SiPixelLorentzAnglePCLWorker") << "i_index" << new_index << " bin name: " << iHists.BPixnewmodulename_[i];
0753 
0754     iHists.h_bySectOccupancy_->setBinLabel(new_index, iHists.BPixnewmodulename_[i]);
0755   }
0756 
0757   //book the 2D histograms
0758   for (int i_layer = 1; i_layer <= iHists.nlay; i_layer++) {
0759     iBooker.setCurrentFolder(fmt::sprintf("%s/BPix/BPixLayer%i", folder_.data(), i_layer));
0760     for (int i_module = 1; i_module <= iHists.nModules_[i_layer - 1]; i_module++) {
0761       unsigned int i_index = i_module + (i_layer - 1) * iHists.nModules_[i_layer - 1];
0762 
0763       sprintf(name, "h_drift_depth_adc_layer%i_module%i", i_layer, i_module);
0764       sprintf(title, "depth vs drift (ADC) layer%i module%i; drift [#mum]; production depth [#mum]", i_layer, i_module);
0765       iHists.h_drift_depth_adc_[i_index] =
0766           iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
0767 
0768       sprintf(name, "h_drift_depth_adc2_layer%i_module%i", i_layer, i_module);
0769       sprintf(
0770           title, "depth vs drift (ADC^{2}) layer%i module%i; drift [#mum]; production depth [#mum]", i_layer, i_module);
0771       iHists.h_drift_depth_adc2_[i_index] =
0772           iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
0773 
0774       sprintf(name, "h_drift_depth_noadc_layer%i_module%i", i_layer, i_module);
0775       sprintf(
0776           title, "depth vs drift (no ADC) layer%i module%i; drift [#mum]; production depth [#mum]", i_layer, i_module);
0777       iHists.h_drift_depth_noadc_[i_index] =
0778           iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
0779 
0780       sprintf(name, "h_drift_depth_layer%i_module%i", i_layer, i_module);
0781       sprintf(title, "depth vs drift layer%i module%i; drift [#mum]; production depth [#mum]", i_layer, i_module);
0782       iHists.h_drift_depth_[i_index] =
0783           iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
0784     }
0785   }
0786 
0787   // book the "new" modules
0788   iBooker.setCurrentFolder(fmt::sprintf("%s/BPix/NewModules", folder_.data()));
0789   for (int i = 0; i < (int)iHists.BPixnewDetIds_.size(); i++) {
0790     int new_index = iHists.nModules_[iHists.nlay - 1] + (iHists.nlay - 1) * iHists.nModules_[iHists.nlay - 1] + 1 + i;
0791 
0792     sprintf(name, "h_BPixnew_drift_depth_adc_%s", iHists.BPixnewmodulename_[i].c_str());
0793     sprintf(
0794         title, "depth vs drift (ADC) %s; drift [#mum]; production depth [#mum]", iHists.BPixnewmodulename_[i].c_str());
0795     iHists.h_drift_depth_adc_[new_index] =
0796         iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
0797 
0798     sprintf(name, "h_BPixnew_drift_depth_adc2_%s", iHists.BPixnewmodulename_[i].c_str());
0799     sprintf(title,
0800             "depth vs drift (ADC^{2}) %s; drift [#mum]; production depth [#mum]",
0801             iHists.BPixnewmodulename_[i].c_str());
0802     iHists.h_drift_depth_adc2_[new_index] =
0803         iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
0804 
0805     sprintf(name, "h_BPixnew_drift_depth_noadc_%s", iHists.BPixnewmodulename_[i].c_str());
0806     sprintf(title,
0807             "depth vs drift (no ADC)%s; drift [#mum]; production depth [#mum]",
0808             iHists.BPixnewmodulename_[i].c_str());
0809     iHists.h_drift_depth_noadc_[new_index] =
0810         iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
0811 
0812     sprintf(name, "h_BPixnew_drift_depth_%s", iHists.BPixnewmodulename_[i].c_str());
0813     sprintf(title, "depth vs drift %s; drift [#mum]; production depth [#mum]", iHists.BPixnewmodulename_[i].c_str());
0814     iHists.h_drift_depth_[new_index] =
0815         iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
0816   }
0817 
0818   // book the track monitoring plots
0819   iBooker.setCurrentFolder(fmt::sprintf("%s/TrackMonitoring", folder_.data()));
0820   iHists.h_tracks_ = iBooker.book1D("h_tracks", ";tracker volume;tracks", 2, -0.5, 1.5);
0821   iHists.h_tracks_->setBinLabel(1, "all tracks", 1);
0822   iHists.h_tracks_->setBinLabel(2, "has pixel hits", 1);
0823   iHists.h_trackEta_ = iBooker.book1D("h_trackEta", ";track #eta; #tracks", 30, -3., 3.);
0824   iHists.h_trackPhi_ = iBooker.book1D("h_trackPhi", ";track #phi; #tracks", 48, -M_PI, M_PI);
0825   iHists.h_trackPt_ = iBooker.book1D("h_trackPt", ";track p_{T} [GeV]; #tracks", 100, 0., 100.);
0826   iHists.h_trackChi2_ = iBooker.book1D("h_trackChi2ndof", ";track #chi^{2}/ndof; #tracks", 100, 0., 10.);
0827 }
0828 
0829 void SiPixelLorentzAnglePCLWorker::dqmEndRun(edm::Run const& run, edm::EventSetup const& iSetup) {
0830   if (notInPCL_) {
0831     hFile_->cd();
0832     hFile_->Write();
0833     hFile_->Close();
0834   }
0835 }
0836 
0837 // method used to fill per pixel info
0838 const Pixinfo SiPixelLorentzAnglePCLWorker::fillPix(const SiPixelCluster& LocPix, const PixelTopology* topol) const {
0839   Pixinfo pixinfo;
0840   const std::vector<SiPixelCluster::Pixel>& pixvector = LocPix.pixels();
0841   pixinfo.npix = 0;
0842   for (std::vector<SiPixelCluster::Pixel>::const_iterator itPix = pixvector.begin(); itPix != pixvector.end();
0843        itPix++) {
0844     pixinfo.row[pixinfo.npix] = itPix->x;
0845     pixinfo.col[pixinfo.npix] = itPix->y;
0846     pixinfo.adc[pixinfo.npix] = itPix->adc;
0847     LocalPoint lp = topol->localPosition(MeasurementPoint(itPix->x + 0.5, itPix->y + 0.5));
0848     pixinfo.x[pixinfo.npix] = lp.x();
0849     pixinfo.y[pixinfo.npix] = lp.y();
0850     pixinfo.npix++;
0851   }
0852   return pixinfo;
0853 }
0854 
0855 // method used to correct for the surface deformation
0856 const std::pair<LocalPoint, LocalPoint> SiPixelLorentzAnglePCLWorker::surface_deformation(
0857     const PixelTopology* topol, TrajectoryStateOnSurface& tsos, const SiPixelRecHit* recHitPix) const {
0858   LocalPoint trackposition = tsos.localPosition();
0859   const LocalTrajectoryParameters& ltp = tsos.localParameters();
0860   const Topology::LocalTrackAngles localTrackAngles(ltp.dxdz(), ltp.dydz());
0861 
0862   std::pair<float, float> pixels_track = topol->pixel(trackposition, localTrackAngles);
0863   std::pair<float, float> pixels_rechit = topol->pixel(recHitPix->localPosition(), localTrackAngles);
0864 
0865   LocalPoint lp_track = topol->localPosition(MeasurementPoint(pixels_track.first, pixels_track.second));
0866 
0867   LocalPoint lp_rechit = topol->localPosition(MeasurementPoint(pixels_rechit.first, pixels_rechit.second));
0868 
0869   std::pair<LocalPoint, LocalPoint> lps = std::make_pair(lp_track, lp_rechit);
0870   return lps;
0871 }
0872 
0873 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0874 void SiPixelLorentzAnglePCLWorker::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0875   edm::ParameterSetDescription desc;
0876   desc.setComment("Worker module of the SiPixel Lorentz Angle PCL monitoring workflow");
0877   desc.add<std::string>("folder", "AlCaReco/SiPixelLorentzAngle")->setComment("directory of PCL Worker output");
0878   desc.add<bool>("notInPCL", false)->setComment("create TTree (true) or not (false)");
0879   desc.add<std::string>("fileName", "testrun.root")->setComment("name of the TTree file if notInPCL = true");
0880   desc.add<std::vector<std::string>>("newmodulelist", {})->setComment("the list of DetIds for new sensors");
0881   desc.add<edm::InputTag>("src", edm::InputTag("TrackRefitter"))->setComment("input track collections");
0882   desc.add<double>("ptMin", 3.)->setComment("minimum pt on tracks");
0883   desc.add<double>("normChi2Max", 2.)->setComment("maximum reduced chi squared");
0884   desc.add<std::vector<int>>("clustSizeYMin", {4, 3, 3, 2})
0885       ->setComment("minimum cluster size on Y axis for all Barrel Layers");
0886   desc.add<int>("clustSizeXMax", 5)->setComment("maximum cluster size on X axis");
0887   desc.add<double>("residualMax", 0.005)->setComment("maximum residual");
0888   desc.add<double>("clustChargeMaxPerLength", 50000)
0889       ->setComment("maximum cluster charge per unit length of pixel depth (z)");
0890   desc.add<int>("binsDepth", 50)->setComment("# bins for electron production depth axis");
0891   desc.add<int>("binsDrift", 100)->setComment("# bins for electron drift axis");
0892   descriptions.addWithDefaultLabel(desc);
0893 }
0894 
0895 // define this as a plug-in
0896 DEFINE_FWK_MODULE(SiPixelLorentzAnglePCLWorker);