CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
SiPixelLorentzAnglePCLWorker Class Reference

#include <CalibTracker/SiPixelLorentzAnglePCLWorker/src/SiPixelLorentzAnglePCLWorker.cc>

Inheritance diagram for SiPixelLorentzAnglePCLWorker:
DQMEDAnalyzer edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >

Public Member Functions

 SiPixelLorentzAnglePCLWorker (const edm::ParameterSet &)
 
 ~SiPixelLorentzAnglePCLWorker () override=default
 
- Public Member Functions inherited from DQMEDAnalyzer
void accumulate (edm::Event const &event, edm::EventSetup const &setup) final
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
void beginStream (edm::StreamID id) final
 
 DQMEDAnalyzer ()
 
void endLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &run, edm::EventSetup const &setup) final
 
virtual bool getCanSaveByLumi ()
 
- Public Member Functions inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from DQMEDAnalyzer
static void globalEndJob (DQMEDAnalyzerGlobalCache const *)
 
static void globalEndLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup, LuminosityBlockContext const *context)
 
static void globalEndRunProduce (edm::Run &run, edm::EventSetup const &setup, RunContext const *context)
 
static std::unique_ptr< DQMEDAnalyzerGlobalCacheinitializeGlobalCache (edm::ParameterSet const &)
 

Private Member Functions

void analyze (edm::Event const &, edm::EventSetup const &) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
LorentzAngleAnalysisTypeEnum convertStringToLorentzAngleAnalysisTypeEnum (std::string type)
 
void dqmBeginRun (edm::Run const &, edm::EventSetup const &) override
 
void dqmEndRun (edm::Run const &, edm::EventSetup const &)
 
const Pixinfo fillPix (const SiPixelCluster &LocPix, const PixelTopology *topol) const
 
const std::pair< LocalPoint, LocalPointsurface_deformation (const PixelTopology *topol, TrajectoryStateOnSurface &tsos, const SiPixelRecHit *recHitPix) const
 

Private Attributes

LorentzAngleAnalysisTypeEnum analysisType_
 
int bladeF_
 
int bx_
 
double chi2_
 
Clust clust_
 
double clustChargeMaxPerLength_
 
Clust clustF_
 
int clustSizeXMax_
 
std::vector< int > clustSizeYMin_
 
int diskF_
 
float eta_
 
long int event_
 
std::string filename_
 
std::string folder_
 
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecordgeomEsToken_
 
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecordgeomPerEventEsToken_
 
std::unique_ptr< TFile > hFile_
 
int hist_depth_
 
int hist_drift_
 
SiPixelLorentzAngleCalibrationHistograms iHists
 
int isflipped_
 
int ladder_
 
int layer_
 
int lumiblock_
 
edm::ESGetToken< MagneticField, IdealMagneticFieldRecordmagneticFieldToken_
 
int module_
 
int moduleF_
 
double ndof_
 
std::vector< std::string > newmodulelist_
 
double normChi2Max_
 
bool notInPCL_
 
int orbit_
 
int panelF_
 
float phi_
 
Pixinfo pixinfo_
 
Pixinfo pixinfoF_
 
float pt_
 
double ptmin_
 
float qScale_
 
float qScaleF_
 
Rechit rechit_
 
Rechit rechitCorr_
 
Rechit rechitCorrF_
 
Rechit rechitF_
 
double residualMax_
 
float rQmQt_
 
float rQmQtF_
 
int run_
 
int sideF_
 
Hit simhit_
 
Hit simhitF_
 
std::unique_ptr< TTree > SiPixelLorentzAngleTreeBarrel_
 
std::unique_ptr< TTree > SiPixelLorentzAngleTreeForward_
 
edm::ESGetToken< SiPixelTemplateDBObject, SiPixelTemplateDBObjectESProducerRcdsiPixelTemplateEsToken_
 
edm::ESGetToken< std::vector< SiPixelTemplateStore >, SiPixelTemplateDBObjectESProducerRcdsiPixelTemplateStoreEsToken_
 
edm::EDGetTokenT< TrajTrackAssociationCollectiont_trajTrack
 
const SiPixelTemplateDBObjecttemplateDBobject_
 
const std::vector< SiPixelTemplateStore > * thePixelTemp_ = nullptr
 
edm::ESGetToken< TrackerTopology, TrackerTopologyRcdtopoEsToken_
 
edm::ESGetToken< TrackerTopology, TrackerTopologyRcdtopoPerEventEsToken_
 
Hit trackhit_
 
float trackhitCorrX_
 
float trackhitCorrXF_
 
float trackhitCorrY_
 
float trackhitCorrYF_
 
Hit trackhitF_
 
edm::ESWatcher< SiPixelTemplateDBObjectESProducerRcdwatchSiPixelTemplateRcd_
 

Additional Inherited Members

- Public Types inherited from DQMEDAnalyzer
typedef dqm::reco::DQMStore DQMStore
 
typedef dqm::reco::MonitorElement MonitorElement
 
- Public Types inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 
- Protected Member Functions inherited from DQMEDAnalyzer
uint64_t meId () const
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 
unsigned int streamId_
 

Detailed Description

Description: generates the intermediate ALCAPROMPT dataset for the measurement of the SiPixel Lorentz Angle in the Prompt Calibration Loop Implementation: Books and fills 2D histograms of the drift vs depth in bins of pixel module rings to be fed into the SiPixelLorentzAnglePCLHarvester

Definition at line 104 of file SiPixelLorentzAnglePCLWorker.cc.

Constructor & Destructor Documentation

◆ SiPixelLorentzAnglePCLWorker()

SiPixelLorentzAnglePCLWorker::SiPixelLorentzAnglePCLWorker ( const edm::ParameterSet iConfig)
explicit

Definition at line 210 of file SiPixelLorentzAnglePCLWorker.cc.

References Pixinfo::adc, bladeF_, bx_, chi2_, clust_, clustF_, Pixinfo::col, diskF_, eta_, event_, filename_, edm::ParameterSet::getParameter(), hFile_, isflipped_, ladder_, layer_, lumiblock_, module_, moduleF_, ndof_, notInPCL_, Pixinfo::npix, orbit_, panelF_, phi_, pixinfo_, pixinfoF_, pt_, qScale_, qScaleF_, rechit_, rechitCorr_, rechitCorrF_, rechitF_, Pixinfo::row, rQmQt_, rQmQtF_, run_, sideF_, SiPixelLorentzAngleTreeBarrel_, SiPixelLorentzAngleTreeForward_, t_trajTrack, trackhit_, trackhitCorrX_, trackhitCorrXF_, trackhitCorrY_, trackhitCorrYF_, trackhitF_, Pixinfo::x, and Pixinfo::y.

212  folder_(iConfig.getParameter<std::string>("folder")),
213  notInPCL_(iConfig.getParameter<bool>("notInPCL")),
214  filename_(iConfig.getParameter<std::string>("fileName")),
215  newmodulelist_(iConfig.getParameter<std::vector<std::string>>("newmodulelist")),
216  ptmin_(iConfig.getParameter<double>("ptMin")),
217  normChi2Max_(iConfig.getParameter<double>("normChi2Max")),
218  clustSizeYMin_(iConfig.getParameter<std::vector<int>>("clustSizeYMin")),
219  clustSizeXMax_(iConfig.getParameter<int>("clustSizeXMax")),
220  residualMax_(iConfig.getParameter<double>("residualMax")),
221  clustChargeMaxPerLength_(iConfig.getParameter<double>("clustChargeMaxPerLength")),
222  hist_depth_(iConfig.getParameter<int>("binsDepth")),
223  hist_drift_(iConfig.getParameter<int>("binsDrift")),
224  geomEsToken_(esConsumes<edm::Transition::BeginRun>()),
225  topoEsToken_(esConsumes<edm::Transition::BeginRun>()),
226  siPixelTemplateEsToken_(esConsumes<edm::Transition::BeginRun>()),
227  siPixelTemplateStoreEsToken_(esConsumes<edm::Transition::BeginRun>()),
231  t_trajTrack = consumes<TrajTrackAssociationCollection>(iConfig.getParameter<edm::InputTag>("src"));
232 
233  // now do what ever initialization is needed
234  int bufsize = 64000;
235  // create tree structure
236  // Barrel pixel
237  if (notInPCL_) {
238  hFile_ = std::make_unique<TFile>(filename_.c_str(), "RECREATE");
240  std::make_unique<TTree>("SiPixelLorentzAngleTreeBarrel_", "SiPixel LorentzAngle tree barrel", bufsize);
241  SiPixelLorentzAngleTreeBarrel_->Branch("run", &run_, "run/I", bufsize);
242  SiPixelLorentzAngleTreeBarrel_->Branch("event", &event_, "event/l", bufsize);
243  SiPixelLorentzAngleTreeBarrel_->Branch("lumiblock", &lumiblock_, "lumiblock/I", bufsize);
244  SiPixelLorentzAngleTreeBarrel_->Branch("bx", &bx_, "bx/I", bufsize);
245  SiPixelLorentzAngleTreeBarrel_->Branch("orbit", &orbit_, "orbit/I", bufsize);
246  SiPixelLorentzAngleTreeBarrel_->Branch("module", &module_, "module/I", bufsize);
247  SiPixelLorentzAngleTreeBarrel_->Branch("ladder", &ladder_, "ladder/I", bufsize);
248  SiPixelLorentzAngleTreeBarrel_->Branch("layer", &layer_, "layer/I", bufsize);
249  SiPixelLorentzAngleTreeBarrel_->Branch("isflipped", &isflipped_, "isflipped/I", bufsize);
250  SiPixelLorentzAngleTreeBarrel_->Branch("pt", &pt_, "pt/F", bufsize);
251  SiPixelLorentzAngleTreeBarrel_->Branch("eta", &eta_, "eta/F", bufsize);
252  SiPixelLorentzAngleTreeBarrel_->Branch("phi", &phi_, "phi/F", bufsize);
253  SiPixelLorentzAngleTreeBarrel_->Branch("chi2", &chi2_, "chi2/D", bufsize);
254  SiPixelLorentzAngleTreeBarrel_->Branch("ndof", &ndof_, "ndof/D", bufsize);
255  SiPixelLorentzAngleTreeBarrel_->Branch("trackhit", &trackhit_, "x/F:y/F:alpha/D:beta/D:gamma_/D", bufsize);
256  SiPixelLorentzAngleTreeBarrel_->Branch("npix", &pixinfo_.npix, "npix/I", bufsize);
257  SiPixelLorentzAngleTreeBarrel_->Branch("rowpix", pixinfo_.row, "row[npix]/F", bufsize);
258  SiPixelLorentzAngleTreeBarrel_->Branch("colpix", pixinfo_.col, "col[npix]/F", bufsize);
259  SiPixelLorentzAngleTreeBarrel_->Branch("adc", pixinfo_.adc, "adc[npix]/F", bufsize);
260  SiPixelLorentzAngleTreeBarrel_->Branch("xpix", pixinfo_.x, "x[npix]/F", bufsize);
261  SiPixelLorentzAngleTreeBarrel_->Branch("ypix", pixinfo_.y, "y[npix]/F", bufsize);
262 
264  "clust",
265  &clust_,
266  "x/F:y/F:charge/F:size_x/I:size_y/I:maxPixelCol/I:maxPixelRow:minPixelCol/I:minPixelRow/I",
267  bufsize);
268  SiPixelLorentzAngleTreeBarrel_->Branch("rechit", &rechit_, "x/F:y/F", bufsize);
269  SiPixelLorentzAngleTreeBarrel_->Branch("rechit_corr", &rechitCorr_, "x/F:y/F", bufsize);
270  SiPixelLorentzAngleTreeBarrel_->Branch("trackhitcorr_x", &trackhitCorrX_, "trackhitcorr_x/F", bufsize);
271  SiPixelLorentzAngleTreeBarrel_->Branch("trackhitcorr_y", &trackhitCorrY_, "trackhitcorr_y/F", bufsize);
272  SiPixelLorentzAngleTreeBarrel_->Branch("qScale", &qScale_, "qScale/F", bufsize);
273  SiPixelLorentzAngleTreeBarrel_->Branch("rQmQt", &rQmQt_, "rQmQt/F", bufsize);
274  // Forward pixel
275 
277  std::make_unique<TTree>("SiPixelLorentzAngleTreeForward_", "SiPixel LorentzAngle tree forward", bufsize);
278  SiPixelLorentzAngleTreeForward_->Branch("run", &run_, "run/I", bufsize);
279  SiPixelLorentzAngleTreeForward_->Branch("event", &event_, "event/l", bufsize);
280  SiPixelLorentzAngleTreeForward_->Branch("lumiblock", &lumiblock_, "lumiblock/I", bufsize);
281  SiPixelLorentzAngleTreeForward_->Branch("bx", &bx_, "bx/I", bufsize);
282  SiPixelLorentzAngleTreeForward_->Branch("orbit", &orbit_, "orbit/I", bufsize);
283  SiPixelLorentzAngleTreeForward_->Branch("side", &sideF_, "side/I", bufsize);
284  SiPixelLorentzAngleTreeForward_->Branch("disk", &diskF_, "disk/I", bufsize);
285  SiPixelLorentzAngleTreeForward_->Branch("blade", &bladeF_, "blade/I", bufsize);
286  SiPixelLorentzAngleTreeForward_->Branch("panel", &panelF_, "panel/I", bufsize);
287  SiPixelLorentzAngleTreeForward_->Branch("module", &moduleF_, "module/I", bufsize);
288  SiPixelLorentzAngleTreeForward_->Branch("pt", &pt_, "pt/F", bufsize);
289  SiPixelLorentzAngleTreeForward_->Branch("eta", &eta_, "eta/F", bufsize);
290  SiPixelLorentzAngleTreeForward_->Branch("phi", &phi_, "phi/F", bufsize);
291  SiPixelLorentzAngleTreeForward_->Branch("chi2", &chi2_, "chi2/D", bufsize);
292  SiPixelLorentzAngleTreeForward_->Branch("ndof", &ndof_, "ndof/D", bufsize);
293  SiPixelLorentzAngleTreeForward_->Branch("trackhit", &trackhitF_, "x/F:y/F:alpha/D:beta/D:gamma_/D", bufsize);
294  SiPixelLorentzAngleTreeForward_->Branch("npix", &pixinfoF_.npix, "npix/I", bufsize);
295  SiPixelLorentzAngleTreeForward_->Branch("rowpix", pixinfoF_.row, "row[npix]/F", bufsize);
296  SiPixelLorentzAngleTreeForward_->Branch("colpix", pixinfoF_.col, "col[npix]/F", bufsize);
297  SiPixelLorentzAngleTreeForward_->Branch("adc", pixinfoF_.adc, "adc[npix]/F", bufsize);
298  SiPixelLorentzAngleTreeForward_->Branch("xpix", pixinfoF_.x, "x[npix]/F", bufsize);
299  SiPixelLorentzAngleTreeForward_->Branch("ypix", pixinfoF_.y, "y[npix]/F", bufsize);
300 
302  "clust",
303  &clustF_,
304  "x/F:y/F:charge/F:size_x/I:size_y/I:maxPixelCol/I:maxPixelRow:minPixelCol/I:minPixelRow/I",
305  bufsize);
306  SiPixelLorentzAngleTreeForward_->Branch("rechit", &rechitF_, "x/F:y/F", bufsize);
307  SiPixelLorentzAngleTreeForward_->Branch("rechit_corr", &rechitCorrF_, "x/F:y/F", bufsize);
308  SiPixelLorentzAngleTreeForward_->Branch("trackhitcorr_x", &trackhitCorrXF_, "trackhitcorr_x/F", bufsize);
309  SiPixelLorentzAngleTreeForward_->Branch("trackhitcorr_y", &trackhitCorrYF_, "trackhitcorr_y/F", bufsize);
310  SiPixelLorentzAngleTreeForward_->Branch("qScale", &qScaleF_, "qScale/F", bufsize);
311  SiPixelLorentzAngleTreeForward_->Branch("rQmQt", &rQmQtF_, "rQmQt/F", bufsize);
312  }
313 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
float x[maxpix]
edm::ESGetToken< std::vector< SiPixelTemplateStore >, SiPixelTemplateDBObjectESProducerRcd > siPixelTemplateStoreEsToken_
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoPerEventEsToken_
float y[maxpix]
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomPerEventEsToken_
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoEsToken_
std::unique_ptr< TTree > SiPixelLorentzAngleTreeBarrel_
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_
std::vector< std::string > newmodulelist_
LorentzAngleAnalysisTypeEnum convertStringToLorentzAngleAnalysisTypeEnum(std::string type)
float adc[maxpix]
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomEsToken_
float col[maxpix]
edm::ESGetToken< SiPixelTemplateDBObject, SiPixelTemplateDBObjectESProducerRcd > siPixelTemplateEsToken_
edm::EDGetTokenT< TrajTrackAssociationCollection > t_trajTrack
LorentzAngleAnalysisTypeEnum analysisType_
float row[maxpix]
std::unique_ptr< TTree > SiPixelLorentzAngleTreeForward_

◆ ~SiPixelLorentzAnglePCLWorker()

SiPixelLorentzAnglePCLWorker::~SiPixelLorentzAnglePCLWorker ( )
overridedefault

Member Function Documentation

◆ analyze()

void SiPixelLorentzAnglePCLWorker::analyze ( edm::Event const &  iEvent,
edm::EventSetup const &  iSetup 
)
overrideprivatevirtual

Reimplemented from DQMEDAnalyzer.

Definition at line 321 of file SiPixelLorentzAnglePCLWorker.cc.

References funct::abs(), Pixinfo::adc, Hit::alpha, analysisType_, edm::AssociationMap< Tag >::begin(), Hit::beta, SiPixelLorentzAngleCalibrationHistograms::betaStartIdx_, mkfit::Config::Bfield, bladeF_, Surface::bounds(), SiPixelLorentzAngleCalibrationHistograms::BPixnewDetIds_, bx_, Clust::charge, chi2_, Trajectory::chiSquared(), clust_, clustChargeMaxPerLength_, SiPixelRecHit::cluster(), clustF_, clustSizeXMax_, clustSizeYMin_, siPixelLACalibration::cmToum, Pixinfo::col, edm::AssociationMap< edm::OneToOne< std::vector< Trajectory >, reco::TrackCollection, unsigned short > >::const_iterator, hcalRecHitTable_cff::depth, hcalRecHitTable_cff::detId, diskF_, HLT_2024v10_cff::distance, shallow::drift(), PVValHelper::dx, PVValHelper::dy, eGrazingAngle, eMinimumClusterSize, edm::AssociationMap< Tag >::empty(), edm::AssociationMap< Tag >::end(), eta_, event_, dqm::impl::MonitorElement::Fill(), fillPix(), spr::find(), Hit::gamma, geomPerEventEsToken_, edm::EventSetup::getData(), SiPixelTemplateDBObject::getTemplateID(), SiPixelLorentzAngleCalibrationHistograms::h_bySectOccupancy_, SiPixelLorentzAngleCalibrationHistograms::h_drift_depth_adc2_, SiPixelLorentzAngleCalibrationHistograms::h_drift_depth_adc_, SiPixelLorentzAngleCalibrationHistograms::h_drift_depth_noadc_, SiPixelLorentzAngleCalibrationHistograms::h_fpixAngleSize_, SiPixelLorentzAngleCalibrationHistograms::h_fpixMagField_, SiPixelLorentzAngleCalibrationHistograms::h_trackChi2_, SiPixelLorentzAngleCalibrationHistograms::h_trackEta_, SiPixelLorentzAngleCalibrationHistograms::h_trackPhi_, SiPixelLorentzAngleCalibrationHistograms::h_trackPt_, SiPixelLorentzAngleCalibrationHistograms::h_tracks_, heavyIonCSV_trainingSettings::idx, iEvent, iHists, SiPixelTemplate::interpolate(), MagneticField::inTesla(), isflipped_, TrajectoryStateOnSurface::isValid(), ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, dqmiolumiharvest::j, ladder_, layer_, TrajectoryStateOnSurface::localDirection(), BaseTrackerRecHit::localPosition(), TrajectoryStateOnSurface::localPosition(), lumiblock_, M_PI, magneticFieldToken_, Clust::maxPixelCol, Clust::maxPixelRow, Trajectory::measurements(), Clust::minPixelCol, Clust::minPixelRow, module_, moduleF_, Trajectory::ndof(), ndof_, SiPixelLorentzAngleCalibrationHistograms::nlay, SiPixelLorentzAngleCalibrationHistograms::nModules_, MagneticField::nominalValue(), normChi2Max_, notInPCL_, SiPixelLorentzAngleCalibrationHistograms::nPanels_, Pixinfo::npix, SiPixelLorentzAngleCalibrationHistograms::nSides_, orbit_, panelF_, PV3DBase< T, PVType, FrameType >::perp(), TrackerGeometry::Ph1PXB, phi_, PixelTopology::pitch(), PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, pixinfo_, pixinfoF_, GloballyPositioned< T >::position(), funct::pow(), pt_, ptmin_, TrackerTopology::pxbLadder(), TrackerTopology::pxbLayer(), TrackerTopology::pxbModule(), TrackerTopology::pxfBlade(), TrackerTopology::pxfDisk(), TrackerTopology::pxfModule(), TrackerTopology::pxfPanel(), TrackerTopology::pxfSide(), SiPixelTemplate::qscale(), qScale_, qScaleF_, SiPixelTemplate::r_qMeas_qTrue(), DetId::rawId(), rpcPointValidation_cfi::recHit, rechit_, rechitCorr_, rechitCorrF_, rechitF_, residualMax_, Pixinfo::row, rQmQt_, rQmQtF_, run_, sideF_, SiPixelLorentzAngleTreeBarrel_, SiPixelLorentzAngleTreeForward_, Clust::size_x, Clust::size_y, PixelGeomDetUnit::specificTopology(), mathSSE::sqrt(), GeomDet::surface(), surface_deformation(), t_trajTrack, funct::tan(), templateDBobject_, thePixelTemp_, Bounds::thickness(), Surface::toGlobal(), GloballyPositioned< T >::toLocal(), topoPerEventEsToken_, HLT_2024v10_cff::track, DetId::Tracker, PbPb_ZMuSkimMuonDPG_cff::tracker, trackhit_, trackhitCorrX_, trackhitCorrXF_, trackhitCorrY_, trackhitCorrYF_, trackhitF_, x, Pixinfo::x, Hit::x, PV3DBase< T, PVType, FrameType >::x(), Clust::x, Rechit::x, y, Pixinfo::y, Hit::y, PV3DBase< T, PVType, FrameType >::y(), Clust::y, Rechit::y, and PV3DBase< T, PVType, FrameType >::z().

321  {
322  // Retrieve tracker topology from geometry
323  const TrackerTopology* const tTopo = &iSetup.getData(topoPerEventEsToken_);
324 
325  // Retrieve track geometry
326  const TrackerGeometry* tracker = &iSetup.getData(geomPerEventEsToken_);
327 
328  // Retrieve magnetic field
329  const MagneticField* magField = &iSetup.getData(magneticFieldToken_);
330 
331  // get the association map between tracks and trajectories
332  edm::Handle<TrajTrackAssociationCollection> trajTrackCollectionHandle;
333  iEvent.getByToken(t_trajTrack, trajTrackCollectionHandle);
334 
335  module_ = -1;
336  layer_ = -1;
337  ladder_ = -1;
338  isflipped_ = -1;
339  pt_ = -999;
340  eta_ = 999;
341  phi_ = 999;
342  pixinfo_.npix = 0;
343 
344  run_ = iEvent.id().run();
345  event_ = iEvent.id().event();
346  lumiblock_ = iEvent.luminosityBlock();
347  bx_ = iEvent.bunchCrossing();
348  orbit_ = iEvent.orbitNumber();
349 
350  if (!trajTrackCollectionHandle->empty()) {
351  for (TrajTrackAssociationCollection::const_iterator it = trajTrackCollectionHandle->begin();
352  it != trajTrackCollectionHandle->end();
353  ++it) {
354  const reco::Track& track = *it->val;
355  const Trajectory& traj = *it->key;
356 
357  // get the trajectory measurements
358  std::vector<TrajectoryMeasurement> tmColl = traj.measurements();
359  pt_ = track.pt();
360  eta_ = track.eta();
361  phi_ = track.phi();
362  chi2_ = traj.chiSquared();
363  ndof_ = traj.ndof();
364 
365  if (pt_ < ptmin_)
366  continue;
367 
372  iHists.h_tracks_->Fill(0);
373  bool pixeltrack = false;
374 
375  // iterate over trajectory measurements
376  for (const auto& itTraj : tmColl) {
377  if (!itTraj.updatedState().isValid())
378  continue;
379  const TransientTrackingRecHit::ConstRecHitPointer& recHit = itTraj.recHit();
380  if (!recHit->isValid() || recHit->geographicalId().det() != DetId::Tracker)
381  continue;
382  unsigned int subDetID = (recHit->geographicalId().subdetId());
383  if (subDetID == PixelSubdetector::PixelBarrel || subDetID == PixelSubdetector::PixelEndcap) {
384  if (!pixeltrack) {
385  iHists.h_tracks_->Fill(1);
386  }
387  pixeltrack = true;
388  }
389 
390  if (subDetID == PixelSubdetector::PixelBarrel) {
391  DetId detIdObj = recHit->geographicalId();
392  const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(tracker->idToDet(detIdObj));
393  if (!theGeomDet)
394  continue;
395 
396  const PixelTopology* topol = &(theGeomDet->specificTopology());
397 
398  float ypitch_ = topol->pitch().second;
399  float width_ = theGeomDet->surface().bounds().thickness();
400 
401  if (!topol)
402  continue;
403 
404  layer_ = tTopo->pxbLayer(detIdObj);
405  ladder_ = tTopo->pxbLadder(detIdObj);
406  module_ = tTopo->pxbModule(detIdObj);
407 
408  float tmp1 = theGeomDet->surface().toGlobal(Local3DPoint(0., 0., 0.)).perp();
409  float tmp2 = theGeomDet->surface().toGlobal(Local3DPoint(0., 0., 1.)).perp();
410 
411  isflipped_ = (tmp2 < tmp1) ? 1 : 0;
412 
413  const SiPixelRecHit* recHitPix = dynamic_cast<const SiPixelRecHit*>((*recHit).hit());
414  if (!recHitPix)
415  continue;
416  rechit_.x = recHitPix->localPosition().x();
417  rechit_.y = recHitPix->localPosition().y();
418  SiPixelRecHit::ClusterRef const& cluster = recHitPix->cluster();
419 
420  pixinfo_ = fillPix(*cluster, topol);
421 
422  // fill entries in clust_
423 
424  clust_.x = (cluster)->x();
425  clust_.y = (cluster)->y();
426  clust_.charge = (cluster->charge()) / 1000.; // clust_.charge: in the unit of 1000e
427  clust_.size_x = cluster->sizeX();
428  clust_.size_y = cluster->sizeY();
429  clust_.maxPixelCol = cluster->maxPixelCol();
430  clust_.maxPixelRow = cluster->maxPixelRow();
431  clust_.minPixelCol = cluster->minPixelCol();
432  clust_.minPixelRow = cluster->minPixelRow();
433 
434  // fill the trackhit info
435  TrajectoryStateOnSurface tsos = itTraj.updatedState();
436  if (!tsos.isValid()) {
437  edm::LogWarning("SiPixelLorentzAnglePCLWorker") << "tsos not valid";
438  continue;
439  }
440  LocalVector trackdirection = tsos.localDirection();
441  LocalPoint trackposition = tsos.localPosition();
442 
443  if (trackdirection.z() == 0)
444  continue;
445  // the local position and direction
446  trackhit_.alpha = atan2(trackdirection.z(), trackdirection.x());
447  trackhit_.beta = atan2(trackdirection.z(), trackdirection.y());
448  trackhit_.gamma = atan2(trackdirection.x(), trackdirection.y());
449  trackhit_.x = trackposition.x();
450  trackhit_.y = trackposition.y();
451 
452  // get qScale_ = templ.qscale() and templ.r_qMeas_qTrue();
453  float cotalpha = trackdirection.x() / trackdirection.z();
454  float cotbeta = trackdirection.y() / trackdirection.z();
455  float cotbeta_min = clustSizeYMin_[layer_ - 1] * ypitch_ / width_;
456  if (std::abs(cotbeta) <= cotbeta_min)
457  continue;
458  double drdz = sqrt(1. + cotalpha * cotalpha + cotbeta * cotbeta);
459  double clusterCharge_cut = clustChargeMaxPerLength_ * drdz;
460 
461  auto detId = detIdObj.rawId();
462  int DetId_index = -1;
463 
464  const auto& newModIt = (std::find(iHists.BPixnewDetIds_.begin(), iHists.BPixnewDetIds_.end(), detId));
465  bool isNewMod = (newModIt != iHists.BPixnewDetIds_.end());
466  if (isNewMod) {
467  DetId_index = std::distance(iHists.BPixnewDetIds_.begin(), newModIt);
468  }
469 
470  if (notInPCL_) {
471  // fill the template from the store (from dqmBeginRun)
472  SiPixelTemplate theTemplate(*thePixelTemp_);
473 
474  float locBx = (cotbeta < 0.) ? -1 : 1.;
475  float locBz = (cotalpha < 0.) ? -locBx : locBx;
476 
477  int TemplID = templateDBobject_->getTemplateID(detId);
478  theTemplate.interpolate(TemplID, cotalpha, cotbeta, locBz, locBx);
479  qScale_ = theTemplate.qscale();
480  rQmQt_ = theTemplate.r_qMeas_qTrue();
481  }
482 
483  // Surface deformation
484  const auto& lp_pair = surface_deformation(topol, tsos, recHitPix);
485 
486  LocalPoint lp_track = lp_pair.first;
487  LocalPoint lp_rechit = lp_pair.second;
488 
489  rechitCorr_.x = lp_rechit.x();
490  rechitCorr_.y = lp_rechit.y();
491  trackhitCorrX_ = lp_track.x();
492  trackhitCorrY_ = lp_track.y();
493 
494  if (notInPCL_) {
496  }
497 
499  continue;
500  // is one pixel in cluster a large pixel ? (hit will be excluded)
501  bool large_pix = false;
502  for (int j = 0; j < pixinfo_.npix; j++) {
503  int colpos = static_cast<int>(pixinfo_.col[j]);
504  if (pixinfo_.row[j] == 0 || pixinfo_.row[j] == 79 || pixinfo_.row[j] == 80 || pixinfo_.row[j] == 159 ||
505  colpos % 52 == 0 || colpos % 52 == 51) {
506  large_pix = true;
507  }
508  }
509 
510  double residualsq = (trackhitCorrX_ - rechitCorr_.x) * (trackhitCorrX_ - rechitCorr_.x) +
512 
513  double xlim1 = trackhitCorrX_ - width_ * cotalpha / 2.;
514  double hypitch_ = ypitch_ / 2.;
515  double ylim1 = trackhitCorrY_ - width_ * cotbeta / 2.;
516  double ylim2 = trackhitCorrY_ + width_ * cotbeta / 2.;
517 
518  int clustSizeY_cut = clustSizeYMin_[layer_ - 1];
519 
520  if (!large_pix && (chi2_ / ndof_) < normChi2Max_ && cluster->sizeY() >= clustSizeY_cut &&
521  residualsq < residualMax_ * residualMax_ && cluster->charge() < clusterCharge_cut &&
522  cluster->sizeX() < clustSizeXMax_) {
523  // iterate over pixels in hit
524  for (int j = 0; j < pixinfo_.npix; j++) {
525  // use trackhits and include bowing correction
526  float ypixlow = pixinfo_.y[j] - hypitch_;
527  float ypixhigh = pixinfo_.y[j] + hypitch_;
528  if (cotbeta > 0.) {
529  if (ylim1 > ypixlow)
530  ypixlow = ylim1;
531  if (ylim2 < ypixhigh)
532  ypixhigh = ylim2;
533  } else {
534  if (ylim2 > ypixlow)
535  ypixlow = ylim2;
536  if (ylim1 < ypixhigh)
537  ypixhigh = ylim1;
538  }
539  float ypixavg = 0.5f * (ypixlow + ypixhigh);
540 
541  float dx = (pixinfo_.x[j] - xlim1) * siPixelLACalibration::cmToum; // dx: in the unit of micrometer
542  float dy = (ypixavg - ylim1) * siPixelLACalibration::cmToum; // dy: in the unit of micrometer
543  float depth = dy * tan(trackhit_.beta);
544  float drift = dx - dy * tan(trackhit_.gamma);
545 
546  if (isNewMod == false) {
547  int i_index = module_ + (layer_ - 1) * iHists.nModules_[layer_ - 1];
548  iHists.h_drift_depth_adc_[i_index]->Fill(drift, depth, pixinfo_.adc[j]);
550  iHists.h_drift_depth_noadc_[i_index]->Fill(drift, depth, 1.);
551  iHists.h_bySectOccupancy_->Fill(i_index - 1); // histogram starts at 0
552 
553  if (tracker->getDetectorType(subDetID) == TrackerGeometry::ModuleType::Ph1PXB) {
554  if ((module_ == 3 || module_ == 5) && (layer_ == 3 || layer_ == 4)) {
555  int i_index_merge = i_index + 1;
556  iHists.h_drift_depth_adc_[i_index_merge]->Fill(drift, depth, pixinfo_.adc[j]);
557  iHists.h_drift_depth_adc2_[i_index_merge]->Fill(drift, depth, pixinfo_.adc[j] * pixinfo_.adc[j]);
558  iHists.h_drift_depth_noadc_[i_index_merge]->Fill(drift, depth, 1.);
559  iHists.h_bySectOccupancy_->Fill(i_index_merge - 1);
560  }
561  if ((module_ == 4 || module_ == 6) && (layer_ == 3 || layer_ == 4)) {
562  int i_index_merge = i_index - 1;
563  iHists.h_drift_depth_adc_[i_index_merge]->Fill(drift, depth, pixinfo_.adc[j]);
564  iHists.h_drift_depth_adc2_[i_index_merge]->Fill(drift, depth, pixinfo_.adc[j] * pixinfo_.adc[j]);
565  iHists.h_drift_depth_noadc_[i_index_merge]->Fill(drift, depth, 1.);
566  iHists.h_bySectOccupancy_->Fill(i_index_merge - 1);
567  }
568  }
569 
570  } else {
571  int new_index = iHists.nModules_[iHists.nlay - 1] +
572  (iHists.nlay - 1) * iHists.nModules_[iHists.nlay - 1] + 1 + DetId_index;
573 
574  iHists.h_drift_depth_adc_[new_index]->Fill(drift, depth, pixinfo_.adc[j]);
575  iHists.h_drift_depth_adc2_[new_index]->Fill(drift, depth, pixinfo_.adc[j] * pixinfo_.adc[j]);
576  iHists.h_drift_depth_noadc_[new_index]->Fill(drift, depth, 1.);
577  iHists.h_bySectOccupancy_->Fill(new_index - 1); // histogram starts at 0
578  }
579  }
580  }
581  } else if (subDetID == PixelSubdetector::PixelEndcap) {
582  DetId detIdObj = recHit->geographicalId();
583  const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(tracker->idToDet(detIdObj));
584  if (!theGeomDet)
585  continue;
586 
587  const PixelTopology* topol = &(theGeomDet->specificTopology());
588 
589  if (!topol)
590  continue;
591 
592  sideF_ = tTopo->pxfSide(detIdObj);
593  diskF_ = tTopo->pxfDisk(detIdObj);
594  bladeF_ = tTopo->pxfBlade(detIdObj);
595  panelF_ = tTopo->pxfPanel(detIdObj);
596  moduleF_ = tTopo->pxfModule(detIdObj);
597 
598  const SiPixelRecHit* recHitPix = dynamic_cast<const SiPixelRecHit*>((*recHit).hit());
599  if (!recHitPix)
600  continue;
601  rechitF_.x = recHitPix->localPosition().x();
602  rechitF_.y = recHitPix->localPosition().y();
603  SiPixelRecHit::ClusterRef const& cluster = recHitPix->cluster();
604 
605  pixinfoF_ = fillPix(*cluster, topol);
606 
607  // fill entries in clust_
608 
609  clustF_.x = (cluster)->x();
610  clustF_.y = (cluster)->y();
611  clustF_.charge = (cluster->charge()) / 1000.; // clustF_.charge: in the unit of 1000e
612  clustF_.size_x = cluster->sizeX();
613  clustF_.size_y = cluster->sizeY();
614  clustF_.maxPixelCol = cluster->maxPixelCol();
615  clustF_.maxPixelRow = cluster->maxPixelRow();
616  clustF_.minPixelCol = cluster->minPixelCol();
617  clustF_.minPixelRow = cluster->minPixelRow();
618 
619  // fill the trackhit info
620  TrajectoryStateOnSurface tsos = itTraj.updatedState();
621  if (!tsos.isValid()) {
622  edm::LogWarning("SiPixelLorentzAnglePCLWorker") << "tsos not valid";
623  continue;
624  }
625  LocalVector trackdirection = tsos.localDirection();
626  LocalPoint trackposition = tsos.localPosition();
627 
628  if (trackdirection.z() == 0)
629  continue;
630  // the local position and direction
631  trackhitF_.alpha = atan2(trackdirection.z(), trackdirection.x());
632  trackhitF_.beta = atan2(trackdirection.z(), trackdirection.y());
633  trackhitF_.gamma = atan2(trackdirection.x(), trackdirection.y());
634  trackhitF_.x = trackposition.x();
635  trackhitF_.y = trackposition.y();
636 
637  float cotalpha = trackdirection.x() / trackdirection.z();
638  float cotbeta = trackdirection.y() / trackdirection.z();
639 
640  auto detId = detIdObj.rawId();
641 
642  if (notInPCL_) {
643  // fill the template from the store (from dqmBeginRun)
644  SiPixelTemplate theTemplate(*thePixelTemp_);
645 
646  float locBx = (cotbeta < 0.) ? -1 : 1.;
647  float locBz = (cotalpha < 0.) ? -locBx : locBx;
648 
649  int TemplID = templateDBobject_->getTemplateID(detId);
650  theTemplate.interpolate(TemplID, cotalpha, cotbeta, locBz, locBx);
651  qScaleF_ = theTemplate.qscale();
652  rQmQtF_ = theTemplate.r_qMeas_qTrue();
653  }
654 
655  // Surface deformation
656  const auto& lp_pair = surface_deformation(topol, tsos, recHitPix);
657 
658  LocalPoint lp_track = lp_pair.first;
659  LocalPoint lp_rechit = lp_pair.second;
660 
661  rechitCorrF_.x = lp_rechit.x();
662  rechitCorrF_.y = lp_rechit.y();
663  trackhitCorrXF_ = lp_track.x();
664  trackhitCorrYF_ = lp_track.y();
665  if (notInPCL_) {
667  }
668 
670  continue;
671 
672  int theMagField = magField->nominalValue();
673  if (theMagField < 37 || theMagField > 39)
674  continue;
675 
676  double chi2_ndof = chi2_ / ndof_;
677  if (chi2_ndof >= normChi2Max_)
678  continue;
679 
680  //--- large pixel cut
681  bool large_pix = false;
682  for (int j = 0; j < pixinfoF_.npix; j++) {
683  int colpos = static_cast<int>(pixinfoF_.col[j]);
684  if (pixinfoF_.row[j] == 0 || pixinfoF_.row[j] == 79 || pixinfoF_.row[j] == 80 || pixinfoF_.row[j] == 159 ||
685  colpos % 52 == 0 || colpos % 52 == 51) {
686  large_pix = true;
687  }
688  }
689 
690  if (large_pix)
691  continue;
692 
693  //--- residual cut
694  double residual = sqrt(pow(trackhitCorrXF_ - rechitCorrF_.x, 2) + pow(trackhitCorrYF_ - rechitCorrF_.y, 2));
695 
696  if (residual > residualMax_)
697  continue;
698 
699  int ringIdx = bladeF_ <= 22 ? 0 : 1;
700  int panelIdx = panelF_ - 1;
701  int sideIdx = sideF_ - 1;
702  int idx = iHists.nSides_ * iHists.nPanels_ * ringIdx + iHists.nSides_ * panelIdx + sideIdx;
703  int idxBeta = iHists.betaStartIdx_ + idx;
704 
705  double cotanAlpha = std::tan(M_PI / 2. - trackhitF_.alpha);
706  double cotanBeta = std::tan(M_PI / 2. - trackhitF_.beta);
707 
708  LocalVector Bfield = theGeomDet->surface().toLocal(magField->inTesla(theGeomDet->surface().position()));
709  iHists.h_fpixMagField_[0][idx]->Fill(Bfield.x());
710  iHists.h_fpixMagField_[1][idx]->Fill(Bfield.y());
711  iHists.h_fpixMagField_[2][idx]->Fill(Bfield.z());
712 
713  if (clustF_.size_y >= 2) {
714  iHists.h_fpixAngleSize_[idx]->Fill(cotanAlpha, clustF_.size_x);
715  }
716 
717  if (clust_.size_x >= 0) {
718  iHists.h_fpixAngleSize_[idxBeta]->Fill(cotanBeta, clustF_.size_y);
719  }
720  }
721  } //end iteration over trajectory measurements
722  } //end iteration over trajectories
723  }
724 }
ClusterRef cluster() const
Definition: SiPixelRecHit.h:47
unsigned int pxbLayer(const DetId &id) const
T perp() const
Definition: PV3DBase.h:69
float x[maxpix]
unsigned int pxfBlade(const DetId &id) const
T z() const
Definition: PV3DBase.h:61
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoPerEventEsToken_
float y[maxpix]
unsigned int pxfModule(const DetId &id) const
LocalVector drift(const StripGeomDetUnit *, const MagneticField &, const SiStripLorentzAngle &)
Definition: ShallowTools.cc:36
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomPerEventEsToken_
float chiSquared() const
Definition: Trajectory.h:241
unsigned int pxbLadder(const DetId &id) const
bool empty() const
return true if empty
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
LocalPoint toLocal(const GlobalPoint &gp) const
const Pixinfo fillPix(const SiPixelCluster &LocPix, const PixelTopology *topol) const
DataContainer const & measurements() const
Definition: Trajectory.h:178
const_iterator end() const
last iterator over the map (read only)
const std::pair< LocalPoint, LocalPoint > surface_deformation(const PixelTopology *topol, TrajectoryStateOnSurface &tsos, const SiPixelRecHit *recHitPix) const
void Fill(long long x)
virtual float thickness() const =0
double beta
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
int iEvent
Definition: GenABIO.cc:224
int ndof(bool bon=true) const
Definition: Trajectory.cc:97
const SiPixelTemplateDBObject * templateDBobject_
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
T sqrt(T t)
Definition: SSEVec.h:19
LocalVector localDirection() const
std::unique_ptr< TTree > SiPixelLorentzAngleTreeBarrel_
constexpr float Bfield
Definition: Config.h:55
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_
unsigned int pxfDisk(const DetId &id) const
const std::vector< SiPixelTemplateStore > * thePixelTemp_
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define M_PI
float adc[maxpix]
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
unsigned int pxfPanel(const DetId &id) const
Definition: DetId.h:17
unsigned int pxfSide(const DetId &id) const
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
const PositionType & position() const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
double alpha
SiPixelLorentzAngleCalibrationHistograms iHists
short getTemplateID(const uint32_t &detid) const
float col[maxpix]
const_iterator begin() const
first iterator over the map (read only)
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
LocalPoint localPosition() const override
edm::EDGetTokenT< TrajTrackAssociationCollection > t_trajTrack
int nominalValue() const
The nominal field value for this map in kGauss.
Definition: MagneticField.h:49
double gamma
virtual std::pair< float, float > pitch() const =0
unsigned int pxbModule(const DetId &id) const
Log< level::Warning, false > LogWarning
LorentzAngleAnalysisTypeEnum analysisType_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
float row[maxpix]
Our base class.
Definition: SiPixelRecHit.h:23
std::unique_ptr< TTree > SiPixelLorentzAngleTreeForward_
const Bounds & bounds() const
Definition: Surface.h:87
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9

◆ bookHistograms()

void SiPixelLorentzAnglePCLWorker::bookHistograms ( DQMStore::IBooker iBooker,
edm::Run const &  run,
edm::EventSetup const &  iSetup 
)
overrideprivatevirtual

Implements DQMEDAnalyzer.

Definition at line 771 of file SiPixelLorentzAnglePCLWorker.cc.

References analysisType_, SiPixelLorentzAngleCalibrationHistograms::betaStartIdx_, dqm::implementation::IBooker::book1D(), dqm::implementation::IBooker::book2D(), SiPixelLorentzAngleCalibrationHistograms::BPixnewDetIds_, SiPixelLorentzAngleCalibrationHistograms::BPixnewmodulename_, ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), eGrazingAngle, folder_, SiPixelLorentzAngleCalibrationHistograms::h_bySectOccupancy_, SiPixelLorentzAngleCalibrationHistograms::h_drift_depth_, SiPixelLorentzAngleCalibrationHistograms::h_drift_depth_adc2_, SiPixelLorentzAngleCalibrationHistograms::h_drift_depth_adc_, SiPixelLorentzAngleCalibrationHistograms::h_drift_depth_noadc_, SiPixelLorentzAngleCalibrationHistograms::h_fpixAngleSize_, SiPixelLorentzAngleCalibrationHistograms::h_fpixMagField_, SiPixelLorentzAngleCalibrationHistograms::h_fpixMean_, SiPixelLorentzAngleCalibrationHistograms::h_mean_, SiPixelLorentzAngleCalibrationHistograms::h_trackChi2_, SiPixelLorentzAngleCalibrationHistograms::h_trackEta_, SiPixelLorentzAngleCalibrationHistograms::h_trackPhi_, SiPixelLorentzAngleCalibrationHistograms::h_trackPt_, SiPixelLorentzAngleCalibrationHistograms::h_tracks_, hist_depth_, hist_drift_, mps_fire::i, heavyIonCSV_trainingSettings::idx, iHists, createfilelist::int, LogDebug, visualization-live-secondInstance_cfg::m, M_PI, Skims_PA_cff::name, SiPixelLorentzAngleCalibrationHistograms::nlay, SiPixelLorentzAngleCalibrationHistograms::nModules_, SiPixelLorentzAngleCalibrationHistograms::nPanels_, SiPixelLorentzAngleCalibrationHistograms::nRings_, SiPixelLorentzAngleCalibrationHistograms::nSides_, AlCaHLTBitMon_ParallelJobs::p, alignCSCRings::r, alignCSCRings::s, dqm::impl::MonitorElement::setBinLabel(), dqm::implementation::NavigatorBase::setCurrentFolder(), AlCaHLTBitMon_QueryRunRegistry::string, and runGCPTkAlMap::title.

773  {
776  if (analysisType_ == eGrazingAngle) {
777  // book the by partition monitoring
778  const auto maxSect = iHists.nlay * iHists.nModules_[iHists.nlay - 1] + (int)iHists.BPixnewDetIds_.size();
779 
780  iBooker.setCurrentFolder(fmt::sprintf("%s/SectorMonitoring", folder_.data()));
781  iHists.h_bySectOccupancy_ = iBooker.book1D(
782  "h_bySectorOccupancy", "hit occupancy by sector;pixel sector;hits on track", maxSect, -0.5, maxSect - 0.5);
783 
784  iBooker.setCurrentFolder(folder_);
785  static constexpr double min_depth_ = -100.;
786  static constexpr double max_depth_ = 400.;
787  static constexpr double min_drift_ = -500.;
788  static constexpr double max_drift_ = 500.;
789 
790  // book the mean values projections and set the bin names of the by sector monitoring
791  for (int i_layer = 1; i_layer <= iHists.nlay; i_layer++) {
792  for (int i_module = 1; i_module <= iHists.nModules_[i_layer - 1]; i_module++) {
793  unsigned int i_index = i_module + (i_layer - 1) * iHists.nModules_[i_layer - 1];
794  std::string binName = fmt::sprintf("BPix Lay%i Mod%i", i_layer, i_module);
795  LogDebug("SiPixelLorentzAnglePCLWorker") << " i_index: " << i_index << " bin name: " << binName
796  << " (i_layer: " << i_layer << " i_module:" << i_module << ")";
797 
798  iHists.h_bySectOccupancy_->setBinLabel(i_index, binName);
799 
800  name = fmt::sprintf("h_mean_layer%i_module%i", i_layer, i_module);
801  title = fmt::sprintf(
802  "average drift vs depth layer%i module%i; production depth [#mum]; #LTdrift#GT [#mum]", i_layer, i_module);
803  iHists.h_mean_[i_index] = iBooker.book1D(name, title, hist_depth_, min_depth_, max_depth_);
804  }
805  }
806  for (int i = 0; i < (int)iHists.BPixnewDetIds_.size(); i++) {
807  name = fmt::sprintf("h_BPixnew_mean_%s", iHists.BPixnewmodulename_[i].c_str());
808  title = fmt::sprintf("average drift vs depth %s; production depth [#mum]; #LTdrift#GT [#mum]",
809  iHists.BPixnewmodulename_[i].c_str());
810  int new_index = iHists.nModules_[iHists.nlay - 1] + (iHists.nlay - 1) * iHists.nModules_[iHists.nlay - 1] + 1 + i;
811  iHists.h_mean_[new_index] = iBooker.book1D(name, title, hist_depth_, min_depth_, max_depth_);
812 
813  LogDebug("SiPixelLorentzAnglePCLWorker")
814  << "i_index" << new_index << " bin name: " << iHists.BPixnewmodulename_[i];
815 
817  }
818 
819  //book the 2D histograms
820  for (int i_layer = 1; i_layer <= iHists.nlay; i_layer++) {
821  iBooker.setCurrentFolder(fmt::sprintf("%s/BPix/BPixLayer%i", folder_.data(), i_layer));
822  for (int i_module = 1; i_module <= iHists.nModules_[i_layer - 1]; i_module++) {
823  unsigned int i_index = i_module + (i_layer - 1) * iHists.nModules_[i_layer - 1];
824 
825  name = fmt::sprintf("h_drift_depth_adc_layer%i_module%i", i_layer, i_module);
826  title = fmt::sprintf(
827  "depth vs drift (ADC) layer%i module%i; drift [#mum]; production depth [#mum]", i_layer, i_module);
828  iHists.h_drift_depth_adc_[i_index] =
829  iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
830 
831  name = fmt::sprintf("h_drift_depth_adc2_layer%i_module%i", i_layer, i_module);
832  title = fmt::sprintf(
833  "depth vs drift (ADC^{2}) layer%i module%i; drift [#mum]; production depth [#mum]", i_layer, i_module);
834  iHists.h_drift_depth_adc2_[i_index] =
835  iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
836 
837  name = fmt::sprintf("h_drift_depth_noadc_layer%i_module%i", i_layer, i_module);
838  title = fmt::sprintf(
839  "depth vs drift (no ADC) layer%i module%i; drift [#mum]; production depth [#mum]", i_layer, i_module);
840  iHists.h_drift_depth_noadc_[i_index] =
841  iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
842 
843  name = fmt::sprintf("h_drift_depth_layer%i_module%i", i_layer, i_module);
844  title =
845  fmt::sprintf("depth vs drift layer%i module%i; drift [#mum]; production depth [#mum]", i_layer, i_module);
846  iHists.h_drift_depth_[i_index] =
847  iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
848  }
849  }
850 
851  // book the "new" modules
852  iBooker.setCurrentFolder(fmt::sprintf("%s/BPix/NewModules", folder_.data()));
853  for (int i = 0; i < (int)iHists.BPixnewDetIds_.size(); i++) {
854  int new_index = iHists.nModules_[iHists.nlay - 1] + (iHists.nlay - 1) * iHists.nModules_[iHists.nlay - 1] + 1 + i;
855 
856  name = fmt::sprintf("h_BPixnew_drift_depth_adc_%s", iHists.BPixnewmodulename_[i].c_str());
857  title = fmt::sprintf("depth vs drift (ADC) %s; drift [#mum]; production depth [#mum]",
858  iHists.BPixnewmodulename_[i].c_str());
859  iHists.h_drift_depth_adc_[new_index] =
860  iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
861 
862  name = fmt::sprintf("h_BPixnew_drift_depth_adc2_%s", iHists.BPixnewmodulename_[i].c_str());
863  title = fmt::sprintf("depth vs drift (ADC^{2}) %s; drift [#mum]; production depth [#mum]",
864  iHists.BPixnewmodulename_[i].c_str());
865  iHists.h_drift_depth_adc2_[new_index] =
866  iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
867 
868  name = fmt::sprintf("h_BPixnew_drift_depth_noadc_%s", iHists.BPixnewmodulename_[i].c_str());
869  title = fmt::sprintf("depth vs drift (no ADC)%s; drift [#mum]; production depth [#mum]",
870  iHists.BPixnewmodulename_[i].c_str());
871  iHists.h_drift_depth_noadc_[new_index] =
872  iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
873 
874  name = fmt::sprintf("h_BPixnew_drift_depth_%s", iHists.BPixnewmodulename_[i].c_str());
875  title = fmt::sprintf("depth vs drift %s; drift [#mum]; production depth [#mum]",
876  iHists.BPixnewmodulename_[i].c_str());
877  iHists.h_drift_depth_[new_index] =
878  iBooker.book2D(name, title, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
879  }
880  } // end if GrazinAngleAnalysis
881  else {
882  iBooker.setCurrentFolder(folder_);
883  std::string baseName;
884  std::string baseTitle;
885 
886  for (int r = 0; r < iHists.nRings_; ++r) {
887  for (int p = 0; p < iHists.nPanels_; ++p) {
888  for (int s = 0; s < iHists.nSides_; ++s) {
889  baseName = fmt::sprintf("R%d_P%d_z%d", r + 1, p + 1, s + 1);
890  if (s == 0)
891  baseTitle = fmt::sprintf("Ring%d_Panel%d_z-", r + 1, p + 1);
892  else
893  baseTitle = fmt::sprintf("Ring%d_Panel%d_z+", r + 1, p + 1);
894 
895  int idx = iHists.nSides_ * iHists.nPanels_ * r + iHists.nSides_ * p + s;
896  int idxBeta = iHists.betaStartIdx_ + idx;
897 
898  name = fmt::sprintf("%s_alphaMean", baseName);
899  title = fmt::sprintf("%s_alphaMean;cot(#alpha); Average cluster size x (pixel)", baseTitle);
900  iHists.h_fpixMean_[idx] = iBooker.book1D(name, title, 60, -3., 3.);
901  name = fmt::sprintf("%s_betaMean", baseName);
902  title = fmt::sprintf("%s_betaMean;cot(#beta); Average cluster size y (pixel)", baseTitle);
903  iHists.h_fpixMean_[idxBeta] = iBooker.book1D(name, title, 60, -3., 3.);
904 
905  } // loop over sides
906  } // loop over panels
907  } // loop over rings
908  iBooker.setCurrentFolder(fmt::sprintf("%s/FPix", folder_.data()));
909  for (int r = 0; r < iHists.nRings_; ++r) {
910  for (int p = 0; p < iHists.nPanels_; ++p) {
911  for (int s = 0; s < iHists.nSides_; ++s) {
912  baseName = fmt::sprintf("R%d_P%d_z%d", r + 1, p + 1, s + 1);
913  if (s == 0)
914  baseTitle = fmt::sprintf("Ring%d_Panel%d_z-", r + 1, p + 1);
915  else
916  baseTitle = fmt::sprintf("Ring%d_Panel%d_z+", r + 1, p + 1);
917 
918  int idx = iHists.nSides_ * iHists.nPanels_ * r + iHists.nSides_ * p + s;
919  int idxBeta = iHists.betaStartIdx_ + idx;
920 
921  name = fmt::sprintf("%s_alpha", baseName);
922  title = fmt::sprintf("%s_alpha;cot(#alpha); Cluster size x (pixel)", baseTitle);
923  iHists.h_fpixAngleSize_[idx] = iBooker.book2D(name, title, 60, -3., 3., 10, 0.5, 10.5);
924  name = fmt::sprintf("%s_beta", baseName);
925  title = fmt::sprintf("%s_beta;cot(#beta); Cluster size y (pixel) ", baseTitle);
926  iHists.h_fpixAngleSize_[idxBeta] = iBooker.book2D(name, title, 60, -3., 3., 10, 0.5, 10.5);
927  for (int m = 0; m < 3; ++m) {
928  name = fmt::sprintf("%s_B%d", baseName, m);
929  char bComp = m == 0 ? 'x' : (m == 1 ? 'y' : 'z');
930  title = fmt::sprintf("%s_magField%d;B_{%c} [T];Entries", baseTitle, m, bComp);
931  iHists.h_fpixMagField_[m][idx] = iBooker.book1D(name, title, 10000, -5., 5.);
932  } // mag. field comps
933  } // loop over sides
934  } // loop over panels
935  } // loop over rings
936  } // if MinimalClusterSize
937 
938  // book the track monitoring plots
939  iBooker.setCurrentFolder(fmt::sprintf("%s/TrackMonitoring", folder_.data()));
940  iHists.h_tracks_ = iBooker.book1D("h_tracks", ";tracker volume;tracks", 2, -0.5, 1.5);
941  iHists.h_tracks_->setBinLabel(1, "all tracks", 1);
942  iHists.h_tracks_->setBinLabel(2, "has pixel hits", 1);
943  iHists.h_trackEta_ = iBooker.book1D("h_trackEta", ";track #eta; #tracks", 30, -3., 3.);
944  iHists.h_trackPhi_ = iBooker.book1D("h_trackPhi", ";track #phi; #tracks", 48, -M_PI, M_PI);
945  iHists.h_trackPt_ = iBooker.book1D("h_trackPt", ";track p_{T} [GeV]; #tracks", 100, 0., 100.);
946  iHists.h_trackChi2_ = iBooker.book1D("h_trackChi2ndof", ";track #chi^{2}/ndof; #tracks", 100, 0., 10.);
947 }
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
virtual void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
#define M_PI
SiPixelLorentzAngleCalibrationHistograms iHists
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:212
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
LorentzAngleAnalysisTypeEnum analysisType_
#define LogDebug(id)

◆ convertStringToLorentzAngleAnalysisTypeEnum()

LorentzAngleAnalysisTypeEnum SiPixelLorentzAnglePCLWorker::convertStringToLorentzAngleAnalysisTypeEnum ( std::string  type)
private

◆ dqmBeginRun()

void SiPixelLorentzAnglePCLWorker::dqmBeginRun ( edm::Run const &  run,
edm::EventSetup const &  iSetup 
)
overrideprivatevirtual

Reimplemented from DQMEDAnalyzer.

Definition at line 726 of file SiPixelLorentzAnglePCLWorker.cc.

References PixelEndcapName::bladeName(), SiPixelLorentzAngleCalibrationHistograms::BPixnewDetIds_, SiPixelLorentzAngleCalibrationHistograms::BPixnewLayer_, SiPixelLorentzAngleCalibrationHistograms::BPixnewModule_, SiPixelLorentzAngleCalibrationHistograms::BPixnewmodulename_, edm::ESWatcher< T >::check(), hcalRecHitTable_cff::detId, PixelEndcapName::diskName(), SiPixelLorentzAngleCalibrationHistograms::FPixnewBlade_, SiPixelLorentzAngleCalibrationHistograms::FPixnewDetIds_, SiPixelLorentzAngleCalibrationHistograms::FPixnewDisk_, SiPixelLorentzAngleCalibrationHistograms::FPixnewmodulename_, relativeConstraints::geom, geomEsToken_, edm::EventSetup::getData(), PixelBarrelName::getDetId(), PixelEndcapName::getDetId(), mps_fire::i, iHists, PixelBarrelName::layerName(), genParticles_cff::map, PixelBarrelName::moduleName(), newmodulelist_, SiPixelLorentzAngleCalibrationHistograms::nlay, SiPixelLorentzAngleCalibrationHistograms::nModules_, notInPCL_, PixelSubdetector::PixelBarrel, siPixelTemplateEsToken_, siPixelTemplateStoreEsToken_, templateDBobject_, thePixelTemp_, topoEsToken_, and watchSiPixelTemplateRcd_.

726  {
727  // geometry
728  const TrackerGeometry* geom = &iSetup.getData(geomEsToken_);
729  const TrackerTopology* tTopo = &iSetup.getData(topoEsToken_);
730 
731  if (notInPCL_) {
732  // Initialize 1D templates
733  if (watchSiPixelTemplateRcd_.check(iSetup)) {
734  templateDBobject_ = &iSetup.getData(siPixelTemplateEsToken_);
735  thePixelTemp_ = &iSetup.getData(siPixelTemplateStoreEsToken_);
736  }
737  }
738 
740  iHists.nlay = geom->numberOfLayers(PixelSubdetector::PixelBarrel);
741  iHists.nModules_.resize(iHists.nlay);
742  for (int i = 0; i < iHists.nlay; i++) {
743  iHists.nModules_[i] = map.getPXBModules(i + 1);
744  }
745 
746  // list of modules already filled, then return (we already entered here)
747  if (!iHists.BPixnewDetIds_.empty() || !iHists.FPixnewDetIds_.empty())
748  return;
749 
750  if (!newmodulelist_.empty()) {
751  for (auto const& modulename : newmodulelist_) {
752  if (modulename.find("BPix_") != std::string::npos) {
753  PixelBarrelName bn(modulename, true);
754  const auto& detId = bn.getDetId(tTopo);
755  iHists.BPixnewmodulename_.push_back(modulename);
756  iHists.BPixnewDetIds_.push_back(detId.rawId());
757  iHists.BPixnewModule_.push_back(bn.moduleName());
758  iHists.BPixnewLayer_.push_back(bn.layerName());
759  } else if (modulename.find("FPix_") != std::string::npos) {
760  PixelEndcapName en(modulename, true);
761  const auto& detId = en.getDetId(tTopo);
762  iHists.FPixnewmodulename_.push_back(modulename);
763  iHists.FPixnewDetIds_.push_back(detId.rawId());
764  iHists.FPixnewDisk_.push_back(en.diskName());
765  iHists.FPixnewBlade_.push_back(en.bladeName());
766  }
767  }
768  }
769 }
edm::ESGetToken< std::vector< SiPixelTemplateStore >, SiPixelTemplateDBObjectESProducerRcd > siPixelTemplateStoreEsToken_
edm::ESWatcher< SiPixelTemplateDBObjectESProducerRcd > watchSiPixelTemplateRcd_
const SiPixelTemplateDBObject * templateDBobject_
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoEsToken_
const std::vector< SiPixelTemplateStore > * thePixelTemp_
std::vector< std::string > newmodulelist_
SiPixelLorentzAngleCalibrationHistograms iHists
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomEsToken_
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
edm::ESGetToken< SiPixelTemplateDBObject, SiPixelTemplateDBObjectESProducerRcd > siPixelTemplateEsToken_

◆ dqmEndRun()

void SiPixelLorentzAnglePCLWorker::dqmEndRun ( edm::Run const &  run,
edm::EventSetup const &  iSetup 
)
private

Definition at line 949 of file SiPixelLorentzAnglePCLWorker.cc.

References hFile_, and notInPCL_.

949  {
950  if (notInPCL_) {
951  hFile_->cd();
952  hFile_->Write();
953  hFile_->Close();
954  }
955 }

◆ fillDescriptions()

void SiPixelLorentzAnglePCLWorker::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 999 of file SiPixelLorentzAnglePCLWorker.cc.

References edm::ConfigurationDescriptions::addWithDefaultLabel(), submitPVResolutionJobs::desc, ProducerED_cfi::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

999  {
1001  desc.setComment("Worker module of the SiPixel Lorentz Angle PCL monitoring workflow");
1002  desc.add<std::string>("analysisType", "GrazingAngle")
1003  ->setComment("analysis type - GrazingAngle (default) or MinimumClusterSize");
1004  desc.add<std::string>("folder", "AlCaReco/SiPixelLorentzAngle")->setComment("directory of PCL Worker output");
1005  desc.add<bool>("notInPCL", false)->setComment("create TTree (true) or not (false)");
1006  desc.add<std::string>("fileName", "testrun.root")->setComment("name of the TTree file if notInPCL = true");
1007  desc.add<std::vector<std::string>>("newmodulelist", {})->setComment("the list of DetIds for new sensors");
1008  desc.add<edm::InputTag>("src", edm::InputTag("TrackRefitter"))->setComment("input track collections");
1009  desc.add<double>("ptMin", 3.)->setComment("minimum pt on tracks");
1010  desc.add<double>("normChi2Max", 2.)->setComment("maximum reduced chi squared");
1011  desc.add<std::vector<int>>("clustSizeYMin", {4, 3, 3, 2})
1012  ->setComment("minimum cluster size on Y axis for all Barrel Layers");
1013  desc.add<int>("clustSizeXMax", 5)->setComment("maximum cluster size on X axis");
1014  desc.add<double>("residualMax", 0.005)->setComment("maximum residual");
1015  desc.add<double>("clustChargeMaxPerLength", 50000)
1016  ->setComment("maximum cluster charge per unit length of pixel depth (z)");
1017  desc.add<int>("binsDepth", 50)->setComment("# bins for electron production depth axis");
1018  desc.add<int>("binsDrift", 100)->setComment("# bins for electron drift axis");
1019  descriptions.addWithDefaultLabel(desc);
1020 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)

◆ fillPix()

const Pixinfo SiPixelLorentzAnglePCLWorker::fillPix ( const SiPixelCluster LocPix,
const PixelTopology topol 
) const
private

Definition at line 958 of file SiPixelLorentzAnglePCLWorker.cc.

References Pixinfo::adc, Pixinfo::col, Topology::localPosition(), Pixinfo::npix, SiPixelCluster::pixels(), Pixinfo::row, Pixinfo::x, PV3DBase< T, PVType, FrameType >::x(), Pixinfo::y, and PV3DBase< T, PVType, FrameType >::y().

Referenced by analyze().

958  {
959  Pixinfo pixinfo;
960  const std::vector<SiPixelCluster::Pixel>& pixvector = LocPix.pixels();
961  pixinfo.npix = 0;
962  for (std::vector<SiPixelCluster::Pixel>::const_iterator itPix = pixvector.begin(); itPix != pixvector.end();
963  itPix++) {
964  pixinfo.row[pixinfo.npix] = itPix->x;
965  pixinfo.col[pixinfo.npix] = itPix->y;
966  pixinfo.adc[pixinfo.npix] = itPix->adc;
967  LocalPoint lp = topol->localPosition(MeasurementPoint(itPix->x + 0.5, itPix->y + 0.5));
968  pixinfo.x[pixinfo.npix] = lp.x();
969  pixinfo.y[pixinfo.npix] = lp.y();
970  pixinfo.npix++;
971  }
972  return pixinfo;
973 }
float x[maxpix]
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
float y[maxpix]
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
const std::vector< Pixel > pixels() const
float adc[maxpix]
float col[maxpix]
float row[maxpix]

◆ surface_deformation()

const std::pair< LocalPoint, LocalPoint > SiPixelLorentzAnglePCLWorker::surface_deformation ( const PixelTopology topol,
TrajectoryStateOnSurface tsos,
const SiPixelRecHit recHitPix 
) const
private

Definition at line 976 of file SiPixelLorentzAnglePCLWorker.cc.

References LocalTrajectoryParameters::dxdz(), LocalTrajectoryParameters::dydz(), TrajectoryStateOnSurface::localParameters(), BaseTrackerRecHit::localPosition(), TrajectoryStateOnSurface::localPosition(), Topology::localPosition(), and PixelTopology::pixel().

Referenced by analyze().

977  {
978  LocalPoint trackposition = tsos.localPosition();
979  const LocalTrajectoryParameters& ltp = tsos.localParameters();
980  const Topology::LocalTrackAngles localTrackAngles(ltp.dxdz(), ltp.dydz());
981 
982  std::pair<float, float> pixels_track = topol->pixel(trackposition, localTrackAngles);
983  std::pair<float, float> pixels_rechit = topol->pixel(recHitPix->localPosition(), localTrackAngles);
984 
985  LocalPoint lp_track = topol->localPosition(MeasurementPoint(pixels_track.first, pixels_track.second));
986 
987  LocalPoint lp_rechit = topol->localPosition(MeasurementPoint(pixels_rechit.first, pixels_rechit.second));
988 
989  std::pair<LocalPoint, LocalPoint> lps = std::make_pair(lp_track, lp_rechit);
990  return lps;
991 }
virtual std::pair< float, float > pixel(const LocalPoint &p) const =0
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
const LocalTrajectoryParameters & localParameters() const
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
LocalPoint localPosition() const override

Member Data Documentation

◆ analysisType_

LorentzAngleAnalysisTypeEnum SiPixelLorentzAnglePCLWorker::analysisType_
private

Definition at line 133 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and bookHistograms().

◆ bladeF_

int SiPixelLorentzAnglePCLWorker::bladeF_
private

Definition at line 167 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ bx_

int SiPixelLorentzAnglePCLWorker::bx_
private

Definition at line 143 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ chi2_

double SiPixelLorentzAnglePCLWorker::chi2_
private

Definition at line 152 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ clust_

Clust SiPixelLorentzAnglePCLWorker::clust_
private

Definition at line 156 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ clustChargeMaxPerLength_

double SiPixelLorentzAnglePCLWorker::clustChargeMaxPerLength_
private

Definition at line 186 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze().

◆ clustF_

Clust SiPixelLorentzAnglePCLWorker::clustF_
private

Definition at line 172 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ clustSizeXMax_

int SiPixelLorentzAnglePCLWorker::clustSizeXMax_
private

Definition at line 184 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze().

◆ clustSizeYMin_

std::vector<int> SiPixelLorentzAnglePCLWorker::clustSizeYMin_
private

Definition at line 183 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze().

◆ diskF_

int SiPixelLorentzAnglePCLWorker::diskF_
private

Definition at line 166 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ eta_

float SiPixelLorentzAnglePCLWorker::eta_
private

Definition at line 150 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ event_

long int SiPixelLorentzAnglePCLWorker::event_
private

Definition at line 141 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ filename_

std::string SiPixelLorentzAnglePCLWorker::filename_
private

Definition at line 136 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by SiPixelLorentzAnglePCLWorker().

◆ folder_

std::string SiPixelLorentzAnglePCLWorker::folder_
private

Definition at line 134 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by bookHistograms().

◆ geomEsToken_

edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> SiPixelLorentzAnglePCLWorker::geomEsToken_
private

Definition at line 195 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by dqmBeginRun().

◆ geomPerEventEsToken_

edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> SiPixelLorentzAnglePCLWorker::geomPerEventEsToken_
private

Definition at line 200 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze().

◆ hFile_

std::unique_ptr<TFile> SiPixelLorentzAnglePCLWorker::hFile_
private

Definition at line 190 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by dqmEndRun(), and SiPixelLorentzAnglePCLWorker().

◆ hist_depth_

int SiPixelLorentzAnglePCLWorker::hist_depth_
private

Definition at line 187 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by bookHistograms().

◆ hist_drift_

int SiPixelLorentzAnglePCLWorker::hist_drift_
private

Definition at line 188 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by bookHistograms().

◆ iHists

SiPixelLorentzAngleCalibrationHistograms SiPixelLorentzAnglePCLWorker::iHists
private

Definition at line 126 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), bookHistograms(), and dqmBeginRun().

◆ isflipped_

int SiPixelLorentzAnglePCLWorker::isflipped_
private

Definition at line 148 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ ladder_

int SiPixelLorentzAnglePCLWorker::ladder_
private

Definition at line 146 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ layer_

int SiPixelLorentzAnglePCLWorker::layer_
private

Definition at line 147 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ lumiblock_

int SiPixelLorentzAnglePCLWorker::lumiblock_
private

Definition at line 142 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ magneticFieldToken_

edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> SiPixelLorentzAnglePCLWorker::magneticFieldToken_
private

Definition at line 201 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze().

◆ module_

int SiPixelLorentzAnglePCLWorker::module_
private

Definition at line 145 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ moduleF_

int SiPixelLorentzAnglePCLWorker::moduleF_
private

Definition at line 169 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ ndof_

double SiPixelLorentzAnglePCLWorker::ndof_
private

Definition at line 153 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ newmodulelist_

std::vector<std::string> SiPixelLorentzAnglePCLWorker::newmodulelist_
private

Definition at line 137 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by dqmBeginRun().

◆ normChi2Max_

double SiPixelLorentzAnglePCLWorker::normChi2Max_
private

Definition at line 182 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze().

◆ notInPCL_

bool SiPixelLorentzAnglePCLWorker::notInPCL_
private

◆ orbit_

int SiPixelLorentzAnglePCLWorker::orbit_
private

Definition at line 144 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ panelF_

int SiPixelLorentzAnglePCLWorker::panelF_
private

Definition at line 168 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ phi_

float SiPixelLorentzAnglePCLWorker::phi_
private

Definition at line 151 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ pixinfo_

Pixinfo SiPixelLorentzAnglePCLWorker::pixinfo_
private

Definition at line 154 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ pixinfoF_

Pixinfo SiPixelLorentzAnglePCLWorker::pixinfoF_
private

Definition at line 170 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ pt_

float SiPixelLorentzAnglePCLWorker::pt_
private

Definition at line 149 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ ptmin_

double SiPixelLorentzAnglePCLWorker::ptmin_
private

Definition at line 181 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze().

◆ qScale_

float SiPixelLorentzAnglePCLWorker::qScale_
private

Definition at line 161 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ qScaleF_

float SiPixelLorentzAnglePCLWorker::qScaleF_
private

Definition at line 177 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ rechit_

Rechit SiPixelLorentzAnglePCLWorker::rechit_
private

Definition at line 157 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ rechitCorr_

Rechit SiPixelLorentzAnglePCLWorker::rechitCorr_
private

Definition at line 158 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ rechitCorrF_

Rechit SiPixelLorentzAnglePCLWorker::rechitCorrF_
private

Definition at line 174 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ rechitF_

Rechit SiPixelLorentzAnglePCLWorker::rechitF_
private

Definition at line 173 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ residualMax_

double SiPixelLorentzAnglePCLWorker::residualMax_
private

Definition at line 185 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze().

◆ rQmQt_

float SiPixelLorentzAnglePCLWorker::rQmQt_
private

Definition at line 162 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ rQmQtF_

float SiPixelLorentzAnglePCLWorker::rQmQtF_
private

Definition at line 178 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ run_

int SiPixelLorentzAnglePCLWorker::run_
private

Definition at line 140 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ sideF_

int SiPixelLorentzAnglePCLWorker::sideF_
private

Definition at line 165 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ simhit_

Hit SiPixelLorentzAnglePCLWorker::simhit_
private

Definition at line 155 of file SiPixelLorentzAnglePCLWorker.cc.

◆ simhitF_

Hit SiPixelLorentzAnglePCLWorker::simhitF_
private

Definition at line 171 of file SiPixelLorentzAnglePCLWorker.cc.

◆ SiPixelLorentzAngleTreeBarrel_

std::unique_ptr<TTree> SiPixelLorentzAnglePCLWorker::SiPixelLorentzAngleTreeBarrel_
private

Definition at line 191 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ SiPixelLorentzAngleTreeForward_

std::unique_ptr<TTree> SiPixelLorentzAnglePCLWorker::SiPixelLorentzAngleTreeForward_
private

Definition at line 192 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ siPixelTemplateEsToken_

edm::ESGetToken<SiPixelTemplateDBObject, SiPixelTemplateDBObjectESProducerRcd> SiPixelLorentzAnglePCLWorker::siPixelTemplateEsToken_
private

Definition at line 197 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by dqmBeginRun().

◆ siPixelTemplateStoreEsToken_

edm::ESGetToken<std::vector<SiPixelTemplateStore>, SiPixelTemplateDBObjectESProducerRcd> SiPixelLorentzAnglePCLWorker::siPixelTemplateStoreEsToken_
private

Definition at line 198 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by dqmBeginRun().

◆ t_trajTrack

edm::EDGetTokenT<TrajTrackAssociationCollection> SiPixelLorentzAnglePCLWorker::t_trajTrack
private

Definition at line 204 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ templateDBobject_

const SiPixelTemplateDBObject* SiPixelLorentzAnglePCLWorker::templateDBobject_
private

Definition at line 130 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and dqmBeginRun().

◆ thePixelTemp_

const std::vector<SiPixelTemplateStore>* SiPixelLorentzAnglePCLWorker::thePixelTemp_ = nullptr
private

Definition at line 131 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and dqmBeginRun().

◆ topoEsToken_

edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> SiPixelLorentzAnglePCLWorker::topoEsToken_
private

Definition at line 196 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by dqmBeginRun().

◆ topoPerEventEsToken_

edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> SiPixelLorentzAnglePCLWorker::topoPerEventEsToken_
private

Definition at line 199 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze().

◆ trackhit_

Hit SiPixelLorentzAnglePCLWorker::trackhit_
private

Definition at line 155 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ trackhitCorrX_

float SiPixelLorentzAnglePCLWorker::trackhitCorrX_
private

Definition at line 159 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ trackhitCorrXF_

float SiPixelLorentzAnglePCLWorker::trackhitCorrXF_
private

Definition at line 175 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ trackhitCorrY_

float SiPixelLorentzAnglePCLWorker::trackhitCorrY_
private

Definition at line 160 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ trackhitCorrYF_

float SiPixelLorentzAnglePCLWorker::trackhitCorrYF_
private

Definition at line 176 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ trackhitF_

Hit SiPixelLorentzAnglePCLWorker::trackhitF_
private

Definition at line 171 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by analyze(), and SiPixelLorentzAnglePCLWorker().

◆ watchSiPixelTemplateRcd_

edm::ESWatcher<SiPixelTemplateDBObjectESProducerRcd> SiPixelLorentzAnglePCLWorker::watchSiPixelTemplateRcd_
private

Definition at line 129 of file SiPixelLorentzAnglePCLWorker.cc.

Referenced by dqmBeginRun().