CMS 3D CMS Logo

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  , pu_binning_ ( getHistoPSet (iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("puPSet")) )
17  , ls_binning_ ( getHistoLSPSet(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("lsPSet")) )
18  , doPixelLumi_ ( iConfig.getParameter<bool>("doPixelLumi") )
19  , pixelClustersToken_ ( doPixelLumi_ ? consumes<edmNew::DetSetVector<SiPixelCluster> >(iConfig.getParameter<edm::InputTag>("pixelClusters") ) : edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster>> () )
20  , useBPixLayer1_ ( doPixelLumi_ ? iConfig.getParameter<bool> ("useBPixLayer1") : false )
21  , minNumberOfPixelsPerCluster_ ( doPixelLumi_ ? iConfig.getParameter<int> ("minNumberOfPixelsPerCluster") : -1 )
22  , minPixelClusterCharge_ ( doPixelLumi_ ? iConfig.getParameter<double> ("minPixelClusterCharge") : -1. )
23  , pixelCluster_binning_ ( doPixelLumi_ ? getHistoPSet (iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("pixelClusterPSet")) : MEbinning {} )
24  , pixellumi_binning_ ( doPixelLumi_ ? getHistoPSet (iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("pixellumiPSet")) : MEbinning {} )
25 {
26 
29  lumiVsLS_ = nullptr;
30  puVsLS_ = nullptr;
31  pixelLumiVsLS_ = nullptr;
32  pixelLumiVsLumi_ = nullptr;
33 
34  if(useBPixLayer1_)
36  else
38 
39 }
40 
42 {
43  return MEbinning{
44  pset.getParameter<int32_t>("nbins"),
45  pset.getParameter<double>("xmin"),
46  pset.getParameter<double>("xmax"),
47  };
48 }
49 
51 {
52  return MEbinning{
53  pset.getParameter<int32_t>("nbins"),
54  -0.5,
55  double(pset.getParameter<int32_t>("nbins")-0.5)
56  };
57 }
58 
60  edm::Run const & iRun,
61  edm::EventSetup const & iSetup)
62 {
63 
64  std::string histname, histtitle;
65 
66  std::string currentFolder = folderName_ ;
67  ibooker.setCurrentFolder(currentFolder);
68 
69  if ( doPixelLumi_ ) {
70  histname = "numberOfPixelClustersVsLS"; histtitle = "number of pixel clusters vs LS";
71  numberOfPixelClustersVsLS_ = ibooker.book1D(histname, histtitle,
73  // numberOfPixelClustersVsLS_->getTH1()->SetCanExtend(TH1::kAllAxes);
75  numberOfPixelClustersVsLS_->setAxisTitle("number of pixel clusters",2);
76 
77  histname = "numberOfPixelClustersVsLumi"; histtitle = "number of pixel clusters vs scal lumi";
78  numberOfPixelClustersVsLumi_ = ibooker.bookProfile(histname, histtitle,
81  numberOfPixelClustersVsLumi_->setAxisTitle("scal inst lumi E30 [Hz cm^{-2}]",1);
82  numberOfPixelClustersVsLumi_->setAxisTitle("number of pixel clusters",2);
83 
84  histname = "pixelLumiVsLS"; histtitle = "pixel-lumi vs LS";
85  pixelLumiVsLS_ = ibooker.bookProfile(histname, histtitle,
88  // pixelLumiVsLS_->getTH1()->SetCanExtend(TH1::kAllAxes);
90  pixelLumiVsLS_->setAxisTitle("pixel-based inst lumi E30 [Hz cm^{-2}]",2);
91 
92  histname = "pixelLumiVsLumi"; histtitle = "pixel-lumi vs scal lumi";
93  pixelLumiVsLumi_ = ibooker.bookProfile(histname, histtitle,
96  pixelLumiVsLumi_->setAxisTitle("scal inst lumi E30 [Hz cm^{-2}]",1);
97  pixelLumiVsLumi_->setAxisTitle("pixel-based inst lumi E30 [Hz cm^{-2}]",2);
98  }
99 
100  histname = "lumiVsLS"; histtitle = "scal lumi vs LS";
101  lumiVsLS_ = ibooker.bookProfile(histname, histtitle,
104  // lumiVsLS_->getTH1()->SetCanExtend(TH1::kAllAxes);
105  lumiVsLS_->setAxisTitle("LS",1);
106  lumiVsLS_->setAxisTitle("scal inst lumi E30 [Hz cm^{-2}]",2);
107 
108  histname = "puVsLS"; histtitle = "scal PU vs LS";
109  puVsLS_ = ibooker.bookProfile(histname, histtitle,
112  // puVsLS_->getTH1()->SetCanExtend(TH1::kAllAxes);
113  puVsLS_->setAxisTitle("LS",1);
114  puVsLS_->setAxisTitle("scal PU",2);
115 
116 }
117 
123 
124  // int bx = iEvent.bunchCrossing();
125  int ls = iEvent.id().luminosityBlock();
126 
127  float scal_lumi = -1.;
128  float scal_pu = -1.;
130  iEvent.getByToken(lumiScalersToken_, lumiScalers);
131  if ( lumiScalers.isValid() && !lumiScalers->empty() ) {
132  LumiScalersCollection::const_iterator scalit = lumiScalers->begin();
133  scal_lumi = scalit->instantLumi();
134  scal_pu = scalit->pileup();
135  } else {
136  scal_lumi = -1.;
137  scal_pu = -1.;
138  }
139  lumiVsLS_ -> Fill(ls, scal_lumi);
140  puVsLS_ -> Fill(ls, scal_pu );
141 
142  if ( doPixelLumi_ ) {
143  size_t pixel_clusters = 0;
144  float pixel_lumi = -1.;
146  iEvent.getByToken(pixelClustersToken_, pixelClusters);
147  if ( pixelClusters.isValid() ) {
148 
149  edm::ESHandle<TrackerTopology> tTopoHandle;
150  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
151  const TrackerTopology* const tTopo = tTopoHandle.product();
152 
153  // Count the number of clusters with at least a minimum
154  // number of pixels per cluster and at least a minimum charge.
155  size_t tot = 0;
156  edmNew::DetSetVector<SiPixelCluster>::const_iterator pixCluDet = pixelClusters->begin();
157  for ( ; pixCluDet!=pixelClusters->end(); ++pixCluDet) {
158 
159  DetId detid = pixCluDet->detId();
160  size_t subdetid = detid.subdetId();
161  if ( subdetid == (int) PixelSubdetector::PixelBarrel )
162  if ( tTopo->layer(detid)==1 )
163  continue;
164 
166  for ( ; pixClu != pixCluDet->end(); ++pixClu ) {
167  ++tot;
168  if ( (pixClu->size() >= minNumberOfPixelsPerCluster_) &&
169  (pixClu->charge() >= minPixelClusterCharge_ ) ) {
170  ++pixel_clusters;
171  }
172  }
173  }
174  pixel_lumi = lumi_factor_per_bx_ * pixel_clusters / GetLumi::CM2_TO_NANOBARN ; // ?!?!
175  } else
176  pixel_lumi = -1.;
177 
178  numberOfPixelClustersVsLS_ -> Fill(ls, pixel_clusters);
179  numberOfPixelClustersVsLumi_ -> Fill(scal_lumi,pixel_clusters);
180  pixelLumiVsLS_ -> Fill(ls, pixel_lumi);
181  pixelLumiVsLumi_ -> Fill(scal_lumi,pixel_lumi);
182  }
183 
184 }
185 
187 {
188  pset.add<int> ( "nbins");
189  pset.add<double>( "xmin" );
190  pset.add<double>( "xmax" );
191 }
192 
194 {
195  pset.add<int> ( "nbins", 2500);
196 }
197 
199 {
201  desc.add<edm::InputTag>( "pixelClusters", edm::InputTag("hltSiPixelClusters") );
202  desc.add<edm::InputTag>( "scalers", edm::InputTag("hltScalersRawToDigi"));
203  desc.add<std::string> ( "FolderName", "HLT/LumiMonitoring" );
204  desc.add<bool> ( "doPixelLumi", false );
205  desc.add<bool> ( "useBPixLayer1", false );
206  desc.add<int> ( "minNumberOfPixelsPerCluster", 2 ); // from DQM/PixelLumi/python/PixelLumiDQM_cfi.py
207  desc.add<double> ( "minPixelClusterCharge", 15000. );
208 
210  edm::ParameterSetDescription pixelClusterPSet;
211  LumiMonitor::fillHistoPSetDescription(pixelClusterPSet);
212  histoPSet.add("pixelClusterPSet", pixelClusterPSet);
213 
215  fillHistoPSetDescription(lumiPSet);
216  histoPSet.add<edm::ParameterSetDescription>("lumiPSet", lumiPSet);
217 
219  fillHistoPSetDescription(puPSet);
220  histoPSet.add<edm::ParameterSetDescription>("puPSet", puPSet);
221 
222  edm::ParameterSetDescription pixellumiPSet;
223  fillHistoPSetDescription(pixellumiPSet);
224  histoPSet.add<edm::ParameterSetDescription>("pixellumiPSet", pixellumiPSet);
225 
228  histoPSet.add<edm::ParameterSetDescription>("lsPSet", lsPSet);
229 
230  desc.add<edm::ParameterSetDescription>("histoPSet",histoPSet);
231 
232  descriptions.add("lumiMonitor", desc);
233 }
234 
235 // Define this as a plug-in
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: LumiMonitor.cc:59
T getParameter(std::string const &) const
bool useBPixLayer1_
Definition: LumiMonitor.h:68
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:508
MEbinning pu_binning_
Definition: LumiMonitor.h:63
float minPixelClusterCharge_
Definition: LumiMonitor.h:70
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
data_type const * const_iterator
Definition: DetSetNew.h:30
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelClustersToken_
Definition: LumiMonitor.h:67
LuminosityBlockNumber_t luminosityBlock() const
Definition: EventID.h:40
static double XSEC_PIXEL_CLUSTER
Definition: GetLumi.h:41
MonitorElement * numberOfPixelClustersVsLumi_
Definition: LumiMonitor.h:78
int iEvent
Definition: GenABIO.cc:230
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: LumiMonitor.cc:198
std::string folderName_
Definition: LumiMonitor.h:59
MonitorElement * pixelLumiVsLumi_
Definition: LumiMonitor.h:82
edm::EDGetTokenT< LumiScalersCollection > lumiScalersToken_
Definition: LumiMonitor.h:61
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
MEbinning pixellumi_binning_
Definition: LumiMonitor.h:72
double xmin
Definition: LumiMonitor.h:30
void analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup) override
Definition: LumiMonitor.cc:122
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
static void fillHistoLSPSetDescription(edm::ParameterSetDescription &pset)
Definition: LumiMonitor.cc:193
static double CM2_TO_NANOBARN
Definition: GetLumi.h:47
static MEbinning getHistoPSet(edm::ParameterSet pset)
Definition: LumiMonitor.cc:41
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
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:79
MonitorElement * pixelLumiVsLS_
Definition: LumiMonitor.h:81
static void fillHistoPSetDescription(edm::ParameterSetDescription &pset)
Definition: LumiMonitor.cc:186
MonitorElement * puVsLS_
Definition: LumiMonitor.h:80
Definition: DetId.h:18
def ls(path, rec=False)
Definition: eostools.py:348
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
const T & get() const
Definition: EventSetup.h:55
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static double rXSEC_PIXEL_CLUSTER
Definition: GetLumi.h:45
bool doPixelLumi_
Definition: LumiMonitor.h:66
MEbinning pixelCluster_binning_
Definition: LumiMonitor.h:71
unsigned int layer(const DetId &id) const
static MEbinning getHistoLSPSet(edm::ParameterSet pset)
Definition: LumiMonitor.cc:50
edm::EventID id() const
Definition: EventBase.h:60
MonitorElement * numberOfPixelClustersVsLS_
Definition: LumiMonitor.h:77
Pixel cluster – collection of neighboring pixels above threshold.
HLT enums.
float lumi_factor_per_bx_
Definition: LumiMonitor.h:84
int minNumberOfPixelsPerCluster_
Definition: LumiMonitor.h:69
std::vector< LumiScalers > LumiScalersCollection
Definition: LumiScalers.h:160
MEbinning lumi_binning_
Definition: LumiMonitor.h:62
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)
T const * product() const
Definition: ESHandle.h:86
static double FREQ_ORBIT
Definition: GetLumi.h:32
const_iterator begin(bool update=false) const
MEbinning ls_binning_
Definition: LumiMonitor.h:64
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