CMS 3D CMS Logo

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

#include <CalibTracker/SiPixelLorentzAngle/src/SiPixelLorentzAnglePCLHarvester.cc>

Inheritance diagram for SiPixelLorentzAnglePCLHarvester:
DQMEDHarvester edm::one::EDProducer< edm::EndLuminosityBlockProducer, edm::EndRunProducer, edm::EndProcessBlockProducer, edm::one::WatchLuminosityBlocks, edm::one::WatchRuns, edm::one::SharedResources, edm::Accumulator > edm::one::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

void beginRun (const edm::Run &, const edm::EventSetup &) override
 
 SiPixelLorentzAnglePCLHarvester (const edm::ParameterSet &)
 
 ~SiPixelLorentzAnglePCLHarvester () override=default
 
- Public Member Functions inherited from DQMEDHarvester
void accumulate (edm::Event const &ev, edm::EventSetup const &es) final
 
void beginJob () override
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &) final
 
void beginRun (edm::Run const &run, edm::EventSetup const &) override
 
virtual void dqmAnalyze (DQMStore::IBooker &, DQMStore::IGetter &, edm::Event const &, edm::EventSetup const &)
 
 DQMEDHarvester (edm::ParameterSet const &iConfig)
 
 DQMEDHarvester ()
 
virtual void dqmEndLuminosityBlock (DQMStore::IBooker &, DQMStore::IGetter &, edm::LuminosityBlock const &, edm::EventSetup const &)
 
virtual void dqmEndRun (DQMStore::IBooker &, DQMStore::IGetter &, edm::Run const &, edm::EventSetup const &)
 
void endLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &) final
 
void endLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &es) final
 
void endProcessBlockProduce (edm::ProcessBlock &) final
 
void endRun (edm::Run const &, edm::EventSetup const &) override
 
void endRunProduce (edm::Run &run, edm::EventSetup const &es) final
 
 ~DQMEDHarvester () override=default
 
- Public Member Functions inherited from edm::one::EDProducer< edm::EndLuminosityBlockProducer, edm::EndRunProducer, edm::EndProcessBlockProducer, edm::one::WatchLuminosityBlocks, edm::one::WatchRuns, edm::one::SharedResources, edm::Accumulator >
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
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
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsInputProcessBlocks () const final
 
bool wantsProcessBlocks () const final
 
- Public Member Functions inherited from edm::one::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
std::vector< bool > const & recordProvenanceList () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
TypeLabelList const & typeLabelList () const
 used by the fwk to register the list of products of this module More...
 
 ~ProducerBase () noexcept(false) override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESProxyIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex > const & esGetTokenRecordIndicesVector (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::array< std::vector< ModuleDescription const *> *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, ProductRegistry const &preg, std::map< std::string, ModuleDescription const *> const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &)
 
- Static Public Member Functions inherited from edm::one::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Private Member Functions

void dqmEndJob (DQMStore::IBooker &, DQMStore::IGetter &) override
 
void endRun (const edm::Run &, const edm::EventSetup &) override
 
void findMean (MonitorElement *h_drift_depth_adc_slice_, int i, int i_ring)
 
SiPixelLAHarvest::fitResults fitAndStore (std::shared_ptr< SiPixelLorentzAngle > theLA, int i_idx, int i_lay, int i_mod)
 
bool isFitGood (SiPixelLAHarvest::fitResults &res)
 

Private Attributes

const SiPixelLorentzAnglecurrentLorentzAngle_
 
const bool doChebyshevFit_
 
const std::string dqmDir_
 
std::unique_ptr< TF1 > f1_
 
const double fitChi2Cut_
 
const std::vector< double > fitRange_
 
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecordgeomEsToken_
 
SiPixelLorentzAngleCalibrationHistograms hists_
 
edm::ESGetToken< MagneticField, IdealMagneticFieldRecordmagneticFieldToken_
 
const int minHitsCut_
 
std::vector< std::string > newmodulelist_
 
const int order_
 
const std::string recordName_
 
edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleRcdsiPixelLAEsToken_
 
std::pair< double, double > theFitRange_ {5., 280.}
 
float theMagField_ {0.f}
 
std::unique_ptr< TrackerTopologytheTrackerTopology_
 
edm::ESGetToken< TrackerTopology, TrackerTopologyRcdtopoEsTokenBR_
 
edm::ESGetToken< TrackerTopology, TrackerTopologyRcdtopoEsTokenER_
 
float width_
 

Static Private Attributes

static constexpr float teslaToInverseGeV_ = 2.99792458e-3f
 

Additional Inherited Members

- Public Types inherited from DQMEDHarvester
typedef dqm::harvesting::DQMStore DQMStore
 
typedef dqm::harvesting::MonitorElement MonitorElement
 
- Public Types inherited from edm::one::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
template<typename T >
using BranchAliasSetterT = ProductRegistryHelper::BranchAliasSetterT< T >
 
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex > >
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Protected Member Functions inherited from edm::ProducerBase
template<Transition Tr = Transition::Event>
auto produces (std::string instanceName) noexcept
 declare what type of product will make and with which optional label More...
 
template<Transition B>
BranchAliasSetter produces (const TypeID &id, std::string instanceName=std::string(), bool recordProvenance=true)
 
template<BranchType B>
BranchAliasSetter produces (const TypeID &id, std::string instanceName=std::string(), bool recordProvenance=true)
 
BranchAliasSetter produces (const TypeID &id, std::string instanceName=std::string(), bool recordProvenance=true)
 
template<typename ProductType , Transition B>
BranchAliasSetterT< ProductType > produces (std::string instanceName)
 
template<class ProductType >
BranchAliasSetterT< ProductType > produces ()
 
template<typename ProductType , BranchType B>
BranchAliasSetterT< ProductType > produces (std::string instanceName)
 
template<typename ProductType , BranchType B>
BranchAliasSetterT< ProductType > produces ()
 
template<class ProductType >
BranchAliasSetterT< ProductType > produces (std::string instanceName)
 
template<typename ProductType , Transition B>
BranchAliasSetterT< ProductType > produces ()
 
template<Transition Tr = Transition::Event>
auto produces () noexcept
 
ProducesCollector producesCollector ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
void resetItemsToGetFrom (BranchType iType)
 
- Protected Attributes inherited from DQMEDHarvester
DQMStoredqmstore_
 
edm::GetterOfProducts< DQMTokenjobmegetter_
 
edm::EDPutTokenT< DQMTokenjobToken_
 
edm::GetterOfProducts< DQMTokenlumimegetter_
 
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::GetterOfProducts< DQMTokenrunmegetter_
 
edm::EDPutTokenT< DQMTokenrunToken_
 

Detailed Description

Description: reads the intermediate ALCAPROMPT DQMIO-like dataset and performs the fitting of the SiPixel Lorentz Angle in the Prompt Calibration Loop Implementation: Reads the 2D histograms of the drift vs depth created by SiPixelLorentzAnglePCLWorker modules and generates 1D profiles which are then fit with a 5th order polinomial. The extracted value of the tan(theta_L)/B are stored in an output sqlite file which is then uploaded to the conditions database

Definition at line 104 of file SiPixelLorentzAnglePCLHarvester.cc.

Constructor & Destructor Documentation

◆ SiPixelLorentzAnglePCLHarvester()

SiPixelLorentzAnglePCLHarvester::SiPixelLorentzAnglePCLHarvester ( const edm::ParameterSet iConfig)

Definition at line 146 of file SiPixelLorentzAnglePCLHarvester.cc.

References Exception, fitRange_, edm::Service< T >::isAvailable(), and theFitRange_.

147  : geomEsToken_(esConsumes<edm::Transition::BeginRun>()),
148  topoEsTokenBR_(esConsumes<edm::Transition::BeginRun>()),
149  topoEsTokenER_(esConsumes<edm::Transition::EndRun>()),
150  siPixelLAEsToken_(esConsumes<edm::Transition::BeginRun>()),
151  magneticFieldToken_(esConsumes<edm::Transition::BeginRun>()),
152  newmodulelist_(iConfig.getParameter<std::vector<std::string>>("newmodulelist")),
153  dqmDir_(iConfig.getParameter<std::string>("dqmDir")),
154  doChebyshevFit_(iConfig.getParameter<bool>("doChebyshevFit")),
155  order_(iConfig.getParameter<int>("order")),
156  fitChi2Cut_(iConfig.getParameter<double>("fitChi2Cut")),
157  fitRange_(iConfig.getParameter<std::vector<double>>("fitRange")),
158  minHitsCut_(iConfig.getParameter<int>("minHitsCut")),
159  recordName_(iConfig.getParameter<std::string>("record")) {
160  // initialize the fit range
161  if (fitRange_.size() == 2) {
162  theFitRange_.first = fitRange_[0];
163  theFitRange_.second = fitRange_[1];
164  } else {
165  throw cms::Exception("SiPixelLorentzAnglePCLHarvester") << "Too many fit range parameters specified";
166  }
167 
168  // first ensure DB output service is available
170  if (!poolDbService.isAvailable())
171  throw cms::Exception("SiPixelLorentzAnglePCLHarvester") << "PoolDBService required";
172 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoEsTokenBR_
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoEsTokenER_
edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleRcd > siPixelLAEsToken_
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomEsToken_
bool isAvailable() const
Definition: Service.h:40
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_

◆ ~SiPixelLorentzAnglePCLHarvester()

SiPixelLorentzAnglePCLHarvester::~SiPixelLorentzAnglePCLHarvester ( )
overridedefault

Member Function Documentation

◆ beginRun()

void SiPixelLorentzAnglePCLHarvester::beginRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)
override

Definition at line 175 of file SiPixelLorentzAnglePCLHarvester.cc.

References PixelEndcapName::bladeName(), Surface::bounds(), SiPixelLorentzAngleCalibrationHistograms::BPixnewDetIds_, SiPixelLorentzAngleCalibrationHistograms::BPixnewLayer_, SiPixelLorentzAngleCalibrationHistograms::BPixnewModule_, SiPixelLorentzAngleCalibrationHistograms::BPixnewmodulename_, submitPVResolutionJobs::count, currentLorentzAngle_, SiPixelLorentzAngleCalibrationHistograms::detIdsList, PixelEndcapName::diskName(), spr::find(), SiPixelLorentzAngleCalibrationHistograms::FPixnewBlade_, SiPixelLorentzAngleCalibrationHistograms::FPixnewDetIds_, SiPixelLorentzAngleCalibrationHistograms::FPixnewDisk_, SiPixelLorentzAngleCalibrationHistograms::FPixnewmodulename_, GeomDet::geographicalId(), relativeConstraints::geom, geomEsToken_, edm::EventSetup::getData(), PixelBarrelName::getDetId(), PixelEndcapName::getDetId(), hists_, mps_fire::i, l1ctLayer2EG_cff::id, MagneticField::inverseBzAtOriginInGeV(), pixelTopology::layer, PixelBarrelName::layerName(), LogDebug, magneticFieldToken_, genParticles_cff::map, callgraph::module, PixelBarrelName::moduleName(), newmodulelist_, SiPixelLorentzAngleCalibrationHistograms::nLadders_, SiPixelLorentzAngleCalibrationHistograms::nlay, SiPixelLorentzAngleCalibrationHistograms::nModules_, PixelSubdetector::PixelBarrel, TrackerTopology::pxbLayer(), TrackerTopology::pxbModule(), DetId::rawId(), siPixelLAEsToken_, GeomDet::surface(), teslaToInverseGeV_, theMagField_, Bounds::thickness(), topoEsTokenBR_, parallelization::uint, and width_.

175  {
176  // geometry
177  const TrackerGeometry* geom = &iSetup.getData(geomEsToken_);
178  const TrackerTopology* tTopo = &iSetup.getData(topoEsTokenBR_);
179 
180  const MagneticField* magField = &iSetup.getData(magneticFieldToken_);
182 
183  // B-field value
184  // inverseBzAtOriginInGeV() returns the inverse of field z component for this map in GeV
185  // for the conversion please consult https://github.com/cms-sw/cmssw/blob/master/MagneticField/Engine/src/MagneticField.cc#L17
186  // theInverseBzAtOriginInGeV = 1.f / (at0z * 2.99792458e-3f);
187  // ==> at0z = 1.f / (theInverseBzAtOriginInGeV * 2.99792458e-3f)
188 
190 
192  hists_.nlay = geom->numberOfLayers(PixelSubdetector::PixelBarrel);
193  hists_.nModules_.resize(hists_.nlay);
194  hists_.nLadders_.resize(hists_.nlay);
195  for (int i = 0; i < hists_.nlay; i++) {
196  hists_.nModules_[i] = map.getPXBModules(i + 1);
197  hists_.nLadders_[i] = map.getPXBLadders(i + 1);
198  }
199 
200  // list of modules already filled, then return (we already entered here)
201  if (!hists_.BPixnewDetIds_.empty() || !hists_.FPixnewDetIds_.empty())
202  return;
203 
204  if (!newmodulelist_.empty()) {
205  for (auto const& modulename : newmodulelist_) {
206  if (modulename.find("BPix_") != std::string::npos) {
207  PixelBarrelName bn(modulename, true);
208  const auto& detId = bn.getDetId(tTopo);
209  hists_.BPixnewmodulename_.push_back(modulename);
210  hists_.BPixnewDetIds_.push_back(detId.rawId());
211  hists_.BPixnewModule_.push_back(bn.moduleName());
212  hists_.BPixnewLayer_.push_back(bn.layerName());
213  } else if (modulename.find("FPix_") != std::string::npos) {
214  PixelEndcapName en(modulename, true);
215  const auto& detId = en.getDetId(tTopo);
216  hists_.FPixnewmodulename_.push_back(modulename);
217  hists_.FPixnewDetIds_.push_back(detId.rawId());
218  hists_.FPixnewDisk_.push_back(en.diskName());
219  hists_.FPixnewBlade_.push_back(en.bladeName());
220  }
221  }
222  }
223 
224  uint count = 0;
225  for (const auto& id : hists_.BPixnewDetIds_) {
226  LogDebug("SiPixelLorentzAnglePCLHarvester") << id;
227  count++;
228  }
229  LogDebug("SiPixelLorentzAnglePCLHarvester") << "Stored a total of " << count << " new detIds.";
230 
231  // list of modules already filled, return (we already entered here)
232  if (!hists_.detIdsList.empty())
233  return;
234 
235  std::vector<uint32_t> treatedIndices;
236 
237  for (const auto& det : geom->detsPXB()) {
238  const PixelGeomDetUnit* pixelDet = dynamic_cast<const PixelGeomDetUnit*>(det);
239  width_ = pixelDet->surface().bounds().thickness();
240  const auto& layer = tTopo->pxbLayer(pixelDet->geographicalId());
241  const auto& module = tTopo->pxbModule(pixelDet->geographicalId());
242  int i_index = module + (layer - 1) * hists_.nModules_[layer - 1];
243 
244  uint32_t rawId = pixelDet->geographicalId().rawId();
245 
246  // if the detId is already accounted for in the special class, do not attach it
247  if (std::find(hists_.BPixnewDetIds_.begin(), hists_.BPixnewDetIds_.end(), rawId) != hists_.BPixnewDetIds_.end())
248  continue;
249 
250  if (std::find(treatedIndices.begin(), treatedIndices.end(), i_index) != treatedIndices.end()) {
251  hists_.detIdsList[i_index].push_back(rawId);
252  } else {
253  hists_.detIdsList.insert(std::pair<uint32_t, std::vector<uint32_t>>(i_index, {rawId}));
254  treatedIndices.push_back(i_index);
255  }
256  }
257 
258  count = 0;
259  for (const auto& i : treatedIndices) {
260  for (const auto& id : hists_.detIdsList[i]) {
261  LogDebug("SiPixelLorentzAnglePCLHarvester") << id;
262  count++;
263  };
264  }
265  LogDebug("SiPixelLorentzAnglePCLHarvester") << "Stored a total of " << count << " detIds.";
266 }
unsigned int pxbLayer(const DetId &id) const
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoEsTokenBR_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
float inverseBzAtOriginInGeV() const
The inverse of field z component for this map in GeV.
Definition: MagneticField.h:52
constexpr std::array< uint8_t, layerIndexSize< TrackerTraits > > layer
virtual float thickness() const =0
edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleRcd > siPixelLAEsToken_
std::unordered_map< uint32_t, std::vector< uint32_t > > detIdsList
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomEsToken_
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
SiPixelLorentzAngleCalibrationHistograms hists_
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
unsigned int pxbModule(const DetId &id) const
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_
#define LogDebug(id)
const Bounds & bounds() const
Definition: Surface.h:87

◆ dqmEndJob()

void SiPixelLorentzAnglePCLHarvester::dqmEndJob ( DQMStore::IBooker iBooker,
DQMStore::IGetter iGetter 
)
overrideprivatevirtual

Implements DQMEDHarvester.

Definition at line 276 of file SiPixelLorentzAnglePCLHarvester.cc.

References newFWLiteAna::bin, dqm::implementation::IBooker::book1D(), dqm::implementation::IBooker::book2D(), SiPixelLorentzAngleCalibrationHistograms::BPixnewDetIds_, SiPixelLorentzAngleCalibrationHistograms::BPixnewLayer_, SiPixelLorentzAngleCalibrationHistograms::BPixnewModule_, SiPixelLorentzAngleCalibrationHistograms::BPixnewmodulename_, dqm::implementation::NavigatorBase::cd(), currentLorentzAngle_, cond::service::PoolDBOutputService::currentTime(), dqmDir_, cppFunctionSkipper::exception, f, dqm::impl::MonitorElement::Fill(), findMean(), fitAndStore(), dqm-mbProfile::format, dqm::implementation::IGetter::get(), SiPixelLorentzAngle::getLorentzAngle(), SiPixelLorentzAngle::getLorentzAngles(), dqm::impl::MonitorElement::getNbinsX(), dqm::impl::MonitorElement::getTH1(), SiPixelLorentzAngleCalibrationHistograms::h2_byLayerDiff_, SiPixelLorentzAngleCalibrationHistograms::h2_byLayerLA_, SiPixelLorentzAngleCalibrationHistograms::h_bySectChi2_, SiPixelLorentzAngleCalibrationHistograms::h_bySectCovMatrixStatus_, SiPixelLorentzAngleCalibrationHistograms::h_bySectDeltaLA_, SiPixelLorentzAngleCalibrationHistograms::h_bySectDriftError_, SiPixelLorentzAngleCalibrationHistograms::h_bySectFitQuality_, SiPixelLorentzAngleCalibrationHistograms::h_bySectFitStatus_, SiPixelLorentzAngleCalibrationHistograms::h_bySectLA_, SiPixelLorentzAngleCalibrationHistograms::h_bySectMeasLA_, SiPixelLorentzAngleCalibrationHistograms::h_bySectOccupancy_, SiPixelLorentzAngleCalibrationHistograms::h_bySectRejectLA_, SiPixelLorentzAngleCalibrationHistograms::h_bySectSetLA_, SiPixelLorentzAngleCalibrationHistograms::h_drift_depth_, SiPixelLorentzAngleCalibrationHistograms::h_drift_depth_adc2_, SiPixelLorentzAngleCalibrationHistograms::h_drift_depth_adc_, SiPixelLorentzAngleCalibrationHistograms::h_drift_depth_noadc_, SiPixelLorentzAngleCalibrationHistograms::h_mean_, hists_, mps_fire::i, caHitNtupletGeneratorKernels::if(), createfilelist::int, edm::Service< T >::isAvailable(), dqmiolumiharvest::j, SummaryClient_cfi::labels, PVValHelper::ladder, pixelTopology::layer, SiStripSimParameters_cfi::LorentzAngle, callgraph::module, SiPixelLorentzAngleCalibrationHistograms::nLadders_, SiPixelLorentzAngleCalibrationHistograms::nlay, SiPixelLorentzAngleCalibrationHistograms::nModules_, PixelSubdetector::PixelBarrel, recordName_, dqm::impl::MonitorElement::setBinLabel(), dqm::implementation::NavigatorBase::setCurrentFolder(), AlCaHLTBitMon_QueryRunRegistry::string, theTrackerTopology_, HcalDetIdTransform::transform(), relativeConstraints::value, cms::Exception::what(), and cond::service::PoolDBOutputService::writeOneIOV().

276  {
277  // go in the right directory
278  iGetter.cd();
279  iGetter.setCurrentFolder(dqmDir_);
280 
281  // fetch the 2D histograms
282  for (int i_layer = 1; i_layer <= hists_.nlay; i_layer++) {
283  const auto& prefix_ = fmt::sprintf("%s/BPix/BPixLayer%i", dqmDir_, i_layer);
284  for (int i_module = 1; i_module <= hists_.nModules_[i_layer - 1]; i_module++) {
285  int i_index = i_module + (i_layer - 1) * hists_.nModules_[i_layer - 1];
286 
287  hists_.h_drift_depth_[i_index] =
288  iGetter.get(fmt::format("{}/h_drift_depth_layer{}_module{}", prefix_, i_layer, i_module));
289 
290  if (hists_.h_drift_depth_[i_index] == nullptr) {
291  edm::LogError("SiPixelLorentzAnglePCLHarvester::dqmEndJob")
292  << "Failed to retrieve electron drift over depth for layer " << i_layer << ", module " << i_module << ".";
293  continue;
294  }
295 
296  hists_.h_drift_depth_adc_[i_index] =
297  iGetter.get(fmt::format("{}/h_drift_depth_adc_layer{}_module{}", prefix_, i_layer, i_module));
298 
299  hists_.h_drift_depth_adc2_[i_index] =
300  iGetter.get(fmt::format("{}/h_drift_depth_adc2_layer{}_module{}", prefix_, i_layer, i_module));
301 
302  hists_.h_drift_depth_noadc_[i_index] =
303  iGetter.get(fmt::format("{}/h_drift_depth_noadc_layer{}_module{}", prefix_, i_layer, i_module));
304 
305  hists_.h_mean_[i_index] = iGetter.get(fmt::format("{}/h_mean_layer{}_module{}", dqmDir_, i_layer, i_module));
306 
307  hists_.h_drift_depth_[i_index]->divide(
308  hists_.h_drift_depth_adc_[i_index], hists_.h_drift_depth_noadc_[i_index], 1., 1., "");
309  }
310  }
311 
312  // fetch the new modules 2D histograms
313  for (int i = 0; i < (int)hists_.BPixnewDetIds_.size(); i++) {
314  int new_index = i + 1 + hists_.nModules_[hists_.nlay - 1] + (hists_.nlay - 1) * hists_.nModules_[hists_.nlay - 1];
315 
316  hists_.h_drift_depth_[new_index] = iGetter.get(
317  fmt::format("{}/h_BPixnew_drift_depth_{}", dqmDir_ + "/BPix/NewModules", hists_.BPixnewmodulename_[i]));
318 
319  if (hists_.h_drift_depth_[new_index] == nullptr) {
320  edm::LogError("SiPixelLorentzAnglePCLHarvester")
321  << "Failed to retrieve electron drift over depth for new module " << hists_.BPixnewmodulename_[i] << ".";
322  continue;
323  }
324 
325  hists_.h_drift_depth_adc_[new_index] = iGetter.get(
326  fmt::format("{}/h_BPixnew_drift_depth_adc_{}", dqmDir_ + "/BPix/NewModules", hists_.BPixnewmodulename_[i]));
327 
328  hists_.h_drift_depth_adc2_[new_index] = iGetter.get(
329  fmt::format("{}/h_BPixnew_drift_depth_adc2_{}", dqmDir_ + "/BPix/NewModules", hists_.BPixnewmodulename_[i]));
330 
331  hists_.h_drift_depth_noadc_[new_index] = iGetter.get(
332  fmt::format("{}/h_BPixnew_drift_depth_noadc_{}", dqmDir_ + "/BPix/NewModules", hists_.BPixnewmodulename_[i]));
333 
334  hists_.h_mean_[new_index] = iGetter.get(fmt::format("{}/h_BPixnew_mean_{}", dqmDir_, hists_.BPixnewmodulename_[i]));
335 
336  hists_.h_drift_depth_[new_index]->divide(
337  hists_.h_drift_depth_adc_[new_index], hists_.h_drift_depth_noadc_[new_index], 1., 1., "");
338  }
339 
340  hists_.h_bySectOccupancy_ = iGetter.get(fmt::format("{}/h_bySectorOccupancy", dqmDir_ + "/SectorMonitoring"));
341  if (hists_.h_bySectOccupancy_ == nullptr) {
342  edm::LogError("SiPixelLorentzAnglePCLHarvester") << "Failed to retrieve the hit on track occupancy.";
343  return;
344  }
345 
346  int hist_drift_;
347  int hist_depth_;
348  double min_drift_;
349  double max_drift_;
350 
351  if (hists_.h_drift_depth_adc_[1] != nullptr) {
352  hist_drift_ = hists_.h_drift_depth_adc_[1]->getNbinsX();
353  hist_depth_ = hists_.h_drift_depth_adc_[1]->getNbinsY();
354  min_drift_ = hists_.h_drift_depth_adc_[1]->getAxisMin(1);
355  max_drift_ = hists_.h_drift_depth_adc_[1]->getAxisMax(1);
356  } else {
357  hist_drift_ = 100;
358  hist_depth_ = 50;
359  min_drift_ = -500.;
360  max_drift_ = 500.;
361  }
362 
363  iBooker.setCurrentFolder(fmt::format("{}Harvesting", dqmDir_));
364  MonitorElement* h_drift_depth_adc_slice_ =
365  iBooker.book1D("h_drift_depth_adc_slice", "slice of adc histogram", hist_drift_, min_drift_, max_drift_);
366 
367  // book histogram of differences
368  MonitorElement* h_diffLA = iBooker.book1D(
369  "h_diffLA", "difference in #mu_{H}; #Delta #mu_{H}/#mu_{H} (old-new)/old [%];n. modules", 100, -150, 150);
370 
371  // retrieve the number of bins from the other monitoring histogram
372  const auto& maxSect = hists_.h_bySectOccupancy_->getNbinsX();
373  const double lo = -0.5;
374  const double hi = maxSect - 0.5;
375 
376  // this will be booked in the Harvesting folder
377  iBooker.setCurrentFolder(fmt::format("{}Harvesting/SectorMonitoring", dqmDir_));
378  std::string repText = "%s tan#theta_{LA}/B by sector;pixel sector;%s tan(#theta_{LA})/B [1/T]";
380  iBooker.book1D("h_LAbySector_Measured", fmt::sprintf(repText, "measured", "measured"), maxSect, lo, hi);
382  iBooker.book1D("h_LAbySector_Accepted", fmt::sprintf(repText, "accepted", "accepted"), maxSect, lo, hi);
384  iBooker.book1D("h_LAbySector_Rejected", fmt::sprintf(repText, "rejected", "rejected"), maxSect, lo, hi);
385  hists_.h_bySectLA_ = iBooker.book1D("h_LAbySector", fmt::sprintf(repText, "payload", "payload"), maxSect, lo, hi);
387  iBooker.book1D("h_deltaLAbySector", fmt::sprintf(repText, "#Delta", "#Delta"), maxSect, lo, hi);
389  iBooker.book1D("h_bySectorChi2", "Fit #chi^{2}/ndf by sector;pixel sector; fit #chi^{2}/ndf", maxSect, lo, hi);
390 
392  iBooker.book1D("h_bySectFitStatus", "Fit Status by sector;pixel sector; fit status", maxSect, lo, hi);
393 
395  "h_bySectorCovMatrixStatus", "Fit Covariance Matrix Status by sector;pixel sector; fit status", maxSect, lo, hi);
397  iBooker.book1D("h_bySectorDriftError",
398  "square error on the measured drift at half-width by sector;pixel sector;#Delta d^{2}(t/2)",
399  maxSect,
400  lo,
401  hi);
402 
403  // set the bin labels for the fit quality dashboard histogram
404  std::vector<std::string> labels = {"#chi^{2}!=0", "#chi^{2} cut", "covStat>=2", "n. entries cut"};
405  hists_.h_bySectFitQuality_ = iBooker.book2D("h_bySectFitQuality",
406  "Fit quality by sector;pixel sector;quality cut",
407  maxSect,
408  lo,
409  hi,
410  labels.size(),
411  -0.5,
412  labels.size() - 0.5);
413 
414  // copy the bin labels from the occupancy histogram
415  for (int bin = 1; bin <= maxSect; bin++) {
416  const auto& binName = hists_.h_bySectOccupancy_->getTH1()->GetXaxis()->GetBinLabel(bin);
420  hists_.h_bySectLA_->setBinLabel(bin, binName);
422  hists_.h_bySectChi2_->setBinLabel(bin, binName);
427  }
428 
429  for (unsigned int binY = 1; binY <= labels.size(); binY++) {
430  hists_.h_bySectFitQuality_->setBinLabel(binY, labels[binY - 1], 2);
431  }
432 
433  // this will be booked in the Harvesting folder
434  iBooker.setCurrentFolder(fmt::format("{}Harvesting/LorentzAngleMaps", dqmDir_));
435  for (int i = 0; i < hists_.nlay; i++) {
436  std::string repName = "h2_byLayerLA_%i";
437  std::string repText = "BPix Layer %i tan#theta_{LA}/B;module number;ladder number;tan#theta_{LA}/B [1/T]";
438 
439  hists_.h2_byLayerLA_.emplace_back(iBooker.book2D(fmt::sprintf(repName, i + 1),
440  fmt::sprintf(repText, i + 1),
441  hists_.nModules_[i],
442  0.5,
443  hists_.nModules_[i] + 0.5,
444  hists_.nLadders_[i],
445  0.5,
446  hists_.nLadders_[i] + 0.5));
447 
448  repName = "h2_byLayerDiff_%i";
449  repText = "BPix Layer %i #Delta#mu_{H}/#mu_{H};module number;ladder number;#Delta#mu_{H}/#mu_{H} [%%]";
450  hists_.h2_byLayerDiff_.emplace_back(iBooker.book2D(fmt::sprintf(repName, i + 1),
451  fmt::sprintf(repText, i + 1),
452  hists_.nModules_[i],
453  0.5,
454  hists_.nModules_[i] + 0.5,
455  hists_.nLadders_[i],
456  0.5,
457  hists_.nLadders_[i] + 0.5));
458  }
459 
460  // clang-format off
461  edm::LogPrint("LorentzAngle") << "module" << "\t" << "layer" << "\t"
462  << "offset" << "\t" << "e0" << "\t"
463  << "slope" << "\t" << "e1" << "\t"
464  << "rel.err" << "\t" << "pull" << "\t"
465  << "p2" << "\t" << "e2" << "\t"
466  << "p3" << "\t" << "e3" << "\t"
467  << "p4" << "\t" << "e4" << "\t"
468  << "p5" << "\t" << "e5" << "\t"
469  << "chi2" << "\t" << "prob" << "\t"
470  << "newDetId" << "\t" << "tan(LA)" << "\t"
471  << "Error(LA)" ;
472  // clang-format on
473 
474  // payload to be written out
475  std::shared_ptr<SiPixelLorentzAngle> LorentzAngle = std::make_shared<SiPixelLorentzAngle>();
476 
477  // fill the map of simulation values
478  double p1_simul_newmodule = 0.294044;
479  double p1_simul[hists_.nlay + 1][hists_.nModules_[hists_.nlay - 1]];
480  for (int i_layer = 1; i_layer <= hists_.nlay; i_layer++) {
481  for (int i_module = 1; i_module <= hists_.nModules_[i_layer - 1]; i_module++) {
482  if (i_layer == 1)
483  p1_simul[i_layer - 1][i_module - 1] = 0.436848;
484  else if (i_layer == 2)
485  p1_simul[i_layer - 1][i_module - 1] = 0.25802;
486  else if (i_layer == 3 && i_module <= 4)
487  p1_simul[i_layer - 1][i_module - 1] = 0.29374;
488  else if (i_layer == 3 && i_module >= 5)
489  p1_simul[i_layer - 1][i_module - 1] = 0.31084;
490  else if (i_layer == 4 && i_module <= 4)
491  p1_simul[i_layer - 1][i_module - 1] = 0.29944;
492  else
493  p1_simul[i_layer - 1][i_module - 1] = 0.31426;
494  }
495  }
496  // fictitious n-th layer to store the values of new modules
497  for (int i_module = 1; i_module <= hists_.nModules_[hists_.nlay - 1]; i_module++) {
498  p1_simul[hists_.nlay][i_module - 1] = p1_simul_newmodule;
499  }
500 
501  // loop over "new" BPix modules
502  for (int j = 0; j < (int)hists_.BPixnewDetIds_.size(); j++) {
503  //uint32_t rawId = hists_.BPixnewDetIds_[j];
504  int new_index = j + 1 + hists_.nModules_[hists_.nlay - 1] + (hists_.nlay - 1) * hists_.nModules_[hists_.nlay - 1];
505  if (hists_.h_drift_depth_adc_[new_index] == nullptr)
506  continue;
507  for (int i = 1; i <= hist_depth_; i++) {
508  findMean(h_drift_depth_adc_slice_, i, new_index);
509  }
510 
511  // fit the distributions and store the LA in the payload
512  const auto& res = fitAndStore(LorentzAngle, new_index, hists_.BPixnewLayer_[j], hists_.BPixnewModule_[j]);
513 
514  edm::LogPrint("SiPixelLorentzAngle") << std::setprecision(4) << hists_.BPixnewModule_[j] << "\t"
515  << hists_.BPixnewLayer_[j] << "\t" << res.p0 << "\t" << res.e0 << "\t"
516  << res.p1 << std::setprecision(3) << "\t" << res.e1 << "\t"
517  << res.e1 / res.p1 * 100. << "\t"
518  << (res.p1 - p1_simul[hists_.nlay][0]) / res.e1 << "\t" << res.p2 << "\t"
519  << res.e2 << "\t" << res.p3 << "\t" << res.e3 << "\t" << res.p4 << "\t"
520  << res.e4 << "\t" << res.p5 << "\t" << res.e5 << "\t" << res.chi2 << "\t"
521  << res.prob << "\t" << hists_.BPixnewDetIds_[j] << "\t" << res.tan_LA << "\t"
522  << res.error_LA;
523  } // loop on BPix new modules
524 
525  //loop over modules and layers to fit the lorentz angle
526  for (int i_layer = 1; i_layer <= hists_.nlay; i_layer++) {
527  for (int i_module = 1; i_module <= hists_.nModules_[i_layer - 1]; i_module++) {
528  int i_index = i_module + (i_layer - 1) * hists_.nModules_[i_layer - 1];
529  if (hists_.h_drift_depth_adc_[i_index] == nullptr)
530  continue;
531  //loop over bins in depth (z-local-coordinate) (in order to fit slices)
532  for (int i = 1; i <= hist_depth_; i++) {
533  findMean(h_drift_depth_adc_slice_, i, i_index);
534  } // end loop over bins in depth
535 
536  // fit the distributions and store the LA in the payload
537  const auto& res = fitAndStore(LorentzAngle, i_index, i_layer, i_module);
538 
539  edm::LogPrint("SiPixelLorentzAngle")
540  << std::setprecision(4) << i_module << "\t" << i_layer << "\t" << res.p0 << "\t" << res.e0 << "\t" << res.p1
541  << std::setprecision(3) << "\t" << res.e1 << "\t" << res.e1 / res.p1 * 100. << "\t"
542  << (res.p1 - p1_simul[i_layer - 1][i_module - 1]) / res.e1 << "\t" << res.p2 << "\t" << res.e2 << "\t"
543  << res.p3 << "\t" << res.e3 << "\t" << res.p4 << "\t" << res.e4 << "\t" << res.p5 << "\t" << res.e5 << "\t"
544  << res.chi2 << "\t" << res.prob << "\t"
545  << "null"
546  << "\t" << res.tan_LA << "\t" << res.error_LA;
547  }
548  } // end loop over modules and layers
549 
550  // fill the rest of DetIds not filled above (for the moment FPix)
551  const auto& currentLAMap = currentLorentzAngle_->getLorentzAngles();
552  const auto& newLAMap = LorentzAngle->getLorentzAngles();
553  std::vector<unsigned int> currentLADets;
554  std::vector<unsigned int> newLADets;
555 
556  std::transform(currentLAMap.begin(),
557  currentLAMap.end(),
558  std::back_inserter(currentLADets),
559  [](const std::map<unsigned int, float>::value_type& pair) { return pair.first; });
560 
561  std::transform(newLAMap.begin(),
562  newLAMap.end(),
563  std::back_inserter(newLADets),
564  [](const std::map<unsigned int, float>::value_type& pair) { return pair.first; });
565 
566  std::vector<unsigned int> notCommon;
567  std::set_symmetric_difference(
568  currentLADets.begin(), currentLADets.end(), newLADets.begin(), newLADets.end(), std::back_inserter(notCommon));
569 
570  for (const auto& id : notCommon) {
571  float fPixLorentzAnglePerTesla_ = currentLorentzAngle_->getLorentzAngle(id);
572  if (!LorentzAngle->putLorentzAngle(id, fPixLorentzAnglePerTesla_)) {
573  edm::LogError("SiPixelLorentzAnglePCLHarvester")
574  << "[SiPixelLorentzAnglePCLHarvester::dqmEndJob] filling rest of payload: detid already exists";
575  }
576  }
577 
578  for (const auto& id : newLADets) {
579  float deltaMuHoverMuH = (currentLorentzAngle_->getLorentzAngle(id) - LorentzAngle->getLorentzAngle(id)) /
581  h_diffLA->Fill(deltaMuHoverMuH * 100.f);
582  }
583 
584  bool isPayloadChanged{false};
585  // fill the 2D output Lorentz Angle maps and check if the payload is different from the input one
586  for (const auto& [id, value] : LorentzAngle->getLorentzAngles()) {
587  DetId ID = DetId(id);
588  if (ID.subdetId() == PixelSubdetector::PixelBarrel) {
589  const auto& layer = theTrackerTopology_->pxbLayer(id);
590  const auto& ladder = theTrackerTopology_->pxbLadder(id);
591  const auto& module = theTrackerTopology_->pxbModule(id);
592  hists_.h2_byLayerLA_[layer - 1]->setBinContent(module, ladder, value);
593 
594  float deltaMuHoverMuH =
596  hists_.h2_byLayerDiff_[layer - 1]->setBinContent(module, ladder, deltaMuHoverMuH * 100.f);
597 
598  if (!isPayloadChanged && (deltaMuHoverMuH != 0.f))
599  isPayloadChanged = true;
600  }
601  }
602 
603  if (isPayloadChanged) {
604  // fill the DB object record
606  if (mydbservice.isAvailable()) {
607  try {
608  mydbservice->writeOneIOV(*LorentzAngle, mydbservice->currentTime(), recordName_);
609  } catch (const cond::Exception& er) {
610  edm::LogError("SiPixelLorentzAngleDB") << er.what();
611  } catch (const std::exception& er) {
612  edm::LogError("SiPixelLorentzAngleDB") << "caught std::exception " << er.what();
613  }
614  } else {
615  edm::LogError("SiPixelLorentzAngleDB") << "Service is unavailable";
616  }
617  } else {
618  edm::LogPrint("SiPixelLorentzAngleDB") << __PRETTY_FUNCTION__ << " there is no new valid measurement to append! ";
619  }
620 }
std::vector< dqm::reco::MonitorElement * > h2_byLayerLA_
Base exception class for the object to relational access.
Definition: Exception.h:11
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
uint32_t ID
Definition: Definitions.h:24
void findMean(MonitorElement *h_drift_depth_adc_slice_, int i, int i_ring)
SiPixelLAHarvest::fitResults fitAndStore(std::shared_ptr< SiPixelLorentzAngle > theLA, int i_idx, int i_lay, int i_mod)
const std::map< unsigned int, float > & getLorentzAngles() const
Log< level::Error, false > LogError
Definition: Electron.h:6
std::unique_ptr< TrackerTopology > theTrackerTopology_
constexpr std::array< uint8_t, layerIndexSize< TrackerTraits > > layer
void Fill(long long x)
Definition: EPCuts.h:4
Hash writeOneIOV(const T &payload, Time_t time, const std::string &recordName)
double f[11][100]
Definition: value.py:1
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)
Log< level::Warning, true > LogPrint
SiPixelLorentzAngleCalibrationHistograms hists_
Definition: DetId.h:17
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
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:697
virtual int getNbinsX() const
get # of bins in X-axis
bool isAvailable() const
Definition: Service.h:40
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
char const * what() const noexcept override
Definition: Exception.cc:103
float getLorentzAngle(const uint32_t &) const
std::vector< dqm::reco::MonitorElement * > h2_byLayerDiff_
unsigned transform(const HcalDetId &id, unsigned transformCode)

◆ endRun()

void SiPixelLorentzAnglePCLHarvester::endRun ( const edm::Run run,
const edm::EventSetup isetup 
)
overrideprivate

Definition at line 269 of file SiPixelLorentzAnglePCLHarvester.cc.

References edm::EventSetup::getData(), theTrackerTopology_, and topoEsTokenER_.

269  {
270  if (!theTrackerTopology_) {
271  theTrackerTopology_ = std::make_unique<TrackerTopology>(isetup.getData(topoEsTokenER_));
272  }
273 }
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoEsTokenER_
std::unique_ptr< TrackerTopology > theTrackerTopology_

◆ fillDescriptions()

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

Definition at line 882 of file SiPixelLorentzAnglePCLHarvester.cc.

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

882  {
884  desc.setComment("Harvester module of the SiPixel Lorentz Angle PCL monitoring workflow");
885  desc.add<std::vector<std::string>>("newmodulelist", {})->setComment("the list of DetIds for new sensors");
886  desc.add<std::string>("dqmDir", "AlCaReco/SiPixelLorentzAngle")->setComment("the directory of PCL Worker output");
887  desc.add<bool>("doChebyshevFit", false)->setComment("use Chebyshev polynomials for the dript vs depth fit");
888  desc.add<int>("order", 5)->setComment("order of the Chebyshev polynomial used for the fit");
889  desc.add<double>("fitChi2Cut", 20.)->setComment("cut on fit chi2/ndof to accept measurement");
890  desc.add<std::vector<double>>("fitRange", {5., 280.})->setComment("range of depths to perform the LA fit");
891  desc.add<int>("minHitsCut", 10000)->setComment("cut on minimum number of on-track hits to accept measurement");
892  desc.add<std::string>("record", "SiPixelLorentzAngleRcd")->setComment("target DB record");
893  descriptions.addWithDefaultLabel(desc);
894 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)

◆ findMean()

void SiPixelLorentzAnglePCLHarvester::findMean ( MonitorElement h_drift_depth_adc_slice_,
int  i,
int  i_ring 
)
private

Definition at line 623 of file SiPixelLorentzAnglePCLHarvester.cc.

References relativeConstraints::error, dqm::impl::MonitorElement::getMean(), dqm::impl::MonitorElement::getNbinsX(), dqm::impl::MonitorElement::getRMS(), SiPixelLorentzAngleCalibrationHistograms::h_drift_depth_, SiPixelLorentzAngleCalibrationHistograms::h_drift_depth_adc2_, SiPixelLorentzAngleCalibrationHistograms::h_drift_depth_adc_, SiPixelLorentzAngleCalibrationHistograms::h_drift_depth_noadc_, SiPixelLorentzAngleCalibrationHistograms::h_mean_, hists_, mps_fire::i, dqmiolumiharvest::j, SiStripPI::mean, dqm::impl::MonitorElement::Reset(), dqm::impl::MonitorElement::setBinContent(), dqm::impl::MonitorElement::setBinError(), and mathSSE::sqrt().

Referenced by dqmEndJob().

623  {
624  double nentries = 0;
625  h_drift_depth_adc_slice_->Reset();
626  int hist_drift_ = h_drift_depth_adc_slice_->getNbinsX();
627 
628  // determine sigma and sigma^2 of the adc counts and average adc counts
629  //loop over bins in drift width
630  for (int j = 1; j <= hist_drift_; j++) {
631  if (hists_.h_drift_depth_noadc_[i_ring]->getBinContent(j, i) >= 1) {
632  double adc_error2 = (hists_.h_drift_depth_adc2_[i_ring]->getBinContent(j, i) -
633  hists_.h_drift_depth_adc_[i_ring]->getBinContent(j, i) *
634  hists_.h_drift_depth_adc_[i_ring]->getBinContent(j, i) /
635  hists_.h_drift_depth_noadc_[i_ring]->getBinContent(j, i)) /
636  hists_.h_drift_depth_noadc_[i_ring]->getBinContent(j, i);
637 
638  hists_.h_drift_depth_adc_[i_ring]->setBinError(j, i, sqrt(adc_error2));
639  double error2 = adc_error2 / (hists_.h_drift_depth_noadc_[i_ring]->getBinContent(j, i) - 1.);
640  hists_.h_drift_depth_[i_ring]->setBinError(j, i, sqrt(error2));
641  } else {
642  hists_.h_drift_depth_[i_ring]->setBinError(j, i, 0);
643  hists_.h_drift_depth_adc_[i_ring]->setBinError(j, i, 0);
644  }
645  h_drift_depth_adc_slice_->setBinContent(j, hists_.h_drift_depth_adc_[i_ring]->getBinContent(j, i));
646  h_drift_depth_adc_slice_->setBinError(j, hists_.h_drift_depth_adc_[i_ring]->getBinError(j, i));
647  nentries += hists_.h_drift_depth_noadc_[i_ring]->getBinContent(j, i);
648  } // end loop over bins in drift width
649 
650  double mean = h_drift_depth_adc_slice_->getMean(1);
651  double error = 0;
652  if (nentries != 0) {
653  error = h_drift_depth_adc_slice_->getRMS(1) / std::sqrt(nentries);
654  }
655  hists_.h_mean_[i_ring]->setBinContent(i, mean);
656  hists_.h_mean_[i_ring]->setBinError(i, error);
657 
658  h_drift_depth_adc_slice_->Reset(); // clear again after extracting the parameters
659 }
virtual void Reset()
Remove all data from the ME, keept the empty histogram with all its settings.
virtual double getRMS(int axis=1) const
get RMS of histogram along x, y or z axis (axis=1, 2, 3 respectively)
T sqrt(T t)
Definition: SSEVec.h:19
SiPixelLorentzAngleCalibrationHistograms hists_
virtual void setBinContent(int binx, double content)
set content of bin (1-D)
virtual double getMean(int axis=1) const
get mean value of histogram along x, y or z axis (axis=1, 2, 3 respectively)
virtual int getNbinsX() const
get # of bins in X-axis
virtual void setBinError(int binx, double error)
set uncertainty on content of bin (1-D)

◆ fitAndStore()

SiPixelLAHarvest::fitResults SiPixelLorentzAnglePCLHarvester::fitAndStore ( std::shared_ptr< SiPixelLorentzAngle theLA,
int  i_idx,
int  i_lay,
int  i_mod 
)
private

Definition at line 662 of file SiPixelLorentzAnglePCLHarvester.cc.

References cms::cuda::assert(), SiPixelLorentzAngleCalibrationHistograms::BPixnewDetIds_, siPixelLACalibration::cmToum, currentLorentzAngle_, SiPixelLorentzAngleCalibrationHistograms::detIdsList, doChebyshevFit_, f1_, dqm::impl::MonitorElement::getBinContent(), SiPixelLorentzAngle::getLorentzAngle(), SiPixelLorentzAngleCalibrationHistograms::h_bySectChi2_, SiPixelLorentzAngleCalibrationHistograms::h_bySectCovMatrixStatus_, SiPixelLorentzAngleCalibrationHistograms::h_bySectDeltaLA_, SiPixelLorentzAngleCalibrationHistograms::h_bySectDriftError_, SiPixelLorentzAngleCalibrationHistograms::h_bySectFitQuality_, SiPixelLorentzAngleCalibrationHistograms::h_bySectFitStatus_, SiPixelLorentzAngleCalibrationHistograms::h_bySectLA_, SiPixelLorentzAngleCalibrationHistograms::h_bySectMeasLA_, SiPixelLorentzAngleCalibrationHistograms::h_bySectOccupancy_, SiPixelLorentzAngleCalibrationHistograms::h_bySectRejectLA_, SiPixelLorentzAngleCalibrationHistograms::h_bySectSetLA_, SiPixelLorentzAngleCalibrationHistograms::h_mean_, hists_, isFitGood(), dqmdumpme::k, SiPixelLAHarvest::kNotCalculated, SiPixelLAHarvest::kUndefined, cmsLHEtoEOSManager::l, LogDebug, SiPixelLorentzAngleCalibrationHistograms::nlay, SiPixelLorentzAngleCalibrationHistograms::nModules_, order_, funct::pow(), mps_fire::result, dqm::impl::MonitorElement::setBinContent(), dqm::impl::MonitorElement::setBinError(), mathSSE::sqrt(), theFitRange_, theMagField_, and width_.

Referenced by dqmEndJob().

663  {
664  // output results
666 
667  double half_width = width_ * siPixelLACalibration::cmToum / 2.f; // pixel half thickness in units of micro meter
668 
669  if (doChebyshevFit_) {
670  const int npar = order_ + 1;
671  auto cheb = std::make_unique<siPixelLACalibration::Chebyshev>(order_, theFitRange_.first, theFitRange_.second);
672  f1_ = std::make_unique<TF1>("f1", cheb.release(), theFitRange_.first, theFitRange_.second, npar, "Chebyshev");
673  } else {
674  f1_ = std::make_unique<TF1>("f1",
675  "[0] + [1]*x + [2]*x*x + [3]*x*x*x + [4]*x*x*x*x + [5]*x*x*x*x*x",
676  theFitRange_.first,
677  theFitRange_.second);
678  }
679 
680  // DO NOT REMOVE
681  // this is needed to ensure it stays synch with the render plugin:
682  // https://github.com/dmwm/deployment/blob/master/dqmgui/style/SiPixelLorentzAnglePCLRenderPlugin.cc#L199
683  // assert(std::string{f1_->GetName()}=="f1");
684  assert(strcmp(f1_->GetName(), "f1") == 0);
685 
686  f1_->SetParName(0, "offset");
687  f1_->SetParName(1, "tan#theta_{LA}");
688  f1_->SetParName(2, "quad term");
689  f1_->SetParName(3, "cubic term");
690  f1_->SetParName(4, "quartic term");
691  f1_->SetParName(5, "quintic term");
692 
693  f1_->SetParameter(0, 0);
694  f1_->SetParError(0, 0);
695  f1_->SetParameter(1, 0.4);
696  f1_->SetParError(1, 0);
697  f1_->SetParameter(2, 0.0);
698  f1_->SetParError(2, 0);
699  f1_->SetParameter(3, 0.0);
700  f1_->SetParError(3, 0);
701  f1_->SetParameter(4, 0.0);
702  f1_->SetParError(4, 0);
703  f1_->SetParameter(5, 0.0);
704  f1_->SetParError(5, 0);
705  f1_->SetChisquare(0);
706 
707  TFitResultPtr result = hists_.h_mean_[i_index]->getTH1()->Fit(f1_.get(), "ERQS");
708  if (result.Get()) {
709  res.fitStatus = result->Status();
710  edm::LogInfo("SiPixelLorentzAnglePCLHarvester") << "Fit status: " << res.fitStatus;
711  } else {
712  edm::LogError("SiPixelLorentzAnglePCLHarvester") << "Could not retrieve the fit status! Setting it to 0 by default";
713  res.fitStatus = -1;
714  }
715 
716  if (res.fitStatus >= 0) {
717  res.covMatrixStatus = result->CovMatrixStatus();
718  } else {
719  res.covMatrixStatus = SiPixelLAHarvest::kUndefined;
720  }
721 
722  // compute the error on the drift-at-half-width parameter (d^2)
723  float dSq{0.};
724  float cov00{0.}; // needed later for the error on the tan(theta_LA)
725  if (!doChebyshevFit_) {
726  if (res.covMatrixStatus > SiPixelLAHarvest::kNotCalculated) {
727  for (int k = 0; k < order_; k++) {
728  for (int l = 0; l < order_; l++) {
729  dSq += (std::pow(half_width, k) * std::pow(half_width, l) * result->CovMatrix(k, l));
730  }
731  }
732  cov00 = result->CovMatrix(0, 0);
733  } // if the covariance matrix is valid
734  } // compute the error on the drift-at-half width only for the regular polynomial fit
735 
736  res.dSq = dSq;
737 
738  res.p0 = f1_->GetParameter(0);
739  res.e0 = f1_->GetParError(0);
740  res.p1 = f1_->GetParameter(1);
741  res.e1 = f1_->GetParError(1);
742  res.p2 = f1_->GetParameter(2);
743  res.e2 = f1_->GetParError(2);
744  res.p3 = f1_->GetParameter(3);
745  res.e3 = f1_->GetParError(3);
746  res.p4 = f1_->GetParameter(4);
747  res.e4 = f1_->GetParError(4);
748  res.p5 = f1_->GetParameter(5);
749  res.e5 = f1_->GetParError(5);
750  res.chi2 = f1_->GetChisquare();
751  res.ndf = f1_->GetNDF();
752  res.prob = f1_->GetProb();
753  res.redChi2 = res.ndf > 0. ? res.chi2 / res.ndf : 0.;
754 
755  double f1_halfwidth = res.p0 + res.p1 * half_width + res.p2 * pow(half_width, 2) + res.p3 * pow(half_width, 3) +
756  res.p4 * pow(half_width, 4) + res.p5 * pow(half_width, 5);
757 
758  double f1_zerowidth = res.p0;
759 
760  // tan_LA = (f1(x = half_width) - f1(x = 0)) / (half_width - 0)
761  res.tan_LA = (f1_halfwidth - f1_zerowidth) / half_width;
762  double errsq_LA =
763  (pow(res.e1, 2) + pow((half_width * res.e2), 2) + pow((half_width * half_width * res.e3), 2) +
764  pow((half_width * half_width * half_width * res.e4), 2) +
765  pow((half_width * half_width * half_width * half_width * res.e5), 2)); // Propagation of uncertainty
766 
767  res.error_LA = doChebyshevFit_ ? sqrt(errsq_LA) : sqrt(res.dSq + cov00) / half_width;
768 
769  hists_.h_bySectMeasLA_->setBinContent(i_index, (res.tan_LA / theMagField_));
770  hists_.h_bySectMeasLA_->setBinError(i_index, (res.error_LA / theMagField_));
771  hists_.h_bySectChi2_->setBinContent(i_index, res.redChi2);
772  hists_.h_bySectChi2_->setBinError(i_index, 0.); // no errors
773 
774  hists_.h_bySectFitStatus_->setBinContent(i_index, res.fitStatus);
775  hists_.h_bySectFitStatus_->setBinError(i_index, 0.); // no errors
776  hists_.h_bySectCovMatrixStatus_->setBinContent(i_index, res.covMatrixStatus);
777  hists_.h_bySectCovMatrixStatus_->setBinError(i_index, 0.); // no errors
779  hists_.h_bySectDriftError_->setBinError(i_index, 0.);
780 
781  res.nentries = hists_.h_bySectOccupancy_->getBinContent(i_index); // number of on track hits in that sector
782 
783  bool isNew = (i_index > hists_.nlay * hists_.nModules_[hists_.nlay - 1]);
784  int shiftIdx = i_index - hists_.nlay * hists_.nModules_[hists_.nlay - 1] - 1;
785 
786  LogDebug("SiPixelLorentzAnglePCLHarvester")
787  << " isNew: " << isNew << " i_index: " << i_index << " shift index: " << shiftIdx;
788 
789  const auto& detIdsToFill =
790  isNew ? std::vector<unsigned int>({hists_.BPixnewDetIds_[shiftIdx]}) : hists_.detIdsList[i_index];
791 
792  LogDebug("SiPixelLorentzAnglePCLHarvester")
793  << "index: " << i_index << " i_module: " << i_module << " i_layer: " << i_layer;
794  for (const auto& id : detIdsToFill) {
795  LogDebug("SiPixelLorentzAnglePCLHarvester") << id << ",";
796  }
797 
798  // no errors on the following MEs
799  hists_.h_bySectSetLA_->setBinError(i_index, 0.);
800  hists_.h_bySectRejectLA_->setBinError(i_index, 0.);
801  hists_.h_bySectLA_->setBinError(i_index, 0.);
802  hists_.h_bySectDeltaLA_->setBinError(i_index, 0.);
803 
804  float LorentzAnglePerTesla_;
805  float currentLA = currentLorentzAngle_->getLorentzAngle(detIdsToFill.front());
806 
807  // fill the fit status dashboard
808  bool goodFit = isFitGood(res);
809  for (unsigned int i_bin = 0; i_bin < res.quality.size(); i_bin++) {
810  if (res.quality[i_bin]) {
811  hists_.h_bySectFitQuality_->setBinContent(i_index, i_bin + 1, 1.);
812  } else {
813  hists_.h_bySectFitQuality_->setBinContent(i_index, i_bin + 1, 0.);
814  }
815  }
816 
817  // if the fit quality is OK
818  if (goodFit) {
819  LorentzAnglePerTesla_ = res.tan_LA / theMagField_;
820  // fill the LA actually written to payload
821  hists_.h_bySectSetLA_->setBinContent(i_index, LorentzAnglePerTesla_);
822  hists_.h_bySectRejectLA_->setBinContent(i_index, 0.);
823  hists_.h_bySectLA_->setBinContent(i_index, LorentzAnglePerTesla_);
824 
825  const auto& deltaLA = (LorentzAnglePerTesla_ - currentLA);
826  hists_.h_bySectDeltaLA_->setBinContent(i_index, deltaLA);
827 
828  for (const auto& id : detIdsToFill) {
829  if (!theLAPayload->putLorentzAngle(id, LorentzAnglePerTesla_)) {
830  edm::LogError("SiPixelLorentzAnglePCLHarvester") << "[SiPixelLorentzAnglePCLHarvester::fitAndStore]: detid ("
831  << i_layer << "," << i_module << ") already exists";
832  }
833  }
834  } else {
835  // just copy the values from the existing payload
836  hists_.h_bySectSetLA_->setBinContent(i_index, 0.);
837  hists_.h_bySectRejectLA_->setBinContent(i_index, (res.tan_LA / theMagField_));
838  hists_.h_bySectLA_->setBinContent(i_index, currentLA);
839  hists_.h_bySectDeltaLA_->setBinContent(i_index, 0.);
840 
841  for (const auto& id : detIdsToFill) {
842  LorentzAnglePerTesla_ = currentLorentzAngle_->getLorentzAngle(id);
843  if (!theLAPayload->putLorentzAngle(id, LorentzAnglePerTesla_)) {
844  edm::LogError("SiPixelLorentzAnglePCLHarvester") << "[SiPixelLorentzAnglePCLHarvester::fitAndStore]: detid ("
845  << i_layer << "," << i_module << ") already exists";
846  }
847  }
848  }
849  // return the struct of fit details
850  return res;
851 }
Log< level::Error, false > LogError
assert(be >=bs)
Definition: Electron.h:6
T sqrt(T t)
Definition: SSEVec.h:19
std::unordered_map< uint32_t, std::vector< uint32_t > > detIdsList
SiPixelLorentzAngleCalibrationHistograms hists_
Log< level::Info, false > LogInfo
virtual void setBinContent(int binx, double content)
set content of bin (1-D)
virtual void setBinError(int binx, double error)
set uncertainty on content of bin (1-D)
float getLorentzAngle(const uint32_t &) const
bool isFitGood(SiPixelLAHarvest::fitResults &res)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
#define LogDebug(id)
virtual double getBinContent(int binx) const
get content of bin (1-D)

◆ isFitGood()

bool SiPixelLorentzAnglePCLHarvester::isFitGood ( SiPixelLAHarvest::fitResults res)
private

Definition at line 854 of file SiPixelLorentzAnglePCLHarvester.cc.

References PixelTripletNoTipGenerator_cfi::chi2Cut, doChebyshevFit_, fitChi2Cut_, SiPixelLAHarvest::kChi2Cut, SiPixelLAHarvest::kCovStatus, SiPixelLAHarvest::kMadePosDef, SiPixelLAHarvest::kNentries, SiPixelLAHarvest::kZeroChi2, and minHitsCut_.

Referenced by fitAndStore().

854  {
855  // check if reduced chi2 is different from 0.
856  const bool notZeroCut = (res.redChi2 != 0.);
857  // check the result of the chi2 only if doing the chebyshev fit
858  const bool chi2Cut = doChebyshevFit_ ? (res.redChi2 < fitChi2Cut_) : true;
859  // check the covariance matrix status
860  const bool covMatrixStatusCut = (res.covMatrixStatus >= SiPixelLAHarvest::kMadePosDef);
861  // check that the number of entries is larger than the minimum amount
862  const bool nentriesCut = (res.nentries > minHitsCut_);
863 
864  if (notZeroCut) {
865  res.setQualityCutBit(SiPixelLAHarvest::kZeroChi2);
866  }
867  if (chi2Cut) {
868  res.setQualityCutBit(SiPixelLAHarvest::kChi2Cut);
869  }
870  if (covMatrixStatusCut) {
871  res.setQualityCutBit(SiPixelLAHarvest::kCovStatus);
872  }
873  if (nentriesCut) {
874  res.setQualityCutBit(SiPixelLAHarvest::kNentries);
875  }
876 
877  // check if all the bits are set
878  return (res.quality.all());
879 }
Definition: Electron.h:6

Member Data Documentation

◆ currentLorentzAngle_

const SiPixelLorentzAngle* SiPixelLorentzAnglePCLHarvester::currentLorentzAngle_
private

Definition at line 141 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by beginRun(), dqmEndJob(), and fitAndStore().

◆ doChebyshevFit_

const bool SiPixelLorentzAnglePCLHarvester::doChebyshevFit_
private

Definition at line 127 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by fitAndStore(), and isFitGood().

◆ dqmDir_

const std::string SiPixelLorentzAnglePCLHarvester::dqmDir_
private

Definition at line 126 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by dqmEndJob().

◆ f1_

std::unique_ptr<TF1> SiPixelLorentzAnglePCLHarvester::f1_
private

Definition at line 133 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by fitAndStore().

◆ fitChi2Cut_

const double SiPixelLorentzAnglePCLHarvester::fitChi2Cut_
private

Definition at line 129 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by isFitGood().

◆ fitRange_

const std::vector<double> SiPixelLorentzAnglePCLHarvester::fitRange_
private

Definition at line 130 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by SiPixelLorentzAnglePCLHarvester().

◆ geomEsToken_

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

Definition at line 120 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by beginRun().

◆ hists_

SiPixelLorentzAngleCalibrationHistograms SiPixelLorentzAnglePCLHarvester::hists_
private

Definition at line 140 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by beginRun(), dqmEndJob(), findMean(), and fitAndStore().

◆ magneticFieldToken_

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

Definition at line 123 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by beginRun().

◆ minHitsCut_

const int SiPixelLorentzAnglePCLHarvester::minHitsCut_
private

Definition at line 131 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by isFitGood().

◆ newmodulelist_

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

Definition at line 125 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by beginRun().

◆ order_

const int SiPixelLorentzAnglePCLHarvester::order_
private

Definition at line 128 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by fitAndStore().

◆ recordName_

const std::string SiPixelLorentzAnglePCLHarvester::recordName_
private

Definition at line 132 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by dqmEndJob().

◆ siPixelLAEsToken_

edm::ESGetToken<SiPixelLorentzAngle, SiPixelLorentzAngleRcd> SiPixelLorentzAnglePCLHarvester::siPixelLAEsToken_
private

Definition at line 122 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by beginRun().

◆ teslaToInverseGeV_

constexpr float SiPixelLorentzAnglePCLHarvester::teslaToInverseGeV_ = 2.99792458e-3f
staticprivate

Definition at line 137 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by beginRun().

◆ theFitRange_

std::pair<double, double> SiPixelLorentzAnglePCLHarvester::theFitRange_ {5., 280.}
private

◆ theMagField_

float SiPixelLorentzAnglePCLHarvester::theMagField_ {0.f}
private

Definition at line 135 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by beginRun(), and fitAndStore().

◆ theTrackerTopology_

std::unique_ptr<TrackerTopology> SiPixelLorentzAnglePCLHarvester::theTrackerTopology_
private

Definition at line 142 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by dqmEndJob(), and endRun().

◆ topoEsTokenBR_

edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> SiPixelLorentzAnglePCLHarvester::topoEsTokenBR_
private

Definition at line 121 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by beginRun().

◆ topoEsTokenER_

edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> SiPixelLorentzAnglePCLHarvester::topoEsTokenER_
private

Definition at line 121 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by endRun().

◆ width_

float SiPixelLorentzAnglePCLHarvester::width_
private

Definition at line 134 of file SiPixelLorentzAnglePCLHarvester.cc.

Referenced by beginRun(), and fitAndStore().