Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:54

0001 #include "DQM/SiStripMonitorSummary/interface/SiStripBackPlaneCorrectionDQM.h"
0002 #include "DataFormats/TrackerCommon/interface/SiStripSubStructure.h"
0003 #include "TCanvas.h"
0004 
0005 // -----
0006 SiStripBackPlaneCorrectionDQM::SiStripBackPlaneCorrectionDQM(
0007     edm::ESGetToken<SiStripBackPlaneCorrection, SiStripBackPlaneCorrectionRcd> token,
0008     edm::RunNumber_t iRun,
0009     edm::ParameterSet const &hPSet,
0010     edm::ParameterSet const &fPSet,
0011     const TrackerTopology *tTopo,
0012     const TkDetMap *tkDetMap)
0013     : SiStripBaseCondObjDQMGet<SiStripBackPlaneCorrection, SiStripBackPlaneCorrectionRcd>{
0014           token, iRun, hPSet, fPSet, tTopo} {
0015   if (HistoMaps_On_) {
0016     Tk_HM_ = std::make_unique<TkHistoMap>(tkDetMap, "SiStrip/Histo_Map", "BP_TkMap", 0.);
0017   }
0018 }
0019 // -----
0020 
0021 // -----
0022 SiStripBackPlaneCorrectionDQM::~SiStripBackPlaneCorrectionDQM() {}
0023 // -----
0024 
0025 // -----
0026 void SiStripBackPlaneCorrectionDQM::getActiveDetIds(const edm::EventSetup &eSetup) {
0027   getConditionObject(eSetup);
0028 
0029   std::map<uint32_t, float>::const_iterator BPMapIter_;
0030   std::map<uint32_t, float> BPMap_ = condObj_->getBackPlaneCorrections();
0031 
0032   for (BPMapIter_ = BPMap_.begin(); BPMapIter_ != BPMap_.end(); BPMapIter_++) {
0033     activeDetIds.push_back((*BPMapIter_).first);
0034   }
0035 }
0036 // -----
0037 
0038 // -----
0039 void SiStripBackPlaneCorrectionDQM::fillSummaryMEs(const std::vector<uint32_t> &selectedDetIds) {
0040   // -----
0041   // BP on layer-level : fill at once all detIds belonging to same layer when
0042   // encountering first detID in the layer
0043 
0044   bool fillNext = true;
0045   for (unsigned int i = 0; i < selectedDetIds.size(); i++) {
0046     int subDetId_ = ((selectedDetIds[i] >> 25) & 0x7);
0047     if (subDetId_ < 3 || subDetId_ > 6) {
0048       edm::LogError("SiStripBackPlaneCorrection")
0049           << "[SiStripBackPlaneCorrection::fillSummaryMEs] WRONG INPUT : no "
0050              "such subdetector type : "
0051           << subDetId_ << " and detId " << selectedDetIds[i] << " therefore no filling!" << std::endl;
0052     } else if (SummaryOnLayerLevel_On_) {
0053       if (fillNext) {
0054         fillMEsForLayer(/*SummaryMEsMap_,*/ selectedDetIds[i]);
0055       }
0056       if (getLayerNameAndId(selectedDetIds[i + 1]) == getLayerNameAndId(selectedDetIds[i])) {
0057         fillNext = false;
0058       } else {
0059         fillNext = true;
0060       }
0061     } else if (SummaryOnStringLevel_On_) {
0062       if (fillNext) {
0063         fillMEsForLayer(/*SummaryMEsMap_,*/ selectedDetIds[i]);
0064       }
0065       if (getStringNameAndId(selectedDetIds[i + 1]) == getStringNameAndId(selectedDetIds[i])) {
0066         fillNext = false;
0067       } else {
0068         fillNext = true;
0069       }
0070     }
0071   }
0072 
0073   for (std::map<uint32_t, ModMEs>::iterator iter = SummaryMEsMap_.begin(); iter != SummaryMEsMap_.end(); iter++) {
0074     ModMEs selME;
0075     selME = iter->second;
0076 
0077     if (SummaryOnStringLevel_On_) {
0078       if (fPSet_.getParameter<bool>("OutputSummaryProfileAtLayerLevelAsImage")) {
0079         TCanvas c1("c1");
0080         selME.SummaryOfProfileDistr->getTProfile()->Draw();
0081         std::string name(selME.SummaryOfProfileDistr->getTProfile()->GetTitle());
0082         name += ".png";
0083         c1.Print(name.c_str());
0084       }
0085 
0086       if (fPSet_.getParameter<bool>("OutputCumulativeSummaryAtLayerLevelAsImage")) {
0087         TCanvas c2("c2");
0088         selME.SummaryOfCumulDistr->getTH1()->Draw();
0089         std::string name2(selME.SummaryOfCumulDistr->getTitle());
0090         name2 += ".png";
0091         c2.Print(name2.c_str());
0092       }
0093 
0094     } else {
0095       if (hPSet_.getParameter<bool>("FillSummaryProfileAtLayerLevel") &&
0096           fPSet_.getParameter<bool>("OutputSummaryProfileAtLayerLevelAsImage")) {
0097         TCanvas c1("c1");
0098         selME.SummaryOfProfileDistr->getTProfile()->Draw();
0099         std::string name(selME.SummaryOfProfileDistr->getTProfile()->GetTitle());
0100         name += ".png";
0101         c1.Print(name.c_str());
0102       }
0103 
0104       if (hPSet_.getParameter<bool>("FillCumulativeSummaryAtLayerLevel") &&
0105           fPSet_.getParameter<bool>("OutputCumulativeSummaryAtLayerLevelAsImage")) {
0106         TCanvas c1("c1");
0107         selME.SummaryOfCumulDistr->getTH1()->Draw();
0108         std::string name(selME.SummaryOfCumulDistr->getTitle());
0109         name += ".png";
0110         c1.Print(name.c_str());
0111       }
0112     }
0113   }
0114 }
0115 // -----
0116 
0117 // -----
0118 void SiStripBackPlaneCorrectionDQM::fillMEsForLayer(
0119     /*std::map<uint32_t, ModMEs> selMEsMap_,*/ uint32_t selDetId_) {
0120   SiStripHistoId hidmanager;
0121 
0122   std::string hSummaryOfProfile_description;
0123   hSummaryOfProfile_description = hPSet_.getParameter<std::string>("SummaryOfProfile_description");
0124 
0125   std::string hSummary_name;
0126 
0127   int subDetId_ = ((selDetId_ >> 25) & 0x7);
0128 
0129   if (subDetId_ < 3 || subDetId_ > 6) {
0130     edm::LogError("SiStripBackPlaneCorrectionDQM")
0131         << "[SiStripBackPlaneCorrectionDQM::fillMEsForLayer] WRONG INPUT : no "
0132            "such subdetector type : "
0133         << subDetId_ << " no folder set!" << std::endl;
0134     return;
0135   }
0136 
0137   uint32_t selSubDetId_ = ((selDetId_ >> 25) & 0x7);
0138 
0139   std::vector<uint32_t> sameLayerDetIds_;
0140   sameLayerDetIds_.clear();
0141 
0142   if (SummaryOnStringLevel_On_) {  // FILLING FOR STRING LEVEL
0143 
0144     hSummary_name =
0145         hidmanager.createHistoLayer(hSummaryOfProfile_description, "layer", getStringNameAndId(selDetId_).first, "");
0146     std::map<uint32_t, ModMEs>::iterator selMEsMapIter_ = SummaryMEsMap_.find(getStringNameAndId(selDetId_).second);
0147 
0148     ModMEs selME_;
0149     if (selMEsMapIter_ != SummaryMEsMap_.end())
0150       selME_ = selMEsMapIter_->second;
0151 
0152     getSummaryMEs(selME_, selDetId_);
0153 
0154     // -----
0155     sameLayerDetIds_.clear();
0156 
0157     if (selSubDetId_ == 3) {  //  TIB
0158       if (tTopo_->tibIsInternalString(selDetId_)) {
0159         SiStripSubStructure::getTIBDetectors(
0160             activeDetIds, sameLayerDetIds_, tTopo_, tTopo_->tibLayer(selDetId_), 0, 1, tTopo_->tibString(selDetId_));
0161       }
0162       if (tTopo_->tibIsExternalString(selDetId_)) {
0163         SiStripSubStructure::getTIBDetectors(
0164             activeDetIds, sameLayerDetIds_, tTopo_, tTopo_->tibLayer(selDetId_), 0, 2, tTopo_->tibString(selDetId_));
0165       }
0166     } else if (selSubDetId_ == 4) {  // TID
0167       SiStripSubStructure::getTIDDetectors(activeDetIds, sameLayerDetIds_, tTopo_, 0, 0, 0, 0);
0168     } else if (selSubDetId_ == 5) {  // TOB
0169       SiStripSubStructure::getTOBDetectors(
0170           activeDetIds, sameLayerDetIds_, tTopo_, tTopo_->tobLayer(selDetId_), 0, tTopo_->tobRod(selDetId_));
0171     } else if (selSubDetId_ == 6) {  // TEC
0172       SiStripSubStructure::getTECDetectors(activeDetIds, sameLayerDetIds_, tTopo_, 0, 0, 0, 0, 0, 0);
0173     }
0174 
0175     // -----
0176 
0177     for (unsigned int i = 0; i < sameLayerDetIds_.size(); i++) {
0178       selME_.SummaryOfProfileDistr->Fill(i + 1, condObj_->getBackPlaneCorrection(sameLayerDetIds_[i]));
0179 
0180       // Fill the Histo_TkMap+TkMap with the BP:
0181       if (HistoMaps_On_)
0182         Tk_HM_->fill(sameLayerDetIds_[i], condObj_->getBackPlaneCorrection(sameLayerDetIds_[i]));
0183 
0184       std::cout << sameLayerDetIds_[i] << "\t" << condObj_->getBackPlaneCorrection(sameLayerDetIds_[i]) << std::endl;
0185 
0186       if (fPSet_.getParameter<bool>("TkMap_On") || hPSet_.getParameter<bool>("TkMap_On")) {
0187         fillTkMap(sameLayerDetIds_[i], condObj_->getBackPlaneCorrection(sameLayerDetIds_[i]));
0188       }
0189     }
0190 
0191     std::string hSummaryOfCumul_description;
0192     hSummaryOfCumul_description = hPSet_.getParameter<std::string>("SummaryOfCumul_description");
0193 
0194     std::string hSummaryOfCumul_name;
0195 
0196     if (subDetId_ < 3 || subDetId_ > 6) {
0197       edm::LogError("SiStripBackPlaneCorrectionDQM")
0198           << "[SiStripBackPlaneCorrectionDQM::fillMEsForLayer] WRONG INPUT : "
0199              "no such subdetector type : "
0200           << subDetId_ << " no folder set!" << std::endl;
0201       return;
0202     }
0203 
0204     hSummaryOfCumul_name =
0205         hidmanager.createHistoLayer(hSummaryOfCumul_description, "layer", getStringNameAndId(selDetId_).first, "");
0206 
0207     for (unsigned int i = 0; i < sameLayerDetIds_.size(); i++) {
0208       selME_.SummaryOfCumulDistr->Fill(condObj_->getBackPlaneCorrection(sameLayerDetIds_[i]));
0209     }
0210   }  // FILLING FOR STRING LEVEL
0211 
0212   else {  // FILLING FOR LAYER LEVEL
0213 
0214     std::map<uint32_t, ModMEs>::iterator selMEsMapIter_ = SummaryMEsMap_.find(getLayerNameAndId(selDetId_).second);
0215 
0216     ModMEs selME_;
0217     if (selMEsMapIter_ != SummaryMEsMap_.end())
0218       selME_ = selMEsMapIter_->second;
0219 
0220     getSummaryMEs(selME_, selDetId_);
0221 
0222     if (hPSet_.getParameter<bool>("FillSummaryProfileAtLayerLevel")) {
0223       hSummary_name =
0224           hidmanager.createHistoLayer(hSummaryOfProfile_description, "layer", getLayerNameAndId(selDetId_).first, "");
0225 
0226       // -----
0227       sameLayerDetIds_.clear();
0228 
0229       sameLayerDetIds_ = GetSameLayerDetId(activeDetIds, selDetId_);
0230 
0231       for (unsigned int i = 0; i < sameLayerDetIds_.size(); i++) {
0232         selME_.SummaryOfProfileDistr->Fill(i + 1, condObj_->getBackPlaneCorrection(sameLayerDetIds_[i]));
0233 
0234         // Fill the Histo_TkMap with BP:
0235         if (HistoMaps_On_)
0236           Tk_HM_->fill(sameLayerDetIds_[i], condObj_->getBackPlaneCorrection(sameLayerDetIds_[i]));
0237 
0238         if (fPSet_.getParameter<bool>("TkMap_On") || hPSet_.getParameter<bool>("TkMap_On")) {
0239           fillTkMap(sameLayerDetIds_[i], condObj_->getBackPlaneCorrection(sameLayerDetIds_[i]));
0240         }
0241       }
0242     }  // if Fill ...
0243 
0244     if (hPSet_.getParameter<bool>("FillCumulativeSummaryAtLayerLevel")) {
0245       std::string hSummaryOfCumul_description;
0246       hSummaryOfCumul_description = hPSet_.getParameter<std::string>("SummaryOfCumul_description");
0247 
0248       std::string hSummaryOfCumul_name;
0249 
0250       if (subDetId_ < 3 || subDetId_ > 6) {
0251         edm::LogError("SiStripBackPlaneCorrectionDQM")
0252             << "[SiStripBackPlaneCorrectionDQM::fillMEsForLayer] WRONG INPUT : "
0253                "no such subdetector type : "
0254             << subDetId_ << " no folder set!" << std::endl;
0255         return;
0256       }
0257 
0258       hSummaryOfCumul_name =
0259           hidmanager.createHistoLayer(hSummaryOfCumul_description, "layer", getLayerNameAndId(selDetId_).first, "");
0260 
0261       for (unsigned int i = 0; i < sameLayerDetIds_.size(); i++) {
0262         selME_.SummaryOfCumulDistr->Fill(condObj_->getBackPlaneCorrection(sameLayerDetIds_[i]));
0263       }
0264     }  // if Fill ...
0265   }    // FILLING FOR LAYER LEVEL
0266 }
0267 // -----