Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-10-25 09:36:37

0001 #include "CondCore/Utilities/interface/PayloadInspectorModule.h"
0002 #include "CondCore/Utilities/interface/PayloadInspector.h"
0003 #include "CondCore/CondDB/interface/Time.h"
0004 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0005 #include "CondCore/HcalPlugins/interface/HcalObjRepresent.h"
0006 
0007 // the data format of the condition to be inspected
0008 #include "CondFormats/HcalObjects/interface/HcalRecoParams.h"
0009 
0010 #include "TH2F.h"
0011 #include "TCanvas.h"
0012 #include "TLine.h"
0013 #include "TStyle.h"
0014 #include "TLatex.h"
0015 #include "TPave.h"
0016 #include "TPaveStats.h"
0017 #include <string>
0018 #include <fstream>
0019 
0020 namespace {
0021 
0022   /********************************************
0023      printing float values of reco paramters
0024   *********************************************/
0025   class HcalRecoParamsSummary : public cond::payloadInspector::PlotImage<HcalRecoParams> {
0026   public:
0027     HcalRecoParamsSummary() : cond::payloadInspector::PlotImage<HcalRecoParams>("HCAL RecoParam Ratios - map ") {
0028       setSingleIov(true);
0029     }
0030 
0031     bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs) override {
0032       auto iov = iovs.front();
0033       std::shared_ptr<HcalRecoParams> payload = fetchPayload(std::get<1>(iov));
0034       if (payload.get()) {
0035         std::string subDetName;
0036         std::vector<HcalRecoParam> itemsVec;
0037         std::pair<std::string, int> valMap;
0038 
0039         //TODO: Abstract into a function that takes valMap as the argument
0040 
0041         TLatex label, val;
0042         TLine* ll;
0043         TLine* lr;
0044         TLine* lt;
0045         TLine* lb;
0046         TCanvas* can = new TCanvas("RecoParamsSummary", "RecoParamsSummary", 1680, 1680);
0047         //can->cd();
0048         //HcalObjRepresent::drawTable(2,2);
0049         can->Divide(2, 2, 0, 0);
0050         int i = 1;
0051         int psID;
0052         label.SetNDC();
0053         label.SetTextAlign(26);
0054         label.SetTextSize(0.05);
0055         label.SetTextColor(2);
0056         label.DrawLatex(0.5, 0.96, Form("Hcal Pulse Shape IDs"));
0057 
0058         for (std::pair<std::string, std::vector<HcalRecoParam> > cont : (*payload).getAllContainers()) {
0059           psID = 0;
0060           subDetName = std::get<0>(cont);
0061           if (subDetName[0] != 'H' || subDetName == "HT")
0062             continue;
0063           itemsVec = std::get<1>(cont);
0064           //valMap.insert(std::make_pair(subDetName,itemsVec.front().pulseShapeID()));
0065           can->cd(i);
0066           ll = new TLine(0, 0, 0, 1);
0067           ll->SetLineWidth(4);
0068           ll->Draw();
0069           lt = new TLine(0, 1, 1, 1);
0070           lt->SetLineWidth(4);
0071           lt->Draw();
0072           lb = new TLine(0, 0, 1, 0);
0073           lb->SetLineWidth(4);
0074           lb->Draw();
0075           lr = new TLine(1, 0, 1, 1);
0076           lr->SetLineWidth(4);
0077           lr->Draw();
0078           label.SetNDC();
0079           label.SetTextAlign(26);
0080           label.SetTextSize(0.15);
0081           //label.SetTextColor(2);
0082           label.DrawLatex(0.5, 0.75, subDetName.c_str());
0083           val.SetNDC();
0084           val.SetTextAlign(26);
0085           val.SetTextSize(0.1);
0086           //        val.SetTextColor(1);
0087           std::vector<HcalRecoParam>::iterator it;
0088           for (it = itemsVec.begin(); it != itemsVec.end(); it++) {
0089             psID = (*it).pulseShapeID();
0090             if (psID != 0) {
0091               psID = (*it).pulseShapeID();
0092               break;
0093             }
0094           }
0095           val.DrawLatex(0.5, 0.25, std::to_string(psID).c_str());
0096           i++;
0097         }
0098 
0099         std::string ImageName(m_imageFileName);
0100         can->SaveAs(ImageName.c_str());
0101         return false;
0102       } else
0103         return false;
0104     }  // fill method
0105   };
0106   //TODO: Add a Change Summary?
0107 
0108 }  // namespace
0109 
0110 // Register the classes as boost python plugin
0111 PAYLOAD_INSPECTOR_MODULE(HcalRecoParams) { PAYLOAD_INSPECTOR_CLASS(HcalRecoParamsSummary); }