CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/CalibTracker/SiStripLorentzAngle/plugins/EnsembleCalibrationLA.cc

Go to the documentation of this file.
00001 #include "CalibTracker/SiStripLorentzAngle/plugins/EnsembleCalibrationLA.h"
00002 #include "CalibTracker/SiStripCommon/interface/Book.h"
00003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00004 #include <TChain.h>
00005 #include <TFile.h>
00006 #include <boost/foreach.hpp>
00007 #include <boost/lexical_cast.hpp>
00008 #include <fstream>
00009 
00010 namespace sistrip {
00011 
00012 EnsembleCalibrationLA::EnsembleCalibrationLA(const edm::ParameterSet& conf) :
00013   inputFiles( conf.getParameter<std::vector<std::string> >("InputFiles") ),
00014   inFileLocation( conf.getParameter<std::string>("InFileLocation")),
00015   Prefix( conf.getUntrackedParameter<std::string>("Prefix","")),
00016   maxEvents( conf.getUntrackedParameter<unsigned>("MaxEvents",0)),
00017   samples( conf.getParameter<unsigned>("Samples")),
00018   nbins( conf.getParameter<unsigned>("NBins")),
00019   lowBin( conf.getParameter<double>("LowBin")),
00020   highBin( conf.getParameter<double>("HighBin")),
00021   vMethods( conf.getParameter<std::vector<int> >("Methods"))
00022 {}
00023 
00024 void EnsembleCalibrationLA::
00025 endJob() 
00026 {
00027   Book book("la_ensemble");
00028   TChain*const chain = new TChain("la_ensemble"); 
00029   BOOST_FOREACH(std::string file, inputFiles) chain->Add((file+inFileLocation).c_str());
00030 
00031   int methods = 0;  BOOST_FOREACH(unsigned method, vMethods) methods|=method;
00032 
00033   LA_Filler_Fitter 
00034     laff(methods,samples,nbins,lowBin,highBin,maxEvents);
00035   laff.fill(chain,book);           
00036   laff.fit(book);                  
00037   laff.summarize_ensembles(book);  
00038 
00039   write_ensembles_text(book);
00040   write_ensembles_plots(book);
00041   write_samples_plots(book);
00042   write_calibrations();
00043 }
00044 
00045 void EnsembleCalibrationLA::
00046 write_ensembles_text(const Book& book) {
00047   std::pair<std::string, std::vector<LA_Filler_Fitter::EnsembleSummary> > ensemble;
00048   BOOST_FOREACH(ensemble, LA_Filler_Fitter::ensemble_summary(book)) {
00049     fstream file((Prefix+ensemble.first+".dat").c_str(),std::ios::out);
00050     BOOST_FOREACH(LA_Filler_Fitter::EnsembleSummary summary, ensemble.second)
00051       file << summary << std::endl;
00052 
00053     const std::pair<std::pair<float,float>,std::pair<float,float> > line = LA_Filler_Fitter::offset_slope(ensemble.second);
00054     const float pull =  LA_Filler_Fitter::pull(ensemble.second);
00055 
00056     unsigned index = 15;
00057     std::string label;
00058     {
00059       std::cout << ensemble.first << std::endl;
00060       boost::regex format(".*(T[IO]B)_layer(\\d)([as])_(.*)");
00061       if(boost::regex_match(ensemble.first,format)) {
00062         const bool TIB = "TIB" == boost::regex_replace(ensemble.first, format, "\\1");
00063         const bool stereo = "s" == boost::regex_replace(ensemble.first, format, "\\3");
00064         const unsigned layer = boost::lexical_cast<unsigned>(boost::regex_replace(ensemble.first, format, "\\2"));
00065         label = boost::regex_replace(ensemble.first, format, "\\4");
00066         index = LA_Filler_Fitter::layer_index(TIB,stereo,layer);
00067 
00068         calibrations[label].slopes[index]=line.second.first;
00069         calibrations[label].offsets[index]=line.first.first;
00070         calibrations[label].pulls[index]=pull;
00071       }
00072     }
00073 
00074     file << std::endl << std::endl
00075          << "# Best Fit Line: "  
00076          << line.first.first <<"("<< line.first.second<<")   +   x* "
00077          << line.second.first<<"("<< line.second.second<<")"           
00078          << std::endl
00079          << "# Pull (average sigma of (x_measure-x_truth)/e_measure): " << pull 
00080          << std::endl
00081          << "LA_Calibration( METHOD_XXXXX , xxx, " << line.second.first << ", " << line.first.first << ", " << pull << ")," << std::endl;
00082     file.close();
00083   } 
00084 }
00085 
00086 void EnsembleCalibrationLA::
00087 write_ensembles_plots(const Book& book) const {
00088   TFile file((Prefix+"sampleFits.root").c_str(),"RECREATE");
00089   for(Book::const_iterator hist = book.begin(".*(profile|ratio|reconstruction|symm|symmchi2|_w\\d)"); hist!=book.end(); ++hist)
00090     hist->second->Write();
00091   file.Close();
00092 }
00093  
00094 void EnsembleCalibrationLA::
00095 write_samples_plots(const Book& book) const {
00096   TFile file((Prefix+"ensembleFits.root").c_str(),"RECREATE");
00097   for(Book::const_iterator hist = book.begin(".*(measure|merr|ensembleReco|pull)"); hist!=book.end(); ++hist)
00098     hist->second->Write();
00099   file.Close();
00100 }
00101 
00102 void EnsembleCalibrationLA::
00103 write_calibrations() const {
00104   fstream file((Prefix+"calibrations.dat").c_str(),std::ios::out);
00105   std::pair<std::string,MethodCalibrations> cal;
00106   BOOST_FOREACH(cal,calibrations) {
00107     file << cal.first << std::endl
00108          << "\t slopes(";    BOOST_FOREACH(float i, cal.second.slopes) file << i<< ","; file << ")" << std::endl
00109          << "\t offsets(";   BOOST_FOREACH(float i, cal.second.offsets) file << i<< ","; file << ")" << std::endl
00110          << "\t pulls(";     BOOST_FOREACH(float i, cal.second.pulls) file << i<< ","; file << ")" << std::endl;
00111   }
00112   file.close();
00113 }
00114   
00115 }
00116