CMS 3D CMS Logo

SiPixelOfflineCalibAnalysisBase.h
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelOfflineCalibAnalysisBase
4 // Class: SiPixelOfflineCalibAnalysisBase
5 //
13 //
14 // Original Author: Evan Klose Friis
15 // additions by: Freya Blekman
16 // Created: Tue Nov 6 17:27:19 CET 2007
17 //
18 //
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
26 
29 //#include "FWCore/Framework/interface/EventSetup.h"
32 
35 
38 
41 
45 
48 
49 #include "TF1.h"
50 
57 
58 #include <map>
59 #include <vector>
60 #include <iostream>
61 #include <string>
62 //
63 // class decleration
64 //
65 
66 
68 public:
71 
72  //no argument constructor only used to throw exception in the case of derived
73  //class constructor not calling SiPixelOfflineCalibAnalysisBase(const edm::ParameterSet&)
75 
76  //pure virtual function, called after each DetID loaded
77  virtual bool doFits(uint32_t detid, std::vector<SiPixelCalibDigi>::const_iterator ipix);
78 
79 
80  //translate DetID to human readable string
81  std::string translateDetIdToString(uint32_t detid);
82 
83  //booking DQM histograms (for dynamic histogram creation)
84 
85  MonitorElement* bookDQMHistogram1D(uint32_t detid, std::string name, std::string title, int nchX, double lowX, double highX);
86  MonitorElement* bookDQMHistogram1D(uint32_t detid, std::string name, std::string title, int nchX, float *xbinsize); //variable size bins
87  MonitorElement* bookDQMHistogram2D(uint32_t detid, std::string name, std::string title, int nchX, double lowX, double highX, int nchY, double lowY, double highY);
88 
89  MonitorElement* bookDQMHistoPlaquetteSummary2D(uint32_t detid, std::string name,std::string title); // take the detid to determine the size of rows and columns, this saves looking up everything in the cabling map by the user.
91 
93  bool setDQMDirectory(uint32_t detID); //automatically create directory hierachy based on DetID
94  static TF1* fitFunction_;
95  static const std::vector<short>* getVcalValues();
96  std::vector<uint32_t> & getRunNumbers() { return runnumbers_;}
97 
98 protected:
99 
100  //calibration parameters from calib.dat/DB
104 
106  short nTriggers_;
107  static std::vector<short> vCalValues_;
108  uint32_t & EventNumber() { return eventCounter_;}
109 
110 
111 private:
112 
118 
119  std::vector<uint32_t> runnumbers_;
120  uint32_t eventCounter_;
121 
122  //store set of detIDs that have been encountered
123  //second int argument can be a derived class result flag
124  std::map<uint32_t, int> detIdsEntered_;
125  std::map<uint32_t, std::string> detIdNames_;
126 
128 
130 
131  //the beginJob is used to load the calib database. It then calls the pure
132  //virtual calibrationSetup() function. Derived classes should put beginJob functionality there.
133  virtual void beginRun(const edm::Run &, const edm::EventSetup &);
134  void beginRun(const edm::EventSetup& iSetup);
135  void beginJob();
136 
137  //calibrationSetup will be used by derived classes
138  virtual void calibrationSetup(const edm::EventSetup& iSetup);
139 
140  // pure virtual function, checks if the calibration analyzer is of the right type wrt the calibration information in the database (as defined in the calibrationMode_ variable). Should be implemented in each analyzer, if false will not do anything in analyze loop. default returns true.
141  virtual bool checkCorrectCalibrationType();
142 
143  //called when new DetID discovered
144  virtual void newDetID(uint32_t detid);
145 
146  void analyze(const edm::Event&, const edm::EventSetup&);
147  // the endJob method is used to save things like histograms etc.
148  // for additional actions derived classes should use the calibrationEnd() method for endJob functionality.
149  void endJob() ;
150  // calibrationEnd() will be used by derived classes
151  virtual void calibrationEnd() ;
152 
153  // checkPixel returns whether a particular pixel is to be expected during the entire run..
154  bool checkPixel(uint32_t detid, short row, short column);
155 };
156 
edm::ESHandle< SiPixelCalibConfiguration > calib_
virtual void beginRun(const edm::Run &, const edm::EventSetup &)
std::map< uint32_t, std::string > detIdNames_
std::string translateDetIdToString(uint32_t detid)
MonitorElement * bookDQMHistoPlaquetteSummary2D(uint32_t detid, std::string name, std::string title)
edm::ESHandle< TrackerGeometry > geom_
virtual void calibrationSetup(const edm::EventSetup &iSetup)
MonitorElement * bookDQMHistogram1D(uint32_t detid, std::string name, std::string title, int nchX, double lowX, double highX)
virtual bool doFits(uint32_t detid, std::vector< SiPixelCalibDigi >::const_iterator ipix)
void analyze(const edm::Event &, const edm::EventSetup &)
static const std::vector< short > * getVcalValues()
bool checkPixel(uint32_t detid, short row, short column)
void addTF1ToDQMMonitoringElement(MonitorElement *ele, TF1 *func)
edm::ESHandle< SiPixelFedCablingMap > theCablingMap_
MonitorElement * bookDQMHistogram2D(uint32_t detid, std::string name, std::string title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
edm::EDGetTokenT< edm::DetSetVector< SiPixelCalibDigi > > tPixelCalibDigi
Definition: Run.h:42