CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/DQMServices/Components/plugins/MEtoMEComparitor.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    MEtoMEComparitor
00004 // Class:      MEtoMEComparitor
00005 // 
00013 //
00014 // Original Author:  jean-roch Vlimant,40 3-A28,+41227671209,
00015 //         Created:  Tue Nov 30 18:55:50 CET 2010
00016 // $Id: MEtoMEComparitor.cc,v 1.5 2010/12/13 16:56:37 vlimant Exp $
00017 //
00018 //
00019 
00020 #include "MEtoMEComparitor.h"
00021 
00022 #include "classlib/utils/StringList.h"
00023 #include "classlib/utils/StringOps.h"
00024 
00025 
00026 MEtoMEComparitor::MEtoMEComparitor(const edm::ParameterSet& iConfig)
00027 
00028 {
00029   _moduleLabel = iConfig.getParameter<std::string>("MEtoEDMLabel");
00030   
00031   _lumiInstance = iConfig.getParameter<std::string>("lumiInstance");
00032   _runInstance = iConfig.getParameter<std::string>("runInstance");
00033 
00034   _process_ref = iConfig.getParameter<std::string>("processRef");
00035   _process_new = iConfig.getParameter<std::string>("processNew");
00036 
00037   if (iConfig.getParameter<bool>("autoProcess")){
00038     //get the last two process from the provenance
00039 
00040   }
00041   _KSgoodness = iConfig.getParameter<double>("KSgoodness");
00042   _diffgoodness = iConfig.getParameter<double>("Diffgoodness");
00043   _dirDepth = iConfig.getParameter<unsigned int>("dirDepth");
00044   _overallgoodness = iConfig.getParameter<double>("OverAllgoodness");
00045   
00046   _dbe = edm::Service<DQMStore>().operator->();
00047 
00048 
00049   /*
00050     produces<MEtoEDM<TH1F>, edm::InLumi>("name");
00051 
00052     product<TH1S,edm::InLumi>();
00053     product<TH1F,edm::InLumi>();
00054     product<TH1D,edm::InLumi>();
00055 
00056     product<TH1S,edm::InRun>();
00057     product<TH1F,edm::InRun>();
00058     product<TH1D,edm::InRun>();
00059   */
00060 }
00061 
00062 template<class T,class W>
00063 void
00064 MEtoMEComparitor::product(){
00065   /*
00066     typedef typename MEtoEDM<T> prod;
00067     produces<prod,W>("ref");
00068     produces<prod,W>("new");
00069     produces<prod,W>("diff");
00070   */
00071 }
00072 
00073 MEtoMEComparitor::~MEtoMEComparitor()
00074 {
00075 
00076 }
00077 
00078 void 
00079 MEtoMEComparitor::endLuminosityBlock(const edm::LuminosityBlock& iLumi, const edm::EventSetup&)
00080 {
00081 
00082   compare<edm::LuminosityBlock,TH1S>(iLumi,_lumiInstance);
00083   compare<edm::LuminosityBlock,TH1F>(iLumi,_lumiInstance);
00084   compare<edm::LuminosityBlock,TH1D>(iLumi,_lumiInstance);
00085 }
00086 
00087 void
00088 MEtoMEComparitor::endRun(const edm::Run& iRun, const edm::EventSetup& iSetup)
00089 {
00090 
00091   compare<edm::Run,TH1S>(iRun,_runInstance);
00092   compare<edm::Run,TH1F>(iRun,_runInstance);
00093   compare<edm::Run,TH1D>(iRun,_runInstance);
00094 
00095 }
00096 
00097 
00098 template <class T>
00099 void MEtoMEComparitor::book(const std::string & directory,const std::string & type, const T * h){
00100   _dbe->setCurrentFolder(type);
00101   std::type_info const & tp = typeid(*h);
00102   if (tp == typeid(TH1S))
00103     _dbe->book1S(h->GetName(),dynamic_cast<TH1S*>(const_cast<T*>(h)));
00104   else if (tp == typeid(TH1F))
00105     _dbe->book1D(h->GetName(),dynamic_cast<TH1F*>(const_cast<T*>(h)));
00106   else if (tp == typeid(TH1D))
00107     _dbe->book1DD(h->GetName(),dynamic_cast<TH1D*>(const_cast<T*>(h)));
00108 }
00109 template <class T> 
00110 void MEtoMEComparitor::keepBadHistograms(const std::string & directory, const T * h_new, const T * h_ref){
00111   //put it in a collection rather.
00112 
00113   
00114   std::string d_n(h_new->GetName());
00115   d_n+="_diff";
00116   T * difference = new T(d_n.c_str(),
00117                          h_new->GetTitle(),
00118                          h_new->GetNbinsX(),
00119                          h_new->GetXaxis()->GetXmin(),
00120                          h_new->GetXaxis()->GetXmax());
00121   difference->Add(h_new);
00122   difference->Add(h_ref,-1.);
00123 
00124   book(directory,"Ref",h_ref);
00125   book(directory,"New",h_new);
00126   book(directory,"Diff",difference);
00127   delete difference;
00128 
00129 }
00130 
00131 
00132 template <class W,
00133           //class Wto,
00134           class T>
00135 void MEtoMEComparitor::compare(const W& where,const std::string & instance){
00136 
00137   edm::Handle<MEtoEDM<T> > metoedm_ref;
00138   edm::Handle<MEtoEDM<T> > metoedm_new;
00139   where.getByLabel(edm::InputTag(_moduleLabel,
00140                                  instance,
00141                                  _process_ref),
00142                    metoedm_ref);
00143   where.getByLabel(edm::InputTag(_moduleLabel,
00144                                  instance,
00145                                  _process_new),
00146                    metoedm_new);
00147 
00148   if (metoedm_ref.failedToGet() || metoedm_new.failedToGet()){
00149     edm::LogError("ProductNotFound")<<"MEtoMEComparitor did not find his products.";
00150     return;
00151   }
00152 
00153   typedef typename MEtoEDM<T>::MEtoEDMObject MEtoEDMObject; 
00154   
00155   const std::vector<MEtoEDMObject> & metoedmobject_ref = metoedm_ref->getMEtoEdmObject();
00156   const std::vector<MEtoEDMObject> & metoedmobject_new = metoedm_new->getMEtoEdmObject();
00157 
00158   typedef std::map<std::string, std::pair<const MEtoEDMObject*, const MEtoEDMObject*> > Mapping;
00159   typedef typename std::map<std::string, std::pair<const MEtoEDMObject*, const MEtoEDMObject*> >::iterator Mapping_iterator;
00160 
00161   Mapping mapping;
00162 
00163   LogDebug("MEtoMEComparitor")<<"going to do the mapping from "<<metoedmobject_ref.size()<<" x "<<metoedmobject_new.size();
00164   unsigned int countMe=0;
00165   for (unsigned int i_new=0; i_new!= metoedmobject_new.size(); ++i_new){
00166     const std::string & pathname = metoedmobject_new[i_new].name;
00167     if (metoedmobject_new[i_new].object.GetEntries()==0 ||
00168         metoedmobject_new[i_new].object.Integral()==0){
00169       countMe--;
00170       continue;
00171     }
00172     mapping[pathname]=std::make_pair(&metoedmobject_new[i_new],(const MEtoEDMObject*)0);
00173   }
00174   for (unsigned int i_ref=0; i_ref!= metoedmobject_ref.size() ; ++i_ref){
00175     const std::string & pathname = metoedmobject_ref[i_ref].name;
00176     Mapping_iterator there = mapping.find(pathname);
00177     if (there != mapping.end()){
00178       there->second.second = &metoedmobject_ref[i_ref];
00179     }
00180   }
00181 
00182   LogDebug("MEtoMEComparitor")<<"found "<<mapping.size()<<" pairs of plots";
00183   countMe=0;
00184 
00185   unsigned int nNoMatch=0;
00186   unsigned int nEmpty=0;
00187   unsigned int nHollow=0;
00188   unsigned int nGoodKS=0;
00189   unsigned int nBadKS=0;
00190   unsigned int nBadDiff=0;
00191   unsigned int nGoodDiff=0;
00192 
00193   typedef std::map<std::string, std::pair<unsigned int,unsigned int> > Subs;
00194   Subs subSystems;
00195 
00196   for (Mapping_iterator it = mapping.begin();
00197        it!=mapping.end(); 
00198        ++it){
00199     if (!it->second.second){
00200       //this is expected by how the map was created
00201       nNoMatch++;
00202       continue;
00203     }
00204     const T * h_ref = &it->second.second->object;
00205     const T * h_new = &it->second.first->object;
00206 
00207     lat::StringList dir = lat::StringOps::split(it->second.second->name,"/");
00208     std::string subsystem = dir[0];
00209     if (dir.size()>=_dirDepth)
00210       for (unsigned int iD=1;iD!=_dirDepth;++iD) subsystem+="/"+dir[iD];
00211     subSystems[subsystem].first++;
00212 
00213     if (h_ref->GetEntries()!=0 && h_ref->Integral()!=0){
00214       double KS=0;
00215       bool cannotComputeKS=false;
00216       try {
00217         KS = h_new->KolmogorovTest(h_ref);
00218       }
00219       catch( cms::Exception& exception ){
00220         cannotComputeKS=true;
00221       }
00222       if (KS<_KSgoodness){
00223 
00224         unsigned int total_ref=0;
00225         unsigned int absdiff=0;
00226         for (unsigned int iBin=0;
00227              iBin!=(unsigned int)h_new->GetNbinsX()+1 ;
00228              ++iBin){
00229           total_ref+=h_ref->GetBinContent(iBin);
00230           absdiff=std::abs(h_new->GetBinContent(iBin) - h_ref->GetBinContent(iBin));
00231         }
00232         double relativediff=1;
00233         if (total_ref!=0){
00234           relativediff=absdiff / (double) total_ref;
00235         }
00236         if (relativediff > _diffgoodness ){
00237           edm::LogWarning("MEtoMEComparitor")<<"for "<<h_new->GetName()
00238                                              <<" in "<<it->first    
00239                                              <<" the KS is "<<KS*100.<<" %"
00240                                              <<" and the relative diff is: "<<relativediff*100.<<" %"
00241                                              <<" KS"<<((cannotComputeKS)?" not valid":" is valid");
00242           //std::string(" KolmogorovTest is not happy on : ")+h_new->GetName() : "";
00243           //there you want to output the plots somewhere
00244           keepBadHistograms(subsystem,h_new,h_ref);
00245           
00246           nBadDiff++;
00247           subSystems[subsystem].second++;
00248         }else{
00249           nGoodDiff++;
00250         }
00251         nBadKS++;
00252       }
00253       else
00254         nGoodKS++;
00255     }
00256     else{
00257       if (h_ref->GetEntries()==0)
00258         nEmpty++;
00259       else
00260         nHollow++;
00261       LogDebug("MEtoMEComparitor")<<h_new->GetName()   <<" in "<<it->first    <<" is empty";
00262       countMe--;
00263     }
00264     
00265   }
00266   
00267   if (mapping.size()!=0){
00268     std::stringstream summary;
00269     summary<<" Summary :"
00270            <<"\n not matched : "<<nNoMatch
00271            <<"\n empty : "<<nEmpty
00272            <<"\n integral zero : "<<nHollow
00273            <<"\n good KS : "<<nGoodKS
00274            <<"\n bad KS : "<<nBadKS
00275            <<"\n bad diff : "<<nBadDiff
00276            <<"\n godd diff : "<<nGoodDiff;
00277     bool tell=false;
00278     for (Subs::iterator iSub=subSystems.begin();
00279          iSub!=subSystems.end();++iSub){
00280       double fraction = 1-(iSub->second.second / (double)iSub->second.first);
00281       summary<<std::endl<<"Subsytem: "<<iSub->first<<" has "<< fraction*100<<" % goodness";
00282       if (fraction < _overallgoodness)
00283         tell=true;
00284     }
00285     if (tell)
00286       edm::LogWarning("MEtoMEComparitor")<<summary.str();
00287     else
00288       edm::LogInfo("MEtoMEComparitor")<<summary.str();
00289   }
00290   
00291 }
00292 
00293 
00294 // ------------ method called once each job just before starting event loop  ------------
00295 void 
00296 MEtoMEComparitor::beginJob()
00297 {
00298 }
00299 
00300 // ------------ method called once each job just after ending the event loop  ------------
00301 void 
00302 MEtoMEComparitor::endJob() {
00303 }
00304 
00305