CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LumiMonitor.cc
Go to the documentation of this file.
2 
4 
6 
7 
8 // -----------------------------
9 // constructors and destructor
10 // -----------------------------
11 
13  folderName_ ( iConfig.getParameter<std::string>("FolderName") )
14  , lumiScalersToken_ ( consumes<LumiScalersCollection>(iConfig.getParameter<edm::InputTag>("scalers") ) )
15  , lumi_binning_ ( getHistoPSet (iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("lumiPSet")) )
16  , ls_binning_ ( getHistoLSPSet(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("lsPSet")) )
17  , doPixelLumi_ ( iConfig.getParameter<bool>("doPixelLumi") )
18  , pixelClustersToken_ ( doPixelLumi_ ? consumes<edmNew::DetSetVector<SiPixelCluster> >(iConfig.getParameter<edm::InputTag>("pixelClusters") ) : edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster>> () )
19  , useBPixLayer1_ ( doPixelLumi_ ? iConfig.getParameter<bool> ("useBPixLayer1") : false )
20  , minNumberOfPixelsPerCluster_ ( doPixelLumi_ ? iConfig.getParameter<int> ("minNumberOfPixelsPerCluster") : -1 )
21  , minPixelClusterCharge_ ( doPixelLumi_ ? iConfig.getParameter<double> ("minPixelClusterCharge") : -1. )
22  , pixelCluster_binning_ ( doPixelLumi_ ? getHistoPSet (iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("pixelClusterPSet")) : MEbinning {} )
23  , pixellumi_binning_ ( doPixelLumi_ ? getHistoPSet (iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("pixellumiPSet")) : MEbinning {} )
24 {
25 
26  numberOfPixelClustersVsLS_ = nullptr;
27  numberOfPixelClustersVsLumi_ = nullptr;
28  lumiVsLS_ = nullptr;
29  pixelLumiVsLS_ = nullptr;
30  pixelLumiVsLumi_ = nullptr;
31 
32  if(useBPixLayer1_)
34  else
36 
37 }
38 
40 {
41  return MEbinning{
42  pset.getParameter<int32_t>("nbins"),
43  pset.getParameter<double>("xmin"),
44  pset.getParameter<double>("xmax"),
45  };
46 }
47 
49 {
50  return MEbinning{
51  pset.getParameter<int32_t>("nbins"),
52  0.,
53  double(pset.getParameter<int32_t>("nbins"))
54  };
55 }
56 
58  edm::Run const & iRun,
59  edm::EventSetup const & iSetup)
60 {
61 
62  std::string histname, histtitle;
63 
64  std::string currentFolder = folderName_ ;
65  ibooker.setCurrentFolder(currentFolder.c_str());
66 
67  if ( doPixelLumi_ ) {
68  histname = "numberOfPixelClustersVsLS"; histtitle = "number of pixel clusters vs LS";
69  numberOfPixelClustersVsLS_ = ibooker.book1D(histname, histtitle,
71  // numberOfPixelClustersVsLS_->getTH1()->SetCanExtend(TH1::kAllAxes);
73  numberOfPixelClustersVsLS_->setAxisTitle("number of pixel clusters",2);
74 
75  histname = "numberOfPixelClustersVsLumi"; histtitle = "number of pixel clusters vs scal lumi";
76  numberOfPixelClustersVsLumi_ = ibooker.bookProfile(histname, histtitle,
79  numberOfPixelClustersVsLumi_->setAxisTitle("scal inst lumi E30 [Hz cm^{-2}]",1);
80  numberOfPixelClustersVsLumi_->setAxisTitle("number of pixel clusters",2);
81 
82  histname = "pixelLumiVsLS"; histtitle = "pixel-lumi vs LS";
83  pixelLumiVsLS_ = ibooker.bookProfile(histname, histtitle,
86  // pixelLumiVsLS_->getTH1()->SetCanExtend(TH1::kAllAxes);
88  pixelLumiVsLS_->setAxisTitle("pixel-based inst lumi E30 [Hz cm^{-2}]",2);
89 
90  histname = "pixelLumiVsLumi"; histtitle = "pixel-lumi vs scal lumi";
91  pixelLumiVsLumi_ = ibooker.bookProfile(histname, histtitle,
94  pixelLumiVsLumi_->setAxisTitle("scal inst lumi E30 [Hz cm^{-2}]",1);
95  pixelLumiVsLumi_->setAxisTitle("pixel-based inst lumi E30 [Hz cm^{-2}]",2);
96  }
97 
98  histname = "lumiVsLS"; histtitle = "scal lumi vs LS";
99  lumiVsLS_ = ibooker.bookProfile(histname, histtitle,
102  // lumiVsLS_->getTH1()->SetCanExtend(TH1::kAllAxes);
103  lumiVsLS_->setAxisTitle("LS",1);
104  lumiVsLS_->setAxisTitle("scal inst lumi E30 [Hz cm^{-2}]",2);
105 
106 }
107 
113 
114  // int bx = iEvent.bunchCrossing();
115  int ls = iEvent.id().luminosityBlock();
116 
117  float scal_lumi = -1.;
119  iEvent.getByToken(lumiScalersToken_, lumiScalers);
120  if ( lumiScalers.isValid() && lumiScalers->size() ) {
121  LumiScalersCollection::const_iterator scalit = lumiScalers->begin();
122  scal_lumi = scalit->instantLumi();
123  } else {
124  scal_lumi = -1.;
125  }
126  lumiVsLS_ -> Fill(ls, scal_lumi);
127 
128  if ( doPixelLumi_ ) {
129  size_t pixel_clusters = 0;
130  float pixel_lumi = -1.;
132  iEvent.getByToken(pixelClustersToken_, pixelClusters);
133  if ( pixelClusters.isValid() ) {
134 
135  edm::ESHandle<TrackerTopology> tTopoHandle;
136  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
137  const TrackerTopology* const tTopo = tTopoHandle.product();
138 
139  // Count the number of clusters with at least a minimum
140  // number of pixels per cluster and at least a minimum charge.
141  size_t tot = 0;
142  edmNew::DetSetVector<SiPixelCluster>::const_iterator pixCluDet = pixelClusters->begin();
143  for ( ; pixCluDet!=pixelClusters->end(); ++pixCluDet) {
144 
145  DetId detid = pixCluDet->detId();
146  size_t subdetid = detid.subdetId();
147  if ( subdetid == (int) PixelSubdetector::PixelBarrel )
148  if ( tTopo->layer(detid)==1 )
149  continue;
150 
152  for ( ; pixClu != pixCluDet->end(); ++pixClu ) {
153  ++tot;
154  if ( (pixClu->size() >= minNumberOfPixelsPerCluster_) &&
155  (pixClu->charge() >= minPixelClusterCharge_ ) ) {
156  ++pixel_clusters;
157  }
158  }
159  }
160  pixel_lumi = lumi_factor_per_bx_ * pixel_clusters / GetLumi::CM2_TO_NANOBARN ; // ?!?!
161  } else
162  pixel_lumi = -1.;
163 
164  numberOfPixelClustersVsLS_ -> Fill(ls, pixel_clusters);
165  numberOfPixelClustersVsLumi_ -> Fill(scal_lumi,pixel_clusters);
166  pixelLumiVsLS_ -> Fill(ls, pixel_lumi);
167  pixelLumiVsLumi_ -> Fill(scal_lumi,pixel_lumi);
168  }
169 
170 }
171 
173 {
174  pset.add<int> ( "nbins");
175  pset.add<double>( "xmin" );
176  pset.add<double>( "xmax" );
177 }
178 
180 {
181  pset.add<int> ( "nbins", 2500);
182 }
183 
185 {
187  desc.add<edm::InputTag>( "pixelClusters", edm::InputTag("hltSiPixelClusters") );
188  desc.add<edm::InputTag>( "scalers", edm::InputTag("hltScalersRawToDigi"));
189  desc.add<std::string> ( "FolderName", "HLT/LumiMonitoring" );
190  desc.add<bool> ( "doPixelLumi", false );
191  desc.add<bool> ( "useBPixLayer1", false );
192  desc.add<int> ( "minNumberOfPixelsPerCluster", 2 ); // from DQM/PixelLumi/python/PixelLumiDQM_cfi.py
193  desc.add<double> ( "minPixelClusterCharge", 15000. );
194 
196  edm::ParameterSetDescription pixelClusterPSet;
197  LumiMonitor::fillHistoPSetDescription(pixelClusterPSet);
198  histoPSet.add("pixelClusterPSet", pixelClusterPSet);
199 
201  fillHistoPSetDescription(lumiPSet);
202  histoPSet.add<edm::ParameterSetDescription>("lumiPSet", lumiPSet);
203 
204  edm::ParameterSetDescription pixellumiPSet;
205  fillHistoPSetDescription(pixellumiPSet);
206  histoPSet.add<edm::ParameterSetDescription>("pixellumiPSet", pixellumiPSet);
207 
210  histoPSet.add<edm::ParameterSetDescription>("lsPSet", lsPSet);
211 
212  desc.add<edm::ParameterSetDescription>("histoPSet",histoPSet);
213 
214  descriptions.add("lumiMonitor", desc);
215 }
216 
217 // Define this as a plug-in
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: LumiMonitor.cc:57
T getParameter(std::string const &) const
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const_iterator end(bool update=false) const
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:157
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
float minPixelClusterCharge_
Definition: LumiMonitor.h:69
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
def ls
Definition: eostools.py:348
data_type const * const_iterator
Definition: DetSetNew.h:30
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelClustersToken_
Definition: LumiMonitor.h:66
LuminosityBlockNumber_t luminosityBlock() const
Definition: EventID.h:40
static double XSEC_PIXEL_CLUSTER
Definition: GetLumi.h:41
MonitorElement * numberOfPixelClustersVsLumi_
Definition: LumiMonitor.h:77
int iEvent
Definition: GenABIO.cc:230
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: LumiMonitor.cc:184
std::string folderName_
Definition: LumiMonitor.h:59
MonitorElement * pixelLumiVsLumi_
Definition: LumiMonitor.h:80
edm::EDGetTokenT< LumiScalersCollection > lumiScalersToken_
Definition: LumiMonitor.h:61
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
MEbinning pixellumi_binning_
Definition: LumiMonitor.h:71
double xmin
Definition: LumiMonitor.h:30
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
static void fillHistoLSPSetDescription(edm::ParameterSetDescription &pset)
Definition: LumiMonitor.cc:179
static double CM2_TO_NANOBARN
Definition: GetLumi.h:47
static MEbinning getHistoPSet(edm::ParameterSet pset)
Definition: LumiMonitor.cc:39
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:75
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
MonitorElement * lumiVsLS_
Definition: LumiMonitor.h:78
MonitorElement * pixelLumiVsLS_
Definition: LumiMonitor.h:79
static void fillHistoPSetDescription(edm::ParameterSetDescription &pset)
Definition: LumiMonitor.cc:172
Definition: DetId.h:18
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:276
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static double rXSEC_PIXEL_CLUSTER
Definition: GetLumi.h:45
bool doPixelLumi_
Definition: LumiMonitor.h:65
MEbinning pixelCluster_binning_
Definition: LumiMonitor.h:70
unsigned int layer(const DetId &id) const
static MEbinning getHistoLSPSet(edm::ParameterSet pset)
Definition: LumiMonitor.cc:48
edm::EventID id() const
Definition: EventBase.h:59
MonitorElement * numberOfPixelClustersVsLS_
Definition: LumiMonitor.h:76
Pixel cluster – collection of neighboring pixels above threshold.
float lumi_factor_per_bx_
Definition: LumiMonitor.h:82
int minNumberOfPixelsPerCluster_
Definition: LumiMonitor.h:68
std::vector< LumiScalers > LumiScalersCollection
Definition: LumiScalers.h:160
void analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup)
Definition: LumiMonitor.cc:112
MEbinning lumi_binning_
Definition: LumiMonitor.h:62
volatile std::atomic< bool > shutdown_flag false
size_type size() const
Definition: DetSetNew.h:87
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
static double FREQ_ORBIT
Definition: GetLumi.h:32
const_iterator begin(bool update=false) const
MEbinning ls_binning_
Definition: LumiMonitor.h:63
Definition: Run.h:43
LumiMonitor(const edm::ParameterSet &)
Definition: LumiMonitor.cc:12
static double SECONDS_PER_LS
Definition: GetLumi.h:33
double xmax
Definition: LumiMonitor.h:31