CMS 3D CMS Logo

CTPPSPixGainCalibsESAnalyzer.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <map>
12 #include "TH2D.h"
13 #include "TFile.h"
14 
16 public:
18  : m_outfilename(p.getUntrackedParameter<std::string>("outputrootfile", "output.root")) {
19  // std::cout<<"CTPPSPixGainCalibsESAnalyzer"<<std::endl;
21  }
23  //std::cout<<"CTPPSPixGainCalibsESAnalyzer "<<i<<std::endl;
25  }
27  //std::cout<<"~CTPPSPixGainCalibsESAnalyzer "<<std::endl;
28  }
29  void analyze(const edm::Event& e, const edm::EventSetup& c) override;
30  void setReadablePlaneNames();
31 
32 private:
33  std::map<uint32_t, std::string> detId_readable;
35 };
36 
38  detId_readable[2014838784] = "Arm_0_Sec_45_St_0_Pot_3_Plane_0";
39  detId_readable[2014904320] = "Arm_0_Sec_45_St_0_Pot_3_Plane_1";
40  detId_readable[2014969856] = "Arm_0_Sec_45_St_0_Pot_3_Plane_2";
41  detId_readable[2015035392] = "Arm_0_Sec_45_St_0_Pot_3_Plane_3";
42  detId_readable[2015100928] = "Arm_0_Sec_45_St_0_Pot_3_Plane_4";
43  detId_readable[2015166464] = "Arm_0_Sec_45_St_0_Pot_3_Plane_5";
44  detId_readable[2023227392] = "Arm_0_Sec_45_St_2_Pot_3_Plane_0";
45  detId_readable[2023292928] = "Arm_0_Sec_45_St_2_Pot_3_Plane_1";
46  detId_readable[2023358464] = "Arm_0_Sec_45_St_2_Pot_3_Plane_2";
47  detId_readable[2023424000] = "Arm_0_Sec_45_St_2_Pot_3_Plane_3";
48  detId_readable[2023489536] = "Arm_0_Sec_45_St_2_Pot_3_Plane_4";
49  detId_readable[2023555072] = "Arm_0_Sec_45_St_2_Pot_3_Plane_5";
50  detId_readable[2031616000] = "Arm_1_Sec_56_St_0_Pot_3_Plane_0";
51  detId_readable[2031681536] = "Arm_1_Sec_56_St_0_Pot_3_Plane_1";
52  detId_readable[2031747072] = "Arm_1_Sec_56_St_0_Pot_3_Plane_2";
53  detId_readable[2031812608] = "Arm_1_Sec_56_St_0_Pot_3_Plane_3";
54  detId_readable[2031878144] = "Arm_1_Sec_56_St_0_Pot_3_Plane_4";
55  detId_readable[2031943680] = "Arm_1_Sec_56_St_0_Pot_3_Plane_5";
56  detId_readable[2040004608] = "Arm_1_Sec_56_St_2_Pot_3_Plane_0";
57  detId_readable[2040070144] = "Arm_1_Sec_56_St_2_Pot_3_Plane_1";
58  detId_readable[2040135680] = "Arm_1_Sec_56_St_2_Pot_3_Plane_2";
59  detId_readable[2040201216] = "Arm_1_Sec_56_St_2_Pot_3_Plane_3";
60  detId_readable[2040266752] = "Arm_1_Sec_56_St_2_Pot_3_Plane_4";
61  detId_readable[2040332288] = "Arm_1_Sec_56_St_2_Pot_3_Plane_5";
62 }
63 
65  edm::LogPrint("CTPPSPixGainCalibsReader") << "###CTPPSPixGainCalibsESAnalyzer::analyze";
66  edm::LogPrint("CTPPSPixGainCalibsReader") << " I AM IN RUN NUMBER " << e.id().run();
67  edm::LogPrint("CTPPSPixGainCalibsReader") << " ---EVENT NUMBER " << e.id().event();
69  edm::eventsetup::EventSetupRecordKey::TypeTag::findType("CTPPSPixelGainCalibrationsRcd"));
71  //record not found
72  edm::LogPrint("CTPPSPixGainCalibsReader") << "Record \"CTPPSPixelGainCalibrationsRcd"
73  << "\" does not exist ";
74  }
75  //this part gets the handle of the event source and the record (i.e. the Database)
77  edm::LogPrint("CTPPSPixGainCalibsReader") << "got eshandle";
78  context.get<CTPPSPixelGainCalibrationsRcd>().get(calhandle);
79  edm::LogPrint("CTPPSPixGainCalibsReader") << "got context";
80  const CTPPSPixelGainCalibrations* pPixelGainCalibrations = calhandle.product();
81  edm::LogPrint("CTPPSPixGainCalibsReader") << "got CTPPSPixelGainCalibrations* ";
82  edm::LogPrint("CTPPSPixGainCalibsReader") << "print pointer address : ";
83  edm::LogPrint("CTPPSPixGainCalibsReader") << pPixelGainCalibrations;
84 
85  TFile myfile(m_outfilename.c_str(), "RECREATE");
86  myfile.cd();
87 
88  // the pPixelGainCalibrations object contains the map of detIds to pixel gains and pedestals for current run
89  // we get the map just to loop over the contents, but from here on it should be as the code (reconstruction etc) needs.
90  // Probably best to check that the key (detid) is in the list before calling the data
91 
92  edm::LogPrint("CTPPSPixGainCalibsReader") << "Size " << pPixelGainCalibrations->size();
93  const CTPPSPixelGainCalibrations::CalibMap& mymap = pPixelGainCalibrations->getCalibMap(); //just to get the keys?
94 
95  for (CTPPSPixelGainCalibrations::CalibMap::const_iterator it = mymap.begin(); it != mymap.end(); ++it) {
96  uint32_t detId = it->first;
97 
98  edm::LogPrint("CTPPSPixGainCalibsReader")
99  << "Address of detId = " << (&detId) << " and of it = " << (&it) << " and of it->first = " << (&(it->first));
100 
101  edm::LogPrint("CTPPSPixGainCalibsReader") << "Content of pPixelGainCalibrations for key: detId= " << detId;
102  CTPPSPixelGainCalibration mycalibs0 = pPixelGainCalibrations->getGainCalibration(detId);
103  const CTPPSPixelGainCalibration& mycalibs = it->second;
104 
105  edm::LogPrint("CTPPSPixGainCalibsReader")
106  << "Address of mycalibs0 = " << (&mycalibs0) << " and of mycalibs = " << (&mycalibs) << " and of it->second "
107  << (&(it->second));
108 
109  std::string namep("pedsFromDB_" + detId_readable[detId]);
110  std::string nameg("gainsFromDB_" + detId_readable[detId]);
111  std::string tlp("Pedestals for " + detId_readable[detId] + "; column; row");
112  std::string tlg("Gains for " + detId_readable[detId] + "; column; row");
113  TH2D mypeds(namep.c_str(), tlp.c_str(), 156, 0., 156., 160, 0., 160.);
114  TH2D mygains(nameg.c_str(), tlg.c_str(), 156, 0., 156., 160, 0., 160.);
115 
116  int ncols = mycalibs.getNCols();
117  int npix = mycalibs.getIEnd();
118  int nrows = mycalibs.getNRows(); //should be == 160
119  edm::LogPrint("CTPPSPixGainCalibsReader") << "Here ncols = " << ncols << " nrows =" << nrows << " npix=" << npix;
120  for (int jrow = 0; jrow < nrows; ++jrow)
121  for (int icol = 0; icol < ncols; ++icol) {
122  if (mycalibs.isDead(icol + jrow * ncols)) {
123  edm::LogPrint("CTPPSPixGainCalibsReader") << "Dead Pixel icol =" << icol << " jrow =" << jrow;
124  continue;
125  }
126  if (mycalibs.isNoisy(icol + jrow * ncols)) {
127  edm::LogPrint("CTPPSPixGainCalibsReader") << "Noisy Pixel icol =" << icol << " jrow =" << jrow;
128  continue;
129  }
130  mygains.Fill(icol, jrow, mycalibs.getGain(icol, jrow));
131  mypeds.Fill(icol, jrow, mycalibs.getPed(icol, jrow));
132  }
133 
134  mypeds.Write();
135  mygains.Write();
136  }
137  myfile.Write();
138  myfile.Close();
139 }
RunNumber_t run() const
Definition: EventID.h:38
EventNumber_t event() const
Definition: EventID.h:40
const CalibMap & getCalibMap() const
CTPPSPixGainCalibsESAnalyzer(edm::ParameterSet const &p)
bool isDead(const uint32_t ipix) const
const CTPPSPixelGainCalibration & getGainCalibration(const uint32_t &detid) const
static HCTypeTag findType(char const *iTypeName)
find a type based on the types name, if not found will return default HCTypeTag
void analyze(const edm::Event &e, const edm::EventSetup &c) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::map< uint32_t, std::string > detId_readable
std::map< uint32_t, CTPPSPixelGainCalibration > CalibMap
heterocontainer::HCTypeTag TypeTag
edm::EventID id() const
Definition: EventBase.h:59
float getPed(const int &col, const int &row) const
T get() const
Definition: EventSetup.h:73
bool isNoisy(const uint32_t ipix) const
T const * product() const
Definition: ESHandle.h:86
float getGain(const int &col, const int &row) const