CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/DQM/SiStripMonitorHardware/plugins/BuildTrackerMap.cc

Go to the documentation of this file.
00001 
00002 // -*- C++ -*-
00003 //
00004 // Package:    DQM/SiStripMonitorHardware
00005 // Class:      BuildTrackerMapPlugin
00006 // 
00011 //
00012 //         Created:  2009/07/22
00013 // $Id: BuildTrackerMap.cc,v 1.11 2010/04/21 10:39:10 amagnan Exp $
00014 //
00015 
00016 #include <sstream>
00017 #include <fstream>
00018 #include <iostream>
00019 #include <memory>
00020 #include <list>
00021 #include <algorithm>
00022 #include <cassert>
00023 
00024 #include "TCanvas.h"
00025 #include "TH1F.h"
00026 #include "TStyle.h"
00027 #include "TPaveStats.h"
00028 
00029 #include "FWCore/Framework/interface/Frameworkfwd.h"
00030 #include "FWCore/Framework/interface/EDAnalyzer.h"
00031 #include "FWCore/Framework/interface/Event.h"
00032 #include "FWCore/Framework/interface/EventSetup.h"
00033 #include "FWCore/Framework/interface/ESHandle.h"
00034 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00035 #include "FWCore/ServiceRegistry/interface/Service.h"
00036 #include "FWCore/Utilities/interface/Exception.h"
00037 
00038 #include "CondFormats/DataRecord/interface/SiStripFedCablingRcd.h"
00039 #include "CondFormats/SiStripObjects/interface/SiStripFedCabling.h"
00040 #include "CondFormats/SiStripObjects/interface/FedChannelConnection.h"
00041 
00042 #include "CommonTools/TrackerMap/interface/TrackerMap.h"
00043 #include "CommonTools/TrackerMap/interface/TmModule.h"
00044 #include "CommonTools/TrackerMap/interface/TmApvPair.h"
00045 
00046 #include "DQM/SiStripCommon/interface/TkHistoMap.h" 
00047 
00048 #include "DQMServices/Core/interface/DQMStore.h"
00049 
00050 //
00051 // Class declaration
00052 //
00053 
00054 class BuildTrackerMapPlugin : public edm::EDAnalyzer
00055 {
00056  public:
00057 
00058   explicit BuildTrackerMapPlugin(const edm::ParameterSet&);
00059   ~BuildTrackerMapPlugin();
00060  private:
00061   virtual void beginJob();
00062   virtual void analyze(const edm::Event&, const edm::EventSetup&);
00063   virtual void endJob();
00064 
00065   void read(bool aMechView,
00066             std::string aFile,
00067             std::vector<TkHistoMap*> & aTkMapVec,
00068             std::vector<bool> & aValidVec);
00069   void subtractMap(TkHistoMap*& aResult,
00070                    TkHistoMap*& aSubtr);
00071 
00072 
00073 
00074   //input file names
00075   std::string fileName_;
00076   //do mechanical view or not
00077   bool mechanicalView_;
00078   //folder name for histograms in DQMStore
00079   std::string folderName_;
00080   //print debug messages when problems are found: 1=error debug, 2=light debug, 3=full debug
00081   unsigned int printDebug_;
00082 
00083   bool doDiff_;
00084   std::string fileNameDiff_;
00085 
00086   std::vector<TkHistoMap*> tkHistoMapVec_;
00087   std::vector<TkHistoMap*> tkHistoMapVecDiff_;
00088   
00089   //name of the tkHistoMap to extract
00090   std::vector<std::string> tkHistoMapNameVec_;
00091   std::vector<double> minVal_;
00092   std::vector<double> maxVal_;
00093 
00094   std::vector<bool> isValidMap_;
00095   std::vector<bool> isValidMapDiff_;
00096 
00097   edm::ParameterSet pset_;
00098   std::vector<TrackerMap*> tkmap_;
00099 
00100 };
00101 
00102 
00103 //
00104 // Constructors and destructor
00105 //
00106 
00107 BuildTrackerMapPlugin::BuildTrackerMapPlugin(const edm::ParameterSet& iConfig)
00108   : fileName_(iConfig.getUntrackedParameter<std::string>("InputFileName","DQMStore.root")),
00109     mechanicalView_(iConfig.getUntrackedParameter<bool>("MechanicalView",true)),
00110     folderName_(iConfig.getUntrackedParameter<std::string>("HistogramFolderName","DQMData/")),
00111     printDebug_(iConfig.getUntrackedParameter<unsigned int>("PrintDebugMessages",1)),
00112     doDiff_(iConfig.getUntrackedParameter<bool>("DoDifference",false)),
00113     fileNameDiff_(iConfig.getUntrackedParameter<std::string>("InputFileNameForDiff","DQMStore.root")),
00114     tkHistoMapNameVec_(iConfig.getUntrackedParameter<std::vector<std::string> >("TkHistoMapNameVec")),
00115     minVal_(iConfig.getUntrackedParameter<std::vector<double> >("MinValueVec")),
00116     maxVal_(iConfig.getUntrackedParameter<std::vector<double> >("MaxValueVec")),
00117     pset_(iConfig.getParameter<edm::ParameterSet>("TkmapParameters"))
00118 {
00119 
00120   read(mechanicalView_,fileName_,tkHistoMapVec_,isValidMap_);
00121   if (doDiff_) read(mechanicalView_,fileNameDiff_,tkHistoMapVecDiff_,isValidMapDiff_);
00122 
00123 
00124 //   for (unsigned int i(0); i<34; i++){
00125 //     if (i<4) histName_[i] << "TIB/layer_" << i+1 << "/" << tkDetMapName_ << "_TIB_L" << i+1;
00126 //     else if (i<7)  histName_[i] << "TID/side_1/wheel_" << i-3 << "/" << tkDetMapName_ << "_TIDM_D" << i-3;
00127 //     else if (i<10)  histName_[i] << "TID/side_2/wheel_" << i-6 << "/" << tkDetMapName_ << "_TIDP_D" << i-6;
00128 //     else if (i<16) histName_[i] << "TOB/layer_" << i-9 << "/" << tkDetMapName_ << "_TOB_L" << i-9;
00129 //     else if (i<25) histName_[i] << "TEC/side_1/wheel_" << i-15 << "/" << tkDetMapName_ << "_TECM_W" << i-15;
00130 //     else if (i<34) histName_[i] << "TEC/side_2/wheel_" << i-24 << "/" << tkDetMapName_ << "_TECP_W" << i-24;
00131 
00132 //     std::cout << "histName[" << i << "] =" << histName_[i] << std::endl;
00133 
00134 //   } 
00135   
00136 }
00137 
00138 BuildTrackerMapPlugin::~BuildTrackerMapPlugin()
00139 {
00140 
00141   tkHistoMapVec_.clear();
00142   if (doDiff_) tkHistoMapVecDiff_.clear();
00143 }
00144 
00145 
00146 //
00147 // Member functions
00148 //
00149 
00150 /*Check that is possible to load in tkhistomaps histograms already stored in a DQM root file (if the folder and name are known)*/
00151 void BuildTrackerMapPlugin::read(bool aMechView,
00152                                  std::string aFile,
00153                                  std::vector<TkHistoMap*> & aTkMapVec,
00154                                  std::vector<bool> & aValidVec)
00155 {
00156   
00157   DQMStore * lDqmStore = edm::Service<DQMStore>().operator->();
00158   lDqmStore->open(aFile);  
00159   std::vector<TkHistoMap *> tkHistoMap;
00160 
00161   unsigned int nHists = tkHistoMapNameVec_.size();
00162   tkHistoMap.resize(nHists,0);
00163   aValidVec.resize(nHists,true);
00164 
00165   std::string dirName = folderName_;
00166   if (dirName == "") {
00167     dirName += "Run ";
00168     dirName += aFile.substr(aFile.find_last_of("_")+5,6);
00169     dirName += "/SiStrip/Run summary";
00170     std::cout << " -- DirName = " << dirName << std::endl;
00171    }
00172 
00173   //lDqmStore->setCurrentFolder(dirName);
00174   //lDqmStore->showDirStructure();
00175   
00176   unsigned int nFailTot=0;
00177   unsigned int nTotTot = 0;
00178   for (unsigned int i(0); i<nHists; i++){
00179 
00180     tkHistoMap[i] = new TkHistoMap();
00181 
00182     tkHistoMap[i]->loadTkHistoMap(dirName,tkHistoMapNameVec_.at(i),aMechView);
00183     
00184     std::vector<MonitorElement*>& lMaps = tkHistoMap[i]->getAllMaps();
00185 
00186     std::cout << " -- map " << i << ", nHistos = " << lMaps.size() << std::endl;
00187     unsigned int nFail=0;
00188     unsigned int nTot=0;
00189  
00190     for (unsigned int im(0); im<lMaps.size(); im++){
00191       if (!lMaps[im]) {
00192         std::cout << " -- Failed to get element " << im << " for map " << i << std::endl;
00193         nFail++;
00194         nFailTot++;
00195       }
00196       nTot++;
00197       nTotTot++;
00198     }
00199 
00200     if (nFail == nTot) aValidVec[i] = false;
00201     aTkMapVec.push_back(tkHistoMap[i]);
00202   }
00203 
00204   if (nFailTot < nTotTot) std::cout << " - " << nTotTot-nFailTot << "/" << nTotTot 
00205                               << " histomaps read with success for file ." << aFile << std::endl;
00206   else {
00207     std::cout << " - Failed to read any map for file " << aFile << ". Exiting line ... " << __LINE__ << std::endl;
00208      exit(1);
00209   }
00210 
00211 //   //get list of detid for which |deltaRMS(APV0-APV1)|>1
00212 //   unsigned int lHistoNumber = 35;
00213 //   TkDetMap lTkDetMap;
00214 //   std::ofstream list,listRms0,listRms1;
00215 //   list.open("./cmBadModuleList.dat",std::ios::out);
00216 //   listRms0.open("./cmBadModuleList_rms0.dat",std::ios::out);
00217 //   listRms1.open("./cmBadModuleList_rms1.dat",std::ios::out);
00218 //   if (!list || !listRms0 || !listRms1) {
00219 //     std::cout << "Warning, can't open output file to write bad module list !" << std::endl;
00220 //     exit(1);
00221 //   }
00222 
00223 //   TCanvas *lCan = new TCanvas("lCan","",1);
00224 //   TH1F *p_deltaMean = new TH1F("p_deltaMean",";CM_{mean}(APV0)-CM_{mean}(APV1)",500,-2,2);
00225 //   TH1F *p_deltaRMS = new TH1F("p_deltaRMS",";CM_{RMS}(APV0)-CM_{RMS}(APV1)",500,0,3);
00226 //   TH1F *p_MeanAPV0 = new TH1F("p_MeanAPV0",";CM_{mean}(APV0)",500,100,140);
00227 //   //TH1F *p_MeanAPV1 = new TH1F("p_MeanAPV1",";CM_{mean}(APV1)",500,100,140);
00228 //   TH1F *p_RMSAPV0 = new TH1F("p_RMSAPV0",";CM_{RMS}(APV0)",500,0,10);
00229 //   //TH1F *p_RMSAPV1 = new TH1F("p_RMSAPV1",";CM_{RMS}(APV1)",500,0,10);
00230 
00231 
00232 
00233 //   gStyle->SetOptStat(1111111);
00234 
00235 //   for(unsigned int layer=1;layer<lHistoNumber;++layer){
00236 //     std::vector<uint32_t> dets;
00237 //     lTkDetMap.getDetsForLayer(layer,dets);
00238 //     for(size_t i=0;i<dets.size();++i){
00239 //       if(dets[i]>0){
00240 //      //if(tkHistoMap[5]->getEntries(dets[i])>0 && tkHistoMap[5]->getValue(dets[i])) {
00241 //      if(nHists > 3){
00242 //        if (tkHistoMap[3]->getValue(dets[i]) > 1) {
00243 //          list << dets[i] << " " << tkHistoMap[3]->getValue(dets[i]) << std::endl;
00244 //        }
00245 //      p_deltaRMS->Fill(tkHistoMap[3]->getValue(dets[i]));
00246 //      }
00247 //      p_MeanAPV0->Fill(tkHistoMap[0]->getValue(dets[i]));
00248 //      //p_MeanAPV1->Fill(tkHistoMap[1]->getValue(dets[i]));
00249 //      p_RMSAPV0->Fill(tkHistoMap[1]->getValue(dets[i]));
00250 //      if (tkHistoMap[1]->getValue(dets[i]) > 2)
00251 //        listRms0 << dets[i] << " " << tkHistoMap[1]->getValue(dets[i]) << std::endl;
00252 //      //p_RMSAPV1->Fill(tkHistoMap[3]->getValue(dets[i]));
00253 //      //if (tkHistoMap[3]->getValue(dets[i]) > 2)
00254 //      //listRms1 << dets[i] << " " << tkHistoMap[3]->getValue(dets[i]) << std::endl;
00255 
00256 //      if(nHists > 2) p_deltaMean->Fill(tkHistoMap[2]->getValue(dets[i]));
00257 //       }
00258 //     }
00259 //   }
00260 //   list.close();
00261 //   listRms0.close();
00262 //   listRms1.close();
00263 
00264 //   lCan->cd();
00265 //   p_deltaRMS->Draw();
00266 //   //lCan->Print("./deltaRMStotal.png");
00267 //   lCan->Print("./deltaRMStotal.C");
00268 
00269 //   p_deltaMean->Draw();
00270 //   lCan->Update();
00271 //   lCan->Print("./deltaMeantotal.C");
00272 
00273 //   TPaveStats *statBox[2] = {0,0};  
00274 //   statBox[0] = (TPaveStats*)p_MeanAPV0->FindObject("stats");
00275 //   //statBox[1] = (TPaveStats*)p_MeanAPV1->FindObject("stats");
00276 
00277 //   p_MeanAPV0->Draw();
00278 //   //p_MeanAPV1->SetLineColor(2);
00279 //   //p_MeanAPV1->Draw("same");
00280 //   if (statBox[0]) statBox[0]->Draw("same");
00281 //   if (statBox[1]) { 
00282 //     statBox[1]->SetLineColor(2);
00283 //     statBox[1]->SetTextColor(2);
00284 //     statBox[1]->Draw("same");
00285 //   }
00286 //   lCan->Update();
00287 //   lCan->Print("./meanAPVstotal.C");
00288 
00289 //   statBox[0] = (TPaveStats*)p_RMSAPV0->FindObject("stats");
00290 //   //statBox[1] = (TPaveStats*)p_RMSAPV1->FindObject("stats");
00291 
00292 //   p_RMSAPV0->Draw();
00293 //   //p_RMSAPV1->SetLineColor(2);
00294 //   //p_RMSAPV1->Draw("same");
00295 //   if (statBox[0]) statBox[0]->Draw("same");
00296 //   if (statBox[1]) { 
00297 //     statBox[1]->SetLineColor(2);
00298 //     statBox[1]->SetTextColor(2);
00299 //     statBox[1]->Draw("same");
00300 //   }
00301 //   lCan->Update();
00302 //   lCan->Print("./rmsAPVstotal.C");
00303 
00304 
00305 }
00306 
00307 
00308 // ------------ method called to for each event  ------------
00309 void
00310 BuildTrackerMapPlugin::analyze(const edm::Event& iEvent, 
00311                                  const edm::EventSetup& iSetup)
00312 {
00313 
00314   static bool firstEvent = true;
00315 
00316   edm::ESHandle<SiStripFedCabling> fedcabling;
00317   iSetup.get<SiStripFedCablingRcd>().get(fedcabling );
00318   
00319   if (firstEvent) {
00320     for (unsigned int i(0); i<tkHistoMapNameVec_.size(); i++){
00321       tkmap_.push_back(new TrackerMap(pset_,fedcabling));
00322     }
00323 
00324   }
00325 
00326   firstEvent = false;
00327 
00328   std::cout << "End of analyze method: tkmap_ size = " << tkmap_.size() << std::endl;
00329 
00330 }//analyze method
00331 
00332 // ------------ method called once each job just before starting event loop  ------------
00333 void 
00334 BuildTrackerMapPlugin::beginJob()
00335 {
00336 
00337 }
00338 
00339 // ------------ method called once each job just after ending the event loop  ------------
00340 void 
00341 BuildTrackerMapPlugin::endJob()
00342 {
00343   //edm::ESHandle<SiStripFedCabling> pDD1;
00344   //iSetup.get<SiStripFedCablingRcd>().get(pDD1);
00345   std::cout << "Processing endjob with " << tkHistoMapNameVec_.size()<< " elements." << std::endl;
00346 
00347   assert (minVal_.size() == tkHistoMapNameVec_.size());
00348   assert (maxVal_.size() == tkHistoMapNameVec_.size());
00349 
00350   for (unsigned int i(0); i<tkHistoMapNameVec_.size(); i++){
00351 
00352     std::cout << "Processing element " << i << ": " << tkHistoMapNameVec_.at(i) << std::endl;
00353     std::cout << "Min, max = " << minVal_.at(i) << " " << maxVal_.at(i) << std::endl;
00354 
00355     TrackerMap* lTkMap = tkmap_.at(i);
00356 
00357     if (!lTkMap) {
00358       std::cout << "tkmap_ is NULL for element " << i << "... continuing ..." << std::endl;
00359       continue;
00360     }
00361 
00362     subtractMap(tkHistoMapVec_.at(i),tkHistoMapVecDiff_.at(i));
00363 
00364 
00365     //(pset_,pDD1); 
00366     lTkMap->setPalette(1);
00367     lTkMap->showPalette(1);
00368     if (!tkHistoMapVec_.at(i) || !isValidMap_.at(i)) {
00369       std::cout << "Warning, tkHistoMap is invalid for element " << i << "... continuing ..." << std::endl;
00370       continue;
00371     }
00372     tkHistoMapVec_.at(i)->dumpInTkMap(lTkMap);
00373 
00374     //to print all figures to create fancy view
00375     //lTkMap->printall(true,0,255,tkHistoMapNameVec_.at(i));
00376     lTkMap->save(true,
00377                  minVal_.at(i),
00378                  maxVal_.at(i),
00379                  tkHistoMapNameVec_.at(i)+std::string(".png"));
00380     lTkMap->save_as_fedtrackermap(true,
00381                                   minVal_.at(i),
00382                                   maxVal_.at(i),
00383                                   tkHistoMapNameVec_.at(i)+std::string("_FED.png"));
00384 
00385   }
00386 
00387 }
00388 
00389 
00390 void BuildTrackerMapPlugin::subtractMap(TkHistoMap *& aResult,
00391                                         TkHistoMap *& aSubtr)
00392 {
00393   
00394   std::vector<MonitorElement*>& lMaps = aResult->getAllMaps();
00395   std::vector<MonitorElement*>& lMapsDiff = aSubtr->getAllMaps();
00396 
00397   assert(lMaps.size() == lMapsDiff.size());  
00398 
00399     for (unsigned int im(0); im<lMaps.size(); im++){
00400       if (!lMaps[im] || !lMapsDiff[im]) {
00401         std::cout << " -- Failed to get element " << im << " for maps." << std::endl;
00402       }
00403       else {
00404         (lMaps[im]->getTProfile2D())->Add(lMapsDiff[im]->getTProfile2D(),-1);
00405       }
00406     }
00407 
00408 }
00409 
00410 
00411 
00412 // Define as a plug-in
00413 //
00414 
00415 #include "FWCore/Framework/interface/MakerMacros.h"
00416 DEFINE_FWK_MODULE(BuildTrackerMapPlugin);