CMS 3D CMS Logo

PostProcessor Class Reference

#include <Validation/Tools/interface/PostProcessor.h>

Inheritance diagram for PostProcessor:

edm::EDAnalyzer

List of all members.

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup)
void computeEfficiency (const std::string &startDir, const std::string &efficMEName, const std::string &efficMETitle, const std::string &recoMEName, const std::string &simMEName, const std::string &type="eff")
void computeResolution (const std::string &, const std::string &fitMEPrefix, const std::string &fitMETitlePrefix, const std::string &srcMEName)
void endJob ()
void limitedFit (MonitorElement *srcME, MonitorElement *meanME, MonitorElement *sigmaME)
 PostProcessor (const edm::ParameterSet &pset)
 ~PostProcessor ()

Private Attributes

std::vector< std::string > effCmds_
bool isWildcardUsed_
std::string outputFileName_
std::vector< std::string > resCmds_
bool resLimitedFit_
std::vector< std::string > subDirs_
DQMStoretheDQM
unsigned int verbose_


Detailed Description

Definition at line 27 of file PostProcessor.h.


Constructor & Destructor Documentation

PostProcessor::PostProcessor ( const edm::ParameterSet pset  ) 

Definition at line 27 of file PostProcessor.cc.

References effCmds_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), isWildcardUsed_, outputFileName_, resCmds_, resLimitedFit_, subDirs_, and verbose_.

00028 {
00029   verbose_ = pset.getUntrackedParameter<unsigned int>("verbose", 0);
00030 
00031   effCmds_ = pset.getParameter<vstring>("efficiency");
00032   resCmds_ = pset.getParameter<vstring>("resolution");
00033 
00034   resLimitedFit_ = pset.getUntrackedParameter<bool>("resolutionLimitedFit",false);
00035 
00036   outputFileName_ = pset.getUntrackedParameter<string>("outputFileName", "");
00037   subDirs_ = pset.getUntrackedParameter<vstring>("subDirs");
00038 
00039   isWildcardUsed_ = false;
00040 }

PostProcessor::~PostProcessor (  )  [inline]

Definition at line 31 of file PostProcessor.h.

00031 {};


Member Function Documentation

void PostProcessor::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
) [inline, virtual]

Implements edm::EDAnalyzer.

Definition at line 33 of file PostProcessor.h.

00033 {};

void PostProcessor::computeEfficiency ( const std::string &  startDir,
const std::string &  efficMEName,
const std::string &  efficMETitle,
const std::string &  recoMEName,
const std::string &  simMEName,
const std::string &  type = "eff" 
)

Referenced by endJob().

void PostProcessor::computeResolution ( const std::string &  ,
const std::string &  fitMEPrefix,
const std::string &  fitMETitlePrefix,
const std::string &  srcMEName 
)

Referenced by endJob().

void PostProcessor::endJob ( void   )  [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 42 of file PostProcessor.cc.

References DQMStore::cd(), computeEfficiency(), computeResolution(), test_cfg::dirName, effCmds_, DQMStore::getSubdirs(), isWildcardUsed_, outputFileName_, resCmds_, DQMStore::save(), DQMStore::showDirStructure(), subDirs_, theDQM, cmsScimarkStop::tokens, and verbose_.

00043 {
00044   theDQM = 0;
00045   theDQM = Service<DQMStore>().operator->();
00046 
00047   if ( ! theDQM ) {
00048     LogError("PostProcessor") << "Cannot create DQMStore instance\n";
00049     return;
00050   }
00051 
00052   // Process wildcard in the sub-directory
00053   set<string> subDirSet;
00054 
00055   for(vstring::const_iterator iSubDir = subDirs_.begin();
00056       iSubDir != subDirs_.end(); ++iSubDir) {
00057     string subDir = *iSubDir;
00058 
00059     if ( subDir[subDir.size()-1] == '/' ) subDir.erase(subDir.size()-1);
00060 
00061     if ( subDir[subDir.size()-1] == '*' ) {
00062       isWildcardUsed_ = true;
00063       const string::size_type shiftPos = subDir.rfind('/');
00064 
00065       const string searchPath = subDir.substr(0, shiftPos);
00066       theDQM->cd(searchPath);
00067 
00068       vector<string> foundDirs = theDQM->getSubdirs();
00069       const string matchStr = subDir.substr(0, subDir.size()-2);
00070 
00071       for(vector<string>::const_iterator iDir = foundDirs.begin();
00072           iDir != foundDirs.end(); ++iDir) {
00073         const string dirPrefix = iDir->substr(0, matchStr.size());
00074 
00075         if ( dirPrefix == matchStr ) {
00076           subDirSet.insert(*iDir);
00077         }
00078       }
00079     }
00080     else {
00081       subDirSet.insert(subDir);
00082     }
00083   }
00084 
00085   for(set<string>::const_iterator iSubDir = subDirSet.begin();
00086       iSubDir != subDirSet.end(); ++iSubDir) {
00087     typedef boost::escaped_list_separator<char> elsc;
00088 
00089     const string& dirName = *iSubDir;
00090 
00091     for(vstring::const_iterator iCmd = effCmds_.begin();
00092         iCmd != effCmds_.end(); ++iCmd) {
00093       if ( iCmd->empty() ) continue;
00094       boost::tokenizer<elsc> tokens(*iCmd, elsc("\\", " \t", "\'"));
00095 
00096       vector<string> args;
00097       for(boost::tokenizer<elsc>::const_iterator iToken = tokens.begin();
00098           iToken != tokens.end(); ++iToken) {
00099         if ( iToken->empty() ) continue;
00100         args.push_back(*iToken);
00101       }
00102 
00103       if ( args.size() < 4 ) {
00104         LogError("PostProcessor") << "Wrong input to effCmds\n";
00105         continue;
00106       }
00107 
00108       if ( args.size() == 4 ) args.push_back("eff");
00109 
00110       computeEfficiency(dirName, args[0], args[1], args[2], args[3], args[4]);
00111     }
00112 
00113     for(vstring::const_iterator iCmd = resCmds_.begin();
00114         iCmd != resCmds_.end(); ++ iCmd) {
00115       if ( iCmd->empty() ) continue;
00116       boost::tokenizer<elsc> tokens(*iCmd, elsc("\\", " \t", "\'"));
00117 
00118       vector<string> args;
00119       for(boost::tokenizer<elsc>::const_iterator iToken = tokens.begin();
00120           iToken != tokens.end(); ++iToken) {
00121         if ( iToken->empty() ) continue;
00122         args.push_back(*iToken);
00123       }
00124 
00125       if ( args.size() != 3 ) {
00126         LogError("PostProcessor") << "Wrong input to resCmds\n";
00127         continue;
00128       }
00129 
00130       computeResolution(dirName, args[0], args[1], args[2]);
00131     }
00132   }
00133 
00134   if ( verbose_ > 0 ) theDQM->showDirStructure();
00135 
00136   if ( ! outputFileName_.empty() ) theDQM->save(outputFileName_);
00137 }

void PostProcessor::limitedFit ( MonitorElement srcME,
MonitorElement meanME,
MonitorElement sigmaME 
)

Definition at line 296 of file PostProcessor.cc.

References cont, MonitorElement::getTH2F(), histo, i, MonitorElement::setBinContent(), and MonitorElement::setBinError().

00297 {
00298   TH2F * histo = srcME->getTH2F();
00299 
00300   static int i = 0;
00301   i++;
00302 
00303   // Fit slices projected along Y from bins in X 
00304   double cont_min = 100;    //Minimum number of entries
00305   Int_t binx =  histo->GetXaxis()->GetNbins();
00306 
00307   for (int i = 1; i <= binx ; i++) {
00308     TString iString(i);
00309     TH1 *histoY =  histo->ProjectionY(" ", i, i);
00310     double cont = histoY->GetEntries();
00311 
00312     if (cont >= cont_min) {
00313       float minfit = histoY->GetMean() - histoY->GetRMS();
00314       float maxfit = histoY->GetMean() + histoY->GetRMS();
00315       
00316       TF1 *fitFcn = new TF1(TString("g")+histo->GetName()+iString,"gaus",minfit,maxfit);
00317       double x1,x2;
00318       fitFcn->GetRange(x1,x2);
00319 
00320       histoY->Fit(fitFcn,"QR0","",x1,x2);
00321 
00322 //      histoY->Fit(fitFcn->GetName(),"RME");
00323       double *par = fitFcn->GetParameters();
00324       double *err = fitFcn->GetParErrors();
00325 
00326       meanME->setBinContent(i,par[1]);
00327       meanME->setBinError(i,err[1]);
00328 
00329       sigmaME->setBinContent(i,par[2]);
00330       sigmaME->setBinError(i,err[2]);
00331       if(fitFcn) delete fitFcn;
00332       if(histoY) delete histoY;
00333     }
00334     else {
00335       if(histoY) delete histoY;
00336       continue;
00337     }
00338   }
00339 }


Member Data Documentation

std::vector<std::string> PostProcessor::effCmds_ [private]

Definition at line 52 of file PostProcessor.h.

Referenced by endJob(), and PostProcessor().

bool PostProcessor::isWildcardUsed_ [private]

Definition at line 47 of file PostProcessor.h.

Referenced by endJob(), and PostProcessor().

std::string PostProcessor::outputFileName_ [private]

Definition at line 51 of file PostProcessor.h.

Referenced by endJob(), and PostProcessor().

std::vector<std::string> PostProcessor::resCmds_ [private]

Definition at line 52 of file PostProcessor.h.

Referenced by endJob(), and PostProcessor().

bool PostProcessor::resLimitedFit_ [private]

Definition at line 53 of file PostProcessor.h.

Referenced by PostProcessor().

std::vector<std::string> PostProcessor::subDirs_ [private]

Definition at line 50 of file PostProcessor.h.

Referenced by endJob(), and PostProcessor().

DQMStore* PostProcessor::theDQM [private]

Definition at line 49 of file PostProcessor.h.

Referenced by endJob().

unsigned int PostProcessor::verbose_ [private]

Definition at line 46 of file PostProcessor.h.

Referenced by endJob(), and PostProcessor().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:30:15 2009 for CMSSW by  doxygen 1.5.4