CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/ElectroWeakAnalysis/Utilities/src/PdfSystematicsAnalyzer.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 PdfSystematicsAnalyzer: public edm::EDFilter {
00006 public:
00007       PdfSystematicsAnalyzer(const edm::ParameterSet& pset);
00008       virtual ~PdfSystematicsAnalyzer();
00009       virtual bool filter(edm::Event &, const edm::EventSetup&) override;
00010       virtual void beginJob() ;
00011       virtual void endJob() ;
00012 private:
00013       std::string selectorPath_;
00014       std::vector<edm::InputTag> pdfWeightTags_;
00015       unsigned int originalEvents_;
00016       unsigned int selectedEvents_;
00017       std::vector<int> pdfStart_;
00018       std::vector<double> weightedSelectedEvents_;
00019       std::vector<double> weighted2SelectedEvents_;
00020       std::vector<double> weightedEvents_;
00021 };
00022 
00024 #include "FWCore/Framework/interface/MakerMacros.h"
00025 #include "FWCore/Framework/interface/Frameworkfwd.h"
00026 #include "FWCore/Framework/interface/Event.h"
00027 #include "FWCore/Framework/interface/EventSetup.h"
00028 #include "FWCore/Framework/interface/ESHandle.h"
00029 #include "DataFormats/Common/interface/Handle.h"
00030 
00031 #include "FWCore/Common/interface/TriggerNames.h"
00032 #include "DataFormats/Common/interface/TriggerResults.h"
00033 
00035 PdfSystematicsAnalyzer::PdfSystematicsAnalyzer(const edm::ParameterSet& pset) :
00036   selectorPath_(pset.getUntrackedParameter<std::string> ("SelectorPath","")),
00037   pdfWeightTags_(pset.getUntrackedParameter<std::vector<edm::InputTag> > ("PdfWeightTags")) { 
00038 }
00039 
00041 PdfSystematicsAnalyzer::~PdfSystematicsAnalyzer(){}
00042 
00044 void PdfSystematicsAnalyzer::beginJob(){
00045       originalEvents_ = 0;
00046       selectedEvents_ = 0;
00047       edm::LogVerbatim("PDFAnalysis") << "PDF uncertainties will be determined for the following sets: ";
00048       for (unsigned int i=0; i<pdfWeightTags_.size(); ++i) {
00049             edm::LogVerbatim("PDFAnalysis") << "\t" << pdfWeightTags_[i].instance();
00050             pdfStart_.push_back(-1);
00051       }
00052 }
00053 
00055 void PdfSystematicsAnalyzer::endJob(){
00056 
00057       if (originalEvents_==0) {
00058             edm::LogVerbatim("PDFAnalysis") << "NO EVENTS => NO RESULTS";
00059             return;
00060       }
00061       if (selectedEvents_==0) {
00062             edm::LogVerbatim("PDFAnalysis") << "NO SELECTED EVENTS => NO RESULTS";
00063             return;
00064       }
00065 
00066       edm::LogVerbatim("PDFAnalysis") << "\n>>>> Begin of PDF weight systematics summary >>>>";
00067       edm::LogVerbatim("PDFAnalysis") << "Total number of analyzed data: " << originalEvents_ << " [events]";
00068       double originalAcceptance = double(selectedEvents_)/originalEvents_;
00069       edm::LogVerbatim("PDFAnalysis") << "Total number of selected data: " << selectedEvents_ << " [events], corresponding to acceptance: [" << originalAcceptance*100 << " +- " << 100*sqrt( originalAcceptance*(1.-originalAcceptance)/originalEvents_) << "] %";
00070 
00071       edm::LogVerbatim("PDFAnalysis") << "\n>>>>> PDF UNCERTAINTIES ON RATE >>>>>>";
00072       for (unsigned int i=0; i<pdfWeightTags_.size(); ++i) {
00073             bool nnpdfFlag = (pdfWeightTags_[i].instance().substr(0,5)=="NNPDF");
00074             unsigned int nmembers = weightedSelectedEvents_.size()-pdfStart_[i];
00075             if (i<pdfWeightTags_.size()-1) nmembers = pdfStart_[i+1] - pdfStart_[i];
00076             unsigned int npairs = (nmembers-1)/2;
00077             edm::LogVerbatim("PDFAnalysis") << "RATE Results for PDF set " << pdfWeightTags_[i].instance() << " ---->";
00078 
00079             double events_central = weightedSelectedEvents_[pdfStart_[i]]; 
00080             edm::LogVerbatim("PDFAnalysis") << "\tEstimate for central PDF member: " << int(events_central) << " [events]";
00081             double events2_central = weighted2SelectedEvents_[pdfStart_[i]];
00082             edm::LogVerbatim("PDFAnalysis") << "\ti.e. [" << std::setprecision(4) << 100*(events_central-selectedEvents_)/selectedEvents_ << " +- " <<
00083                 100*sqrt(events2_central-events_central+selectedEvents_*(1-originalAcceptance))/selectedEvents_ 
00084             << "] % relative variation with respect to original PDF";
00085 
00086             if (npairs>0) {
00087                   edm::LogVerbatim("PDFAnalysis") << "\tNumber of eigenvectors for uncertainty estimation: " << npairs;
00088               double wplus = 0.;
00089               double wminus = 0.;
00090               unsigned int nplus = 0;
00091               unsigned int nminus = 0;
00092               for (unsigned int j=0; j<npairs; ++j) {
00093                   double wa = weightedSelectedEvents_[pdfStart_[i]+2*j+1]/events_central-1.;
00094                   double wb = weightedSelectedEvents_[pdfStart_[i]+2*j+2]/events_central-1.; 
00095                   if (nnpdfFlag) {
00096                         if (wa>0.) {
00097                               wplus += wa*wa; 
00098                               nplus++;
00099                         } else {
00100                               wminus += wa*wa;
00101                               nminus++;
00102                         }
00103                         if (wb>0.) {
00104                               wplus += wb*wb; 
00105                               nplus++;
00106                         } else {
00107                               wminus += wb*wb;
00108                               nminus++;
00109                         }
00110                   } else {
00111                         if (wa>wb) {
00112                               if (wa<0.) wa = 0.;
00113                               if (wb>0.) wb = 0.;
00114                               wplus += wa*wa;
00115                               wminus += wb*wb;
00116                         } else {
00117                               if (wb<0.) wb = 0.;
00118                               if (wa>0.) wa = 0.;
00119                               wplus += wb*wb;
00120                               wminus += wa*wa;
00121                         }
00122                   }
00123               }
00124               if (wplus>0) wplus = sqrt(wplus);
00125               if (wminus>0) wminus = sqrt(wminus);
00126               if (nnpdfFlag) {
00127                   if (nplus>0) wplus /= sqrt(nplus);
00128                   if (nminus>0) wminus /= sqrt(nminus);
00129               }
00130               edm::LogVerbatim("PDFAnalysis") << "\tRelative uncertainty with respect to central member: +" << std::setprecision(4) << 100.*wplus << " / -" << std::setprecision(4) << 100.*wminus << " [%]";
00131             } else {
00132                   edm::LogVerbatim("PDFAnalysis") << "\tNO eigenvectors for uncertainty estimation";
00133             }
00134       }
00135 
00136       edm::LogVerbatim("PDFAnalysis") << "\n>>>>> PDF UNCERTAINTIES ON ACCEPTANCE >>>>>>";
00137       for (unsigned int i=0; i<pdfWeightTags_.size(); ++i) {
00138             bool nnpdfFlag = (pdfWeightTags_[i].instance().substr(0,5)=="NNPDF");
00139             unsigned int nmembers = weightedEvents_.size()-pdfStart_[i];
00140             if (i<pdfWeightTags_.size()-1) nmembers = pdfStart_[i+1] - pdfStart_[i];
00141             unsigned int npairs = (nmembers-1)/2;
00142             edm::LogVerbatim("PDFAnalysis") << "ACCEPTANCE Results for PDF set " << pdfWeightTags_[i].instance() << " ---->";
00143 
00144             double acc_central = 0.;
00145             double acc2_central = 0.;
00146             if (weightedEvents_[pdfStart_[i]]>0) {
00147                   acc_central = weightedSelectedEvents_[pdfStart_[i]]/weightedEvents_[pdfStart_[i]]; 
00148                   acc2_central = weighted2SelectedEvents_[pdfStart_[i]]/weightedEvents_[pdfStart_[i]]; 
00149             }
00150             double waverage = weightedEvents_[pdfStart_[i]]/originalEvents_;
00151             edm::LogVerbatim("PDFAnalysis") << "\tEstimate for central PDF member acceptance: [" << acc_central*100 << " +- " << 
00152             100*sqrt((acc2_central/waverage-acc_central*acc_central)/originalEvents_)
00153             << "] %";
00154             double xi = acc_central-originalAcceptance;
00155             double deltaxi = (acc2_central-(originalAcceptance+2*xi+xi*xi))/originalEvents_;
00156             if (deltaxi>0) deltaxi = sqrt(deltaxi); //else deltaxi = 0.;
00157             edm::LogVerbatim("PDFAnalysis") << "\ti.e. [" << std::setprecision(4) << 100*xi/originalAcceptance << " +- " << std::setprecision(4) << 100*deltaxi/originalAcceptance << "] % relative variation with respect to the original PDF";
00158 
00159             if (npairs>0) {
00160                   edm::LogVerbatim("PDFAnalysis") << "\tNumber of eigenvectors for uncertainty estimation: " << npairs;
00161               double wplus = 0.;
00162               double wminus = 0.;
00163               unsigned int nplus = 0;
00164               unsigned int nminus = 0;
00165               for (unsigned int j=0; j<npairs; ++j) {
00166                   double wa = 0.;
00167                   if (weightedEvents_[pdfStart_[i]+2*j+1]>0) wa = (weightedSelectedEvents_[pdfStart_[i]+2*j+1]/weightedEvents_[pdfStart_[i]+2*j+1])/acc_central-1.;
00168                   double wb = 0.;
00169                   if (weightedEvents_[pdfStart_[i]+2*j+2]>0) wb = (weightedSelectedEvents_[pdfStart_[i]+2*j+2]/weightedEvents_[pdfStart_[i]+2*j+2])/acc_central-1.;
00170                   if (nnpdfFlag) {
00171                         if (wa>0.) {
00172                               wplus += wa*wa; 
00173                               nplus++;
00174                         } else {
00175                               wminus += wa*wa;
00176                               nminus++;
00177                         }
00178                         if (wb>0.) {
00179                               wplus += wb*wb; 
00180                               nplus++;
00181                         } else {
00182                               wminus += wb*wb;
00183                               nminus++;
00184                         }
00185                   } else {
00186                         if (wa>wb) {
00187                               if (wa<0.) wa = 0.;
00188                               if (wb>0.) wb = 0.;
00189                               wplus += wa*wa;
00190                               wminus += wb*wb;
00191                         } else {
00192                               if (wb<0.) wb = 0.;
00193                               if (wa>0.) wa = 0.;
00194                               wplus += wb*wb;
00195                               wminus += wa*wa;
00196                         }
00197                   }
00198               }
00199               if (wplus>0) wplus = sqrt(wplus);
00200               if (wminus>0) wminus = sqrt(wminus);
00201               if (nnpdfFlag) {
00202                   if (nplus>0) wplus /= sqrt(nplus);
00203                   if (nminus>0) wminus /= sqrt(nminus);
00204               }
00205               edm::LogVerbatim("PDFAnalysis") << "\tRelative uncertainty with respect to central member: +" << std::setprecision(4) << 100.*wplus << " / -" << std::setprecision(4) << 100.*wminus << " [%]";
00206             } else {
00207                   edm::LogVerbatim("PDFAnalysis") << "\tNO eigenvectors for uncertainty estimation";
00208             }
00209       }
00210       edm::LogVerbatim("PDFAnalysis") << ">>>> End of PDF weight systematics summary >>>>";
00211 
00212 }
00213 
00215 bool PdfSystematicsAnalyzer::filter(edm::Event & ev, const edm::EventSetup&){
00216 
00217       edm::Handle<std::vector<double> > weightHandle;
00218       for (unsigned int i=0; i<pdfWeightTags_.size(); ++i) {
00219             if (!ev.getByLabel(pdfWeightTags_[i], weightHandle)) {
00220                   if (originalEvents_==0) {
00221                         edm::LogError("PDFAnalysis") << ">>> WARNING: some weights not found!";
00222                         edm::LogError("PDFAnalysis") << ">>> But maybe OK, if you are prefiltering!";
00223                         edm::LogError("PDFAnalysis") << ">>> If things are OK, this warning should disappear after a while!";
00224                   }
00225                   return false;
00226             }
00227       }
00228 
00229       originalEvents_++;
00230 
00231       bool selectedEvent = false;
00232       edm::Handle<edm::TriggerResults> triggerResults;
00233       if (!ev.getByLabel(edm::InputTag("TriggerResults"), triggerResults)) {
00234             edm::LogError("PDFAnalysis") << ">>> TRIGGER collection does not exist !!!";
00235             return false;
00236       }
00237 
00238 
00239       const edm::TriggerNames & trigNames = ev.triggerNames(*triggerResults);
00240       unsigned int pathIndex = trigNames.triggerIndex(selectorPath_);
00241       bool pathFound = (pathIndex<trigNames.size()); // pathIndex >= 0, since pathIndex is unsigned
00242       if (pathFound) {
00243             if (triggerResults->accept(pathIndex)) selectedEvent = true;
00244       }
00245       //edm::LogVerbatim("PDFAnalysis") << ">>>> Path Name: " << selectorPath_ << ", selected? " << selectedEvent;
00246 
00247       if (selectedEvent) selectedEvents_++;
00248 
00249       for (unsigned int i=0; i<pdfWeightTags_.size(); ++i) {
00250             if (!ev.getByLabel(pdfWeightTags_[i], weightHandle)) return false;
00251             std::vector<double> weights = (*weightHandle);
00252             unsigned int nmembers = weights.size();
00253             // Set up arrays the first time wieghts are read
00254             if (pdfStart_[i]<0) {
00255                   pdfStart_[i] = weightedEvents_.size();
00256                   for (unsigned int j=0; j<nmembers; ++j) {
00257                         weightedEvents_.push_back(0.);
00258                         weightedSelectedEvents_.push_back(0.);
00259                         weighted2SelectedEvents_.push_back(0.);
00260                   }
00261             }
00262             
00263             for (unsigned int j=0; j<nmembers; ++j) {
00264                   weightedEvents_[pdfStart_[i]+j] += weights[j];
00265                   if (selectedEvent) {
00266                         weightedSelectedEvents_[pdfStart_[i]+j] += weights[j];
00267                         weighted2SelectedEvents_[pdfStart_[i]+j] += weights[j]*weights[j];
00268                   }
00269             }
00270 
00271             /*
00272             printf("\n>>>>>>>>> Run %8d Event %d, members %3d PDF set %s : Weights >>>> \n", ev.id().run(), ev.id().event(), nmembers, pdfWeightTags_[i].instance().data());
00273             for (unsigned int i=0; i<nmembers; i+=5) {
00274                   for (unsigned int j=0; ((j<5)&&(i+j<nmembers)); ++j) {
00275                         printf(" %2d: %7.4f", i+j, weights[i+j]);
00276                   }
00277                   safe_printf("\n");
00278             }
00279             */
00280 
00281       }
00282 
00283       return true;
00284 }
00285 
00286 DEFINE_FWK_MODULE(PdfSystematicsAnalyzer);