CMS 3D CMS Logo

Public Member Functions | Private Attributes

PdfSystematicsAnalyzer Class Reference

Inheritance diagram for PdfSystematicsAnalyzer:
edm::EDFilter edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

virtual void beginJob ()
virtual void endJob ()
virtual bool filter (edm::Event &, const edm::EventSetup &)
 PdfSystematicsAnalyzer (const edm::ParameterSet &pset)
virtual ~PdfSystematicsAnalyzer ()

Private Attributes

unsigned int originalEvents_
std::vector< int > pdfStart_
std::vector< edm::InputTagpdfWeightTags_
unsigned int selectedEvents_
std::string selectorPath_
std::vector< double > weighted2SelectedEvents_
std::vector< double > weightedEvents_
std::vector< double > weightedSelectedEvents_

Detailed Description

Definition at line 5 of file PdfSystematicsAnalyzer.cc.


Constructor & Destructor Documentation

PdfSystematicsAnalyzer::PdfSystematicsAnalyzer ( const edm::ParameterSet pset)

Definition at line 35 of file PdfSystematicsAnalyzer.cc.

                                                                          :
  selectorPath_(pset.getUntrackedParameter<std::string> ("SelectorPath","")),
  pdfWeightTags_(pset.getUntrackedParameter<std::vector<edm::InputTag> > ("PdfWeightTags")) { 
}
PdfSystematicsAnalyzer::~PdfSystematicsAnalyzer ( ) [virtual]

Definition at line 41 of file PdfSystematicsAnalyzer.cc.

{}

Member Function Documentation

void PdfSystematicsAnalyzer::beginJob ( void  ) [virtual]

Reimplemented from edm::EDFilter.

Definition at line 44 of file PdfSystematicsAnalyzer.cc.

References i, originalEvents_, pdfStart_, pdfWeightTags_, and selectedEvents_.

                                     {
      originalEvents_ = 0;
      selectedEvents_ = 0;
      edm::LogVerbatim("PDFAnalysis") << "PDF uncertainties will be determined for the following sets: ";
      for (unsigned int i=0; i<pdfWeightTags_.size(); ++i) {
            edm::LogVerbatim("PDFAnalysis") << "\t" << pdfWeightTags_[i].instance();
            pdfStart_.push_back(-1);
      }
}
void PdfSystematicsAnalyzer::endJob ( void  ) [virtual]

Reimplemented from edm::EDFilter.

Definition at line 55 of file PdfSystematicsAnalyzer.cc.

References i, j, originalEvents_, pdfStart_, pdfWeightTags_, selectedEvents_, mathSSE::sqrt(), weighted2SelectedEvents_, weightedEvents_, and weightedSelectedEvents_.

                                   {

      if (originalEvents_==0) {
            edm::LogVerbatim("PDFAnalysis") << "NO EVENTS => NO RESULTS";
            return;
      }
      if (selectedEvents_==0) {
            edm::LogVerbatim("PDFAnalysis") << "NO SELECTED EVENTS => NO RESULTS";
            return;
      }

      edm::LogVerbatim("PDFAnalysis") << "\n>>>> Begin of PDF weight systematics summary >>>>";
      edm::LogVerbatim("PDFAnalysis") << "Total number of analyzed data: " << originalEvents_ << " [events]";
      double originalAcceptance = double(selectedEvents_)/originalEvents_;
      edm::LogVerbatim("PDFAnalysis") << "Total number of selected data: " << selectedEvents_ << " [events], corresponding to acceptance: [" << originalAcceptance*100 << " +- " << 100*sqrt( originalAcceptance*(1.-originalAcceptance)/originalEvents_) << "] %";

      edm::LogVerbatim("PDFAnalysis") << "\n>>>>> PDF UNCERTAINTIES ON RATE >>>>>>";
      for (unsigned int i=0; i<pdfWeightTags_.size(); ++i) {
            bool nnpdfFlag = (pdfWeightTags_[i].instance().substr(0,5)=="NNPDF");
            unsigned int nmembers = weightedSelectedEvents_.size()-pdfStart_[i];
            if (i<pdfWeightTags_.size()-1) nmembers = pdfStart_[i+1] - pdfStart_[i];
            unsigned int npairs = (nmembers-1)/2;
            edm::LogVerbatim("PDFAnalysis") << "RATE Results for PDF set " << pdfWeightTags_[i].instance() << " ---->";

            double events_central = weightedSelectedEvents_[pdfStart_[i]]; 
            edm::LogVerbatim("PDFAnalysis") << "\tEstimate for central PDF member: " << int(events_central) << " [events]";
            double events2_central = weighted2SelectedEvents_[pdfStart_[i]];
            edm::LogVerbatim("PDFAnalysis") << "\ti.e. [" << std::setprecision(4) << 100*(events_central-selectedEvents_)/selectedEvents_ << " +- " <<
                100*sqrt(events2_central-events_central+selectedEvents_*(1-originalAcceptance))/selectedEvents_ 
            << "] % relative variation with respect to original PDF";

            if (npairs>0) {
                  edm::LogVerbatim("PDFAnalysis") << "\tNumber of eigenvectors for uncertainty estimation: " << npairs;
              double wplus = 0.;
              double wminus = 0.;
              unsigned int nplus = 0;
              unsigned int nminus = 0;
              for (unsigned int j=0; j<npairs; ++j) {
                  double wa = weightedSelectedEvents_[pdfStart_[i]+2*j+1]/events_central-1.;
                  double wb = weightedSelectedEvents_[pdfStart_[i]+2*j+2]/events_central-1.; 
                  if (nnpdfFlag) {
                        if (wa>0.) {
                              wplus += wa*wa; 
                              nplus++;
                        } else {
                              wminus += wa*wa;
                              nminus++;
                        }
                        if (wb>0.) {
                              wplus += wb*wb; 
                              nplus++;
                        } else {
                              wminus += wb*wb;
                              nminus++;
                        }
                  } else {
                        if (wa>wb) {
                              if (wa<0.) wa = 0.;
                              if (wb>0.) wb = 0.;
                              wplus += wa*wa;
                              wminus += wb*wb;
                        } else {
                              if (wb<0.) wb = 0.;
                              if (wa>0.) wa = 0.;
                              wplus += wb*wb;
                              wminus += wa*wa;
                        }
                  }
              }
              if (wplus>0) wplus = sqrt(wplus);
              if (wminus>0) wminus = sqrt(wminus);
              if (nnpdfFlag) {
                  if (nplus>0) wplus /= sqrt(nplus);
                  if (nminus>0) wminus /= sqrt(nminus);
              }
              edm::LogVerbatim("PDFAnalysis") << "\tRelative uncertainty with respect to central member: +" << std::setprecision(4) << 100.*wplus << " / -" << std::setprecision(4) << 100.*wminus << " [%]";
            } else {
                  edm::LogVerbatim("PDFAnalysis") << "\tNO eigenvectors for uncertainty estimation";
            }
      }

      edm::LogVerbatim("PDFAnalysis") << "\n>>>>> PDF UNCERTAINTIES ON ACCEPTANCE >>>>>>";
      for (unsigned int i=0; i<pdfWeightTags_.size(); ++i) {
            bool nnpdfFlag = (pdfWeightTags_[i].instance().substr(0,5)=="NNPDF");
            unsigned int nmembers = weightedEvents_.size()-pdfStart_[i];
            if (i<pdfWeightTags_.size()-1) nmembers = pdfStart_[i+1] - pdfStart_[i];
            unsigned int npairs = (nmembers-1)/2;
            edm::LogVerbatim("PDFAnalysis") << "ACCEPTANCE Results for PDF set " << pdfWeightTags_[i].instance() << " ---->";

            double acc_central = 0.;
            double acc2_central = 0.;
            if (weightedEvents_[pdfStart_[i]]>0) {
                  acc_central = weightedSelectedEvents_[pdfStart_[i]]/weightedEvents_[pdfStart_[i]]; 
                  acc2_central = weighted2SelectedEvents_[pdfStart_[i]]/weightedEvents_[pdfStart_[i]]; 
            }
            double waverage = weightedEvents_[pdfStart_[i]]/originalEvents_;
            edm::LogVerbatim("PDFAnalysis") << "\tEstimate for central PDF member acceptance: [" << acc_central*100 << " +- " << 
            100*sqrt((acc2_central/waverage-acc_central*acc_central)/originalEvents_)
            << "] %";
            double xi = acc_central-originalAcceptance;
            double deltaxi = (acc2_central-(originalAcceptance+2*xi+xi*xi))/originalEvents_;
            if (deltaxi>0) deltaxi = sqrt(deltaxi); //else deltaxi = 0.;
            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";

            if (npairs>0) {
                  edm::LogVerbatim("PDFAnalysis") << "\tNumber of eigenvectors for uncertainty estimation: " << npairs;
              double wplus = 0.;
              double wminus = 0.;
              unsigned int nplus = 0;
              unsigned int nminus = 0;
              for (unsigned int j=0; j<npairs; ++j) {
                  double wa = 0.;
                  if (weightedEvents_[pdfStart_[i]+2*j+1]>0) wa = (weightedSelectedEvents_[pdfStart_[i]+2*j+1]/weightedEvents_[pdfStart_[i]+2*j+1])/acc_central-1.;
                  double wb = 0.;
                  if (weightedEvents_[pdfStart_[i]+2*j+2]>0) wb = (weightedSelectedEvents_[pdfStart_[i]+2*j+2]/weightedEvents_[pdfStart_[i]+2*j+2])/acc_central-1.;
                  if (nnpdfFlag) {
                        if (wa>0.) {
                              wplus += wa*wa; 
                              nplus++;
                        } else {
                              wminus += wa*wa;
                              nminus++;
                        }
                        if (wb>0.) {
                              wplus += wb*wb; 
                              nplus++;
                        } else {
                              wminus += wb*wb;
                              nminus++;
                        }
                  } else {
                        if (wa>wb) {
                              if (wa<0.) wa = 0.;
                              if (wb>0.) wb = 0.;
                              wplus += wa*wa;
                              wminus += wb*wb;
                        } else {
                              if (wb<0.) wb = 0.;
                              if (wa>0.) wa = 0.;
                              wplus += wb*wb;
                              wminus += wa*wa;
                        }
                  }
              }
              if (wplus>0) wplus = sqrt(wplus);
              if (wminus>0) wminus = sqrt(wminus);
              if (nnpdfFlag) {
                  if (nplus>0) wplus /= sqrt(nplus);
                  if (nminus>0) wminus /= sqrt(nminus);
              }
              edm::LogVerbatim("PDFAnalysis") << "\tRelative uncertainty with respect to central member: +" << std::setprecision(4) << 100.*wplus << " / -" << std::setprecision(4) << 100.*wminus << " [%]";
            } else {
                  edm::LogVerbatim("PDFAnalysis") << "\tNO eigenvectors for uncertainty estimation";
            }
      }
      edm::LogVerbatim("PDFAnalysis") << ">>>> End of PDF weight systematics summary >>>>";

}
bool PdfSystematicsAnalyzer::filter ( edm::Event ev,
const edm::EventSetup  
) [virtual]

Implements edm::EDFilter.

Definition at line 215 of file PdfSystematicsAnalyzer.cc.

References edm::Event::getByLabel(), i, j, originalEvents_, pdfStart_, pdfWeightTags_, selectedEvents_, selectorPath_, edm::TriggerNames::size(), edm::TriggerNames::triggerIndex(), edm::Event::triggerNames(), patRefSel_triggerSelection_cff::triggerResults, trigNames, weighted2SelectedEvents_, weightedEvents_, and weightedSelectedEvents_.

                                                                      {

      edm::Handle<std::vector<double> > weightHandle;
      for (unsigned int i=0; i<pdfWeightTags_.size(); ++i) {
            if (!ev.getByLabel(pdfWeightTags_[i], weightHandle)) {
                  if (originalEvents_==0) {
                        edm::LogError("PDFAnalysis") << ">>> WARNING: some weights not found!";
                        edm::LogError("PDFAnalysis") << ">>> But maybe OK, if you are prefiltering!";
                        edm::LogError("PDFAnalysis") << ">>> If things are OK, this warning should disappear after a while!";
                  }
                  return false;
            }
      }

      originalEvents_++;

      bool selectedEvent = false;
      edm::Handle<edm::TriggerResults> triggerResults;
      if (!ev.getByLabel(edm::InputTag("TriggerResults"), triggerResults)) {
            edm::LogError("PDFAnalysis") << ">>> TRIGGER collection does not exist !!!";
            return false;
      }


      const edm::TriggerNames & trigNames = ev.triggerNames(*triggerResults);
      unsigned int pathIndex = trigNames.triggerIndex(selectorPath_);
      bool pathFound = (pathIndex<trigNames.size()); // pathIndex >= 0, since pathIndex is unsigned
      if (pathFound) {
            if (triggerResults->accept(pathIndex)) selectedEvent = true;
      }
      //edm::LogVerbatim("PDFAnalysis") << ">>>> Path Name: " << selectorPath_ << ", selected? " << selectedEvent;

      if (selectedEvent) selectedEvents_++;

      for (unsigned int i=0; i<pdfWeightTags_.size(); ++i) {
            if (!ev.getByLabel(pdfWeightTags_[i], weightHandle)) return false;
            std::vector<double> weights = (*weightHandle);
            unsigned int nmembers = weights.size();
            // Set up arrays the first time wieghts are read
            if (pdfStart_[i]<0) {
                  pdfStart_[i] = weightedEvents_.size();
                  for (unsigned int j=0; j<nmembers; ++j) {
                        weightedEvents_.push_back(0.);
                        weightedSelectedEvents_.push_back(0.);
                        weighted2SelectedEvents_.push_back(0.);
                  }
            }
            
            for (unsigned int j=0; j<nmembers; ++j) {
                  weightedEvents_[pdfStart_[i]+j] += weights[j];
                  if (selectedEvent) {
                        weightedSelectedEvents_[pdfStart_[i]+j] += weights[j];
                        weighted2SelectedEvents_[pdfStart_[i]+j] += weights[j]*weights[j];
                  }
            }

            /*
            printf("\n>>>>>>>>> Run %8d Event %d, members %3d PDF set %s : Weights >>>> \n", ev.id().run(), ev.id().event(), nmembers, pdfWeightTags_[i].instance().data());
            for (unsigned int i=0; i<nmembers; i+=5) {
                  for (unsigned int j=0; ((j<5)&&(i+j<nmembers)); ++j) {
                        printf(" %2d: %7.4f", i+j, weights[i+j]);
                  }
                  safe_printf("\n");
            }
            */

      }

      return true;
}

Member Data Documentation

Definition at line 15 of file PdfSystematicsAnalyzer.cc.

Referenced by beginJob(), endJob(), and filter().

std::vector<int> PdfSystematicsAnalyzer::pdfStart_ [private]

Definition at line 17 of file PdfSystematicsAnalyzer.cc.

Referenced by beginJob(), endJob(), and filter().

Definition at line 14 of file PdfSystematicsAnalyzer.cc.

Referenced by beginJob(), endJob(), and filter().

Definition at line 16 of file PdfSystematicsAnalyzer.cc.

Referenced by beginJob(), endJob(), and filter().

Definition at line 13 of file PdfSystematicsAnalyzer.cc.

Referenced by filter().

std::vector<double> PdfSystematicsAnalyzer::weighted2SelectedEvents_ [private]

Definition at line 19 of file PdfSystematicsAnalyzer.cc.

Referenced by endJob(), and filter().

std::vector<double> PdfSystematicsAnalyzer::weightedEvents_ [private]

Definition at line 20 of file PdfSystematicsAnalyzer.cc.

Referenced by endJob(), and filter().

std::vector<double> PdfSystematicsAnalyzer::weightedSelectedEvents_ [private]

Definition at line 18 of file PdfSystematicsAnalyzer.cc.

Referenced by endJob(), and filter().