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