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