CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/ElectroWeakAnalysis/Utilities/src/SimpleSystematicsAnalyzer.cc

Go to the documentation of this file.
00001 
00002 #include "FWCore/Framework/interface/EDFilter.h"
00003 #include "FWCore/Utilities/interface/InputTag.h"
00004 
00005 class SimpleSystematicsAnalyzer: public edm::EDFilter {
00006 public:
00007       SimpleSystematicsAnalyzer(const edm::ParameterSet& pset);
00008       virtual ~SimpleSystematicsAnalyzer();
00009       virtual bool filter(edm::Event &, const edm::EventSetup&);
00010       virtual void beginJob() ;
00011       virtual void endJob() ;
00012 private:
00013       std::string selectorPath_;
00014       std::vector<edm::InputTag> weightTags_;
00015       unsigned int originalEvents_;
00016       std::vector<double> weightedEvents_;
00017       unsigned int selectedEvents_;
00018       std::vector<double> weightedSelectedEvents_;
00019       std::vector<double> weighted2SelectedEvents_;
00020 };
00021 
00023 #include "FWCore/Framework/interface/MakerMacros.h"
00024 #include "FWCore/Framework/interface/Frameworkfwd.h"
00025 #include "FWCore/Framework/interface/Event.h"
00026 #include "FWCore/Framework/interface/EventSetup.h"
00027 #include "FWCore/Framework/interface/ESHandle.h"
00028 #include "DataFormats/Common/interface/Handle.h"
00029 
00030 #include "FWCore/Common/interface/TriggerNames.h"
00031 #include "DataFormats/Common/interface/TriggerResults.h"
00032 
00034 SimpleSystematicsAnalyzer::SimpleSystematicsAnalyzer(const edm::ParameterSet& pset) :
00035   selectorPath_(pset.getUntrackedParameter<std::string> ("SelectorPath","")),
00036   weightTags_(pset.getUntrackedParameter<std::vector<edm::InputTag> > ("WeightTags")) { 
00037 }
00038 
00040 SimpleSystematicsAnalyzer::~SimpleSystematicsAnalyzer(){}
00041 
00043 void SimpleSystematicsAnalyzer::beginJob(){
00044       originalEvents_ = 0;
00045       selectedEvents_ = 0;
00046       edm::LogVerbatim("SimpleSystematicsAnalysis") << "Uncertainties will be determined for the following tags: ";
00047       for (unsigned int i=0; i<weightTags_.size(); ++i) {
00048             edm::LogVerbatim("SimpleSystematicsAnalysis") << "\t" << weightTags_[i].encode();
00049             weightedEvents_.push_back(0.);
00050             weightedSelectedEvents_.push_back(0.);
00051             weighted2SelectedEvents_.push_back(0.);
00052       }
00053 }
00054 
00056 void SimpleSystematicsAnalyzer::endJob(){
00057       if (originalEvents_==0) {
00058             edm::LogVerbatim("SimpleSystematicsAnalysis") << "NO EVENTS => NO RESULTS";
00059             return;
00060       }
00061       if (selectedEvents_==0) {
00062             edm::LogVerbatim("SimpleSystematicsAnalysis") << "NO SELECTED EVENTS => NO RESULTS";
00063             return;
00064       }
00065 
00066       edm::LogVerbatim("SimpleSystematicsAnalysis") << "\n>>>> Begin of Weight systematics summary >>>>";
00067       edm::LogVerbatim("SimpleSystematicsAnalysis") << "Total number of analyzed data: " << originalEvents_ << " [events]";
00068       double originalAcceptance = double(selectedEvents_)/originalEvents_;
00069       edm::LogVerbatim("SimpleSystematicsAnalysis") << "Total number of selected data: " << selectedEvents_ << " [events], corresponding to acceptance: [" << originalAcceptance*100 << " +- " << 100*sqrt( originalAcceptance*(1.-originalAcceptance)/originalEvents_) << "] %";
00070       
00071       for (unsigned int i=0; i<weightTags_.size(); ++i) {
00072             edm::LogVerbatim("SimpleSystematicsAnalysis") << "Results for Weight Tag: " << weightTags_[i].encode() << " ---->";
00073 
00074             double acc_central = 0.;
00075             double acc2_central = 0.;
00076             if (weightedEvents_[i]>0) {
00077                   acc_central = weightedSelectedEvents_[i]/weightedEvents_[i]; 
00078                   acc2_central = weighted2SelectedEvents_[i]/weightedEvents_[i]; 
00079             }
00080             double waverage = weightedEvents_[i]/originalEvents_;
00081             edm::LogVerbatim("SimpleSystematicsAnalysis") << "\tTotal Events after reweighting: " << weightedEvents_[i] << " [events]";
00082             edm::LogVerbatim("SimpleSystematicsAnalysis") << "\tEvents selected after reweighting: " << weightedSelectedEvents_[i] << " [events]";
00083             edm::LogVerbatim("SimpleSystematicsAnalysis") << "\tAcceptance after reweighting: [" << acc_central*100 << " +- " << 
00084             100*sqrt((acc2_central/waverage-acc_central*acc_central)/originalEvents_)
00085             << "] %";
00086             double xi = acc_central-originalAcceptance;
00087             double deltaxi = (acc2_central-(originalAcceptance+2*xi+xi*xi))/originalEvents_;
00088             if (deltaxi>0) deltaxi = sqrt(deltaxi); else deltaxi = 0.;
00089             edm::LogVerbatim("SimpleSystematicsAnalysis") << "\ti.e. [" << std::setprecision(4) << 100*xi/originalAcceptance << " +- " << std::setprecision(4) << 100*deltaxi/originalAcceptance << "] % relative variation with respect to the original acceptance";
00090 
00091       }
00092       edm::LogVerbatim("SimpleSystematicsAnalysis") << ">>>> End of Weight systematics summary >>>>";
00093 
00094 }
00095 
00097 bool SimpleSystematicsAnalyzer::filter(edm::Event & ev, const edm::EventSetup&){
00098       originalEvents_++;
00099 
00100       bool selectedEvent = false;
00101       edm::Handle<edm::TriggerResults> triggerResults;
00102       if (!ev.getByLabel(edm::InputTag("TriggerResults"), triggerResults)) {
00103             edm::LogError("SimpleSystematicsAnalysis") << ">>> TRIGGER collection does not exist !!!";
00104             return false;
00105       }
00106 
00107       const edm::TriggerNames & trigNames = ev.triggerNames(*triggerResults);
00108       unsigned int pathIndex = trigNames.triggerIndex(selectorPath_);
00109       bool pathFound = (pathIndex<trigNames.size()); // pathIndex >= 0, since pathIndex is unsigned
00110       if (pathFound) {
00111             if (triggerResults->accept(pathIndex)) selectedEvent = true;
00112       }
00113       //edm::LogVerbatim("SimpleSystematicsAnalysis") << ">>>> Path Name: " << selectorPath_ << ", selected? " << selectedEvent;
00114 
00115       if (selectedEvent) selectedEvents_++;
00116 
00117       for (unsigned int i=0; i<weightTags_.size(); ++i) {
00118             edm::Handle<double> weightHandle;
00119             ev.getByLabel(weightTags_[i], weightHandle);
00120             weightedEvents_[i] += (*weightHandle);
00121             if (selectedEvent) {
00122                   weightedSelectedEvents_[i] += (*weightHandle);
00123                   weighted2SelectedEvents_[i] += pow((*weightHandle),2);
00124             }
00125       }
00126 
00127       return true;
00128 }
00129 
00130 DEFINE_FWK_MODULE(SimpleSystematicsAnalyzer);