CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/CondFormats/JetMETObjects/src/JetCorrectorParameters.cc

Go to the documentation of this file.
00001 //
00002 // Original Author:  Fedor Ratnikov Nov 9, 2007
00003 // $Id: JetCorrectorParameters.cc,v 1.19 2011/01/27 12:14:13 kkousour Exp $
00004 //
00005 // Generic parameters for Jet corrections
00006 //
00007 #include "CondFormats/JetMETObjects/interface/JetCorrectorParameters.h"
00008 #include "CondFormats/JetMETObjects/src/Utilities.cc"
00009 #include <iostream>
00010 #include <iomanip>
00011 #include <fstream>
00012 #include <sstream>
00013 #include <stdlib.h>
00014 #include <algorithm>
00015 #include <cmath>
00016 #include <iterator>
00017 
00018 //------------------------------------------------------------------------ 
00019 //--- JetCorrectorParameters::Definitions constructor --------------------
00020 //--- takes specific arguments for the member variables ------------------
00021 //------------------------------------------------------------------------
00022 JetCorrectorParameters::Definitions::Definitions(const std::vector<std::string>& fBinVar, const std::vector<std::string>& fParVar, const std::string& fFormula, bool fIsResponse)
00023 {
00024   for(unsigned i=0;i<fBinVar.size();i++)
00025     mBinVar.push_back(fBinVar[i]);
00026   for(unsigned i=0;i<fParVar.size();i++)
00027     mParVar.push_back(fParVar[i]);
00028   mFormula    = fFormula;
00029   mIsResponse = fIsResponse;
00030   mLevel      = "";
00031 }
00032 //------------------------------------------------------------------------
00033 //--- JetCorrectorParameters::Definitions constructor --------------------
00034 //--- reads the member variables from a string ---------------------------
00035 //------------------------------------------------------------------------
00036 JetCorrectorParameters::Definitions::Definitions(const std::string& fLine)
00037 {
00038   std::vector<std::string> tokens = getTokens(fLine); 
00039   if (!tokens.empty())
00040     { 
00041       if (tokens.size() < 6) 
00042         {
00043           std::stringstream sserr;
00044           sserr<<"(line "<<fLine<<"): less than 6 expected tokens:"<<tokens.size();
00045           handleError("JetCorrectorParameters::Definitions",sserr.str());
00046         }
00047       unsigned nvar = getUnsigned(tokens[0]);
00048       unsigned npar = getUnsigned(tokens[nvar+1]);
00049       for(unsigned i=0;i<nvar;i++)
00050         mBinVar.push_back(tokens[i+1]);
00051       for(unsigned i=0;i<npar;i++)
00052         mParVar.push_back(tokens[nvar+2+i]);
00053       mFormula = tokens[npar+nvar+2];
00054       std::string ss = tokens[npar+nvar+3]; 
00055       if (ss == "Response")
00056         mIsResponse = true;
00057       else if (ss == "Correction")
00058         mIsResponse = false;
00059       else if (ss == "Resolution")
00060         mIsResponse = false;
00061       else if (ss.find("PAR")==0)
00062         mIsResponse = false;
00063       else
00064         {
00065           std::stringstream sserr;
00066           sserr<<"unknown option ("<<ss<<")"; 
00067           handleError("JetCorrectorParameters::Definitions",sserr.str());
00068         }
00069       mLevel = tokens[npar+nvar+4]; 
00070     }
00071 }
00072 //------------------------------------------------------------------------
00073 //--- JetCorrectorParameters::Record constructor -------------------------
00074 //--- reads the member variables from a string ---------------------------
00075 //------------------------------------------------------------------------
00076 JetCorrectorParameters::Record::Record(const std::string& fLine,unsigned fNvar) : mMin(0),mMax(0)
00077 {
00078   mNvar = fNvar;
00079   // quckly parse the line
00080   std::vector<std::string> tokens = getTokens(fLine);
00081   if (!tokens.empty())
00082     { 
00083       if (tokens.size() < 3) 
00084         {
00085           std::stringstream sserr;
00086           sserr<<"(line "<<fLine<<"): "<<"three tokens expected, "<<tokens.size()<<" provided.";
00087           handleError("JetCorrectorParameters::Record",sserr.str());
00088         }
00089       for(unsigned i=0;i<mNvar;i++)
00090         {
00091           mMin.push_back(getFloat(tokens[i*mNvar]));
00092           mMax.push_back(getFloat(tokens[i*mNvar+1])); 
00093         }
00094       unsigned nParam = getUnsigned(tokens[2*mNvar]);
00095       if (nParam != tokens.size()-(2*mNvar+1)) 
00096         {
00097           std::stringstream sserr;
00098           sserr<<"(line "<<fLine<<"): "<<tokens.size()-(2*mNvar+1)<<" parameters, but nParam="<<nParam<<".";
00099           handleError("JetCorrectorParameters::Record",sserr.str());
00100         }
00101       for (unsigned i = (2*mNvar+1); i < tokens.size(); ++i)
00102         mParameters.push_back(getFloat(tokens[i]));
00103     } 
00104 }
00105 //------------------------------------------------------------------------
00106 //--- JetCorrectorParameters constructor ---------------------------------
00107 //--- reads the member variables from a string ---------------------------
00108 //------------------------------------------------------------------------
00109 JetCorrectorParameters::JetCorrectorParameters(const std::string& fFile, const std::string& fSection) 
00110 {
00111   std::ifstream input(fFile.c_str());
00112   std::string currentSection = "";
00113   std::string line;
00114   std::string currentDefinitions = "";
00115   while (std::getline(input,line)) 
00116     {
00117       std::string section = getSection(line);
00118       std::string tmp = getDefinitions(line);
00119       if (!section.empty() && tmp.empty()) 
00120         {
00121           currentSection = section;
00122           continue;
00123         }
00124       if (currentSection == fSection) 
00125         {
00126           if (!tmp.empty()) 
00127             {
00128               currentDefinitions = tmp;
00129               continue; 
00130             }
00131           Definitions definitions(currentDefinitions);
00132           if (!(definitions.nBinVar()==0 && definitions.formula()==""))
00133             mDefinitions = definitions;
00134           Record record(line,mDefinitions.nBinVar());
00135           bool check(true);
00136           for(unsigned i=0;i<mDefinitions.nBinVar();++i)
00137             if (record.xMin(i)==0 && record.xMax(i)==0)
00138               check = false;
00139           if (record.nParameters() == 0)
00140             check = false;  
00141           if (check)
00142             mRecords.push_back(record);
00143         } 
00144     }
00145   if (currentDefinitions=="")
00146     handleError("JetCorrectorParameters","No definitions found!!!");
00147   if (mRecords.empty() && currentSection == "") mRecords.push_back(Record());
00148   if (mRecords.empty() && currentSection != "") 
00149     {
00150       std::stringstream sserr; 
00151       sserr<<"the requested section "<<fSection<<" doesn't exist!";
00152       handleError("JetCorrectorParameters",sserr.str()); 
00153     }
00154   std::sort(mRecords.begin(), mRecords.end());
00155   valid_ = true;
00156 }
00157 //------------------------------------------------------------------------
00158 //--- returns the index of the record defined by fX ----------------------
00159 //------------------------------------------------------------------------
00160 int JetCorrectorParameters::binIndex(const std::vector<float>& fX) const 
00161 {
00162   int result = -1;
00163   unsigned N = mDefinitions.nBinVar();
00164   if (N != fX.size()) 
00165     {
00166       std::stringstream sserr; 
00167       sserr<<"# bin variables "<<N<<" doesn't correspont to requested #: "<<fX.size();
00168       handleError("JetCorrectorParameters",sserr.str());
00169     }
00170   unsigned tmp;
00171   for (unsigned i = 0; i < size(); ++i) 
00172     {
00173       tmp = 0;
00174       for (unsigned j=0;j<N;j++)
00175         if (fX[j] >= record(i).xMin(j) && fX[j] < record(i).xMax(j))
00176           tmp+=1;
00177       if (tmp==N)
00178         { 
00179           result = i;
00180           break;
00181         }
00182     } 
00183   return result;
00184 }
00185 //------------------------------------------------------------------------
00186 //--- returns the neighbouring bins of fIndex in the direction of fVar ---
00187 //------------------------------------------------------------------------
00188 int JetCorrectorParameters::neighbourBin(unsigned fIndex, unsigned fVar, bool fNext) const 
00189 {
00190   int result = -1;
00191   unsigned N = mDefinitions.nBinVar();
00192   if (fVar >= N) 
00193     {
00194       std::stringstream sserr; 
00195       sserr<<"# of bin variables "<<N<<" doesn't correspond to requested #: "<<fVar;
00196       handleError("JetCorrectorParameters",sserr.str()); 
00197     }
00198   unsigned tmp;
00199   for (unsigned i = 0; i < size(); ++i) 
00200     {
00201       tmp = 0;
00202       for (unsigned j=0;j<fVar;j++)
00203         if (fabs(record(i).xMin(j)-record(fIndex).xMin(j))<0.0001)
00204           tmp+=1;
00205       for (unsigned j=fVar+1;j<N;j++)
00206         if (fabs(record(i).xMin(j)-record(fIndex).xMin(j))<0.0001)
00207           tmp+=1;
00208       if (tmp<N-1)
00209         continue; 
00210       if (tmp==N-1)
00211         {
00212           if (fNext)
00213             if (fabs(record(i).xMin(fVar)-record(fIndex).xMax(fVar))<0.0001)
00214               tmp+=1;
00215           if (!fNext)
00216             if (fabs(record(i).xMax(fVar)-record(fIndex).xMin(fVar))<0.0001)
00217               tmp+=1;
00218         } 
00219       if (tmp==N)
00220         { 
00221           result = i;
00222           break;
00223         }
00224     } 
00225   return result;
00226 }
00227 //------------------------------------------------------------------------
00228 //--- returns the number of bins in the direction of fVar ----------------
00229 //------------------------------------------------------------------------
00230 unsigned JetCorrectorParameters::size(unsigned fVar) const
00231 {
00232   if (fVar >= mDefinitions.nBinVar()) 
00233     { 
00234       std::stringstream sserr; 
00235       sserr<<"requested bin variable index "<<fVar<<" is greater than number of variables "<<mDefinitions.nBinVar();
00236       handleError("JetCorrectorParameters",sserr.str()); 
00237     }    
00238   unsigned result = 0;
00239   float tmpMin(-9999),tmpMax(-9999);
00240   for (unsigned i = 0; i < size(); ++i)
00241     if (record(i).xMin(fVar) > tmpMin && record(i).xMax(fVar) > tmpMax)
00242       { 
00243         result++;
00244         tmpMin = record(i).xMin(fVar);
00245         tmpMax = record(i).xMax(fVar);
00246       }
00247   return result; 
00248 }
00249 //------------------------------------------------------------------------
00250 //--- returns the vector of bin centers of fVar --------------------------
00251 //------------------------------------------------------------------------
00252 std::vector<float> JetCorrectorParameters::binCenters(unsigned fVar) const 
00253 {
00254   std::vector<float> result;
00255   for (unsigned i = 0; i < size(); ++i)
00256     result.push_back(record(i).xMiddle(fVar));
00257   return result;
00258 }
00259 //------------------------------------------------------------------------
00260 //--- prints parameters on screen ----------------------------------------
00261 //------------------------------------------------------------------------
00262 void JetCorrectorParameters::printScreen() const
00263 {
00264   std::cout<<"--------------------------------------------"<<std::endl;
00265   std::cout<<"////////  PARAMETERS: //////////////////////"<<std::endl;
00266   std::cout<<"--------------------------------------------"<<std::endl;
00267   std::cout<<"Number of binning variables:   "<<definitions().nBinVar()<<std::endl;
00268   std::cout<<"Names of binning variables:    ";
00269   for(unsigned i=0;i<definitions().nBinVar();i++)
00270     std::cout<<definitions().binVar(i)<<" ";
00271   std::cout<<std::endl;
00272   std::cout<<"--------------------------------------------"<<std::endl;
00273   std::cout<<"Number of parameter variables: "<<definitions().nParVar()<<std::endl;
00274   std::cout<<"Names of parameter variables:  ";
00275   for(unsigned i=0;i<definitions().nParVar();i++)
00276     std::cout<<definitions().parVar(i)<<" ";
00277   std::cout<<std::endl;
00278   std::cout<<"--------------------------------------------"<<std::endl;
00279   std::cout<<"Parametrization Formula:       "<<definitions().formula()<<std::endl;
00280   if (definitions().isResponse())
00281     std::cout<<"Type (Response or Correction): "<<"Response"<<std::endl;
00282   else
00283     std::cout<<"Type (Response or Correction): "<<"Correction"<<std::endl;
00284   std::cout<<"Correction Level:              "<<definitions().level()<<std::endl;
00285   std::cout<<"--------------------------------------------"<<std::endl;
00286   std::cout<<"------- Bin contents -----------------------"<<std::endl;
00287   for(unsigned i=0;i<size();i++)
00288     {
00289       for(unsigned j=0;j<definitions().nBinVar();j++)
00290         std::cout<<record(i).xMin(j)<<" "<<record(i).xMax(j)<<" ";
00291       std::cout<<record(i).nParameters()<<" ";
00292       for(unsigned j=0;j<record(i).nParameters();j++)
00293         std::cout<<record(i).parameter(j)<<" ";
00294       std::cout<<std::endl;
00295     }  
00296 }
00297 //------------------------------------------------------------------------
00298 //--- prints parameters on file ----------------------------------------
00299 //------------------------------------------------------------------------
00300 void JetCorrectorParameters::printFile(const std::string& fFileName) const
00301 {
00302   std::ofstream txtFile;
00303   txtFile.open(fFileName.c_str());
00304   txtFile.setf(std::ios::right);
00305   txtFile<<"{"<<definitions().nBinVar()<<std::setw(15);
00306   for(unsigned i=0;i<definitions().nBinVar();i++)
00307     txtFile<<definitions().binVar(i)<<std::setw(15);
00308   txtFile<<definitions().nParVar()<<std::setw(15);
00309   for(unsigned i=0;i<definitions().nParVar();i++)
00310     txtFile<<definitions().parVar(i)<<std::setw(15);
00311   txtFile<<std::setw(definitions().formula().size()+15)<<definitions().formula()<<std::setw(15);
00312   if (definitions().isResponse())
00313     txtFile<<"Response"<<std::setw(15);
00314   else
00315     txtFile<<"Correction"<<std::setw(15);
00316   txtFile<<definitions().level()<<"}"<<"\n";
00317   for(unsigned i=0;i<size();i++)
00318     {
00319       for(unsigned j=0;j<definitions().nBinVar();j++)
00320         txtFile<<record(i).xMin(j)<<std::setw(15)<<record(i).xMax(j)<<std::setw(15);
00321       txtFile<<record(i).nParameters()<<std::setw(15);
00322       for(unsigned j=0;j<record(i).nParameters();j++)
00323         txtFile<<record(i).parameter(j)<<std::setw(15);
00324       txtFile<<"\n";
00325     }
00326   txtFile.close();
00327 }
00328 
00329 
00330 
00331 const char * 
00332 JetCorrectorParametersCollection::labelsArray_[JetCorrectorParametersCollection::N_LEVELS] = 
00333   {
00334     "L1Offset",
00335     "L2Relative",
00336     "L3Absolute",
00337     "L4EMF",
00338     "L5Flavor",
00339     "L6UE",
00340     "L7Parton",
00341     "L1JPTOffset",
00342     "L2L3Residual",
00343     "Uncertainty",
00344     "L1FastJet" 
00345   }; 
00346 
00347 const char *
00348 JetCorrectorParametersCollection::l5FlavorArray_[JetCorrectorParametersCollection::N_L5_SPECIES] = 
00349   {
00350     "L5Flavor_bJ",
00351     "L5Flavor_cJ",
00352     "L5Flavor_qJ",
00353     "L5Flavor_gJ",
00354     "L5Flavor_bT",
00355     "L5Flavor_cT",
00356     "L5Flavor_qT",
00357     "L5Flavor_gT"
00358   };
00359 
00360 const char *
00361 JetCorrectorParametersCollection::l7PartonArray_[JetCorrectorParametersCollection::N_L7_SPECIES] = 
00362   {
00363     "L7Parton_gJ",
00364     "L7Parton_qJ",
00365     "L7Parton_cJ",
00366     "L7Parton_bJ",
00367     "L7Parton_jJ",
00368     "L7Parton_qT",
00369     "L7Parton_cT",
00370     "L7Parton_bT",
00371     "L7Parton_tT"
00372   };
00373 
00374 
00375 std::vector<std::string>
00376 JetCorrectorParametersCollection::labels_(labelsArray_, 
00377                                           labelsArray_ + sizeof(labelsArray_)/sizeof(*labelsArray_) );
00378 
00379 std::vector<std::string>
00380 JetCorrectorParametersCollection::l5Flavors_(l5FlavorArray_, 
00381                                              l5FlavorArray_ + sizeof(l5FlavorArray_)/sizeof(*l5FlavorArray_) );
00382 
00383 std::vector<std::string>
00384 JetCorrectorParametersCollection::l7Partons_(l7PartonArray_, 
00385                                              l7PartonArray_ + sizeof(l7PartonArray_)/sizeof(*l7PartonArray_) );
00386 
00387 
00388 void JetCorrectorParametersCollection::getSections( std::string inputFile,
00389                                                     std::vector<std::string> & outputs )
00390 {
00391   outputs.clear();
00392   std::ifstream input( inputFile.c_str() );
00393   while( !input.eof() ) {
00394     char buff[10000];
00395     input.getline(buff,10000);
00396     std::string in(buff); 
00397     if ( in[0] == '[' ) {
00398       std::string tok = getSection(in);
00399       if ( tok != "" ) {
00400         outputs.push_back( tok );
00401       }
00402     }
00403   }
00404   std::cout << "Found these sections for file: " << std::endl;
00405   copy(outputs.begin(),outputs.end(), std::ostream_iterator<std::string>(std::cout, "\n") );
00406 } 
00407 
00408 
00409 // Add a JetCorrectorParameter object, possibly with flavor. 
00410 void JetCorrectorParametersCollection::push_back( key_type i, value_type const & j, label_type const & flav) { 
00411   std::cout << "i    = " << i << std::endl;  
00412   std::cout << "flav = " << flav << std::endl;
00413   if ( isL5(i) ) {
00414     std::cout << "This is L5, getL5Bin = " << getL5Bin(flav) << std::endl;
00415     correctionsL5_.push_back( pair_type(getL5Bin(flav),j) ); 
00416   }
00417   else if ( isL7(i) ) {
00418     std::cout << "This is L7, getL7Bin = " << getL7Bin(flav) << std::endl;
00419     correctionsL7_.push_back( pair_type(getL7Bin(flav),j) );
00420   }
00421   else if ( flav == "" ) {
00422     corrections_.push_back( pair_type(i,j) );
00423   } else {
00424     std::cout << "***** NOT ADDING " << flav << ", corresponding position in JetCorrectorParameters is not found." << std::endl;
00425   }
00426 }
00427 
00428 
00429 // Access the JetCorrectorParameter via the key k.
00430 // key_type is hashed to deal with the three collections
00431 JetCorrectorParameters const & JetCorrectorParametersCollection::operator[]( key_type k ) const {
00432   collection_type::const_iterator ibegin, iend, i;
00433   if ( isL5(k) ) {
00434     ibegin = correctionsL5_.begin();
00435     iend = correctionsL5_.end();
00436     i = ibegin;
00437   } else if ( isL7(k) ) {
00438     ibegin = correctionsL7_.begin();
00439     iend = correctionsL7_.end();
00440     i = ibegin;      
00441   } else { 
00442     ibegin = corrections_.begin();
00443     iend = corrections_.end();
00444     i = ibegin;
00445   }
00446   for ( ; i != iend; ++i ) {
00447     if ( k == i->first ) return i->second;
00448   }
00449   throw cms::Exception("InvalidInput") << " cannot find key " << static_cast<int>(k) << std::endl;
00450 }
00451 
00452 // Get a list of valid keys. These will contain hashed keys
00453 // that are aware of all three collections. 
00454 void JetCorrectorParametersCollection::validKeys(std::vector<key_type> & keys ) const {
00455   keys.clear();
00456   for ( collection_type::const_iterator ibegin = corrections_.begin(),
00457           iend = corrections_.end(), i = ibegin; i != iend; ++i ) {
00458     keys.push_back( i->first );
00459   }
00460   for ( collection_type::const_iterator ibegin = correctionsL5_.begin(),
00461           iend = correctionsL5_.end(), i = ibegin; i != iend; ++i ) {
00462     keys.push_back( i->first );
00463   }
00464   for ( collection_type::const_iterator ibegin = correctionsL7_.begin(),
00465           iend = correctionsL7_.end(), i = ibegin; i != iend; ++i ) {
00466     keys.push_back( i->first );
00467   }
00468 }
00469 
00470 
00471 // Find the L5 bin for hashing
00472 JetCorrectorParametersCollection::key_type
00473 JetCorrectorParametersCollection::getL5Bin( std::string const & flav ){
00474   std::vector<std::string>::const_iterator found = 
00475     find( l5Flavors_.begin(), l5Flavors_.end(), flav );
00476   if ( found != l5Flavors_.end() ) {
00477     return (found - l5Flavors_.begin() + 1) * 100;
00478   }
00479   else return L5Flavor;
00480 }
00481 // Find the L7 bin for hashing
00482 JetCorrectorParametersCollection::key_type
00483 JetCorrectorParametersCollection::getL7Bin( std::string const & flav ){
00484   std::vector<std::string>::const_iterator found = 
00485     find( l7Partons_.begin(), l7Partons_.end(), flav );
00486   if ( found != l7Partons_.end() ) {
00487     return (found - l7Partons_.begin() + 1) * 1000;
00488   }
00489   else return L7Parton;
00490 }
00491 
00492 // Check if this is an L5 hashed value
00493 bool JetCorrectorParametersCollection::isL5( key_type k ) {
00494   return k == L5Flavor ||
00495     ( k / 100 > 0 && k / 1000 == 0 );
00496 }
00497 // Check if this is an L7 hashed value
00498 bool JetCorrectorParametersCollection::isL7( key_type k ) {
00499   return k == L7Parton ||
00500     ( k / 1000 > 0 );
00501 }
00502 
00503 
00504 // Find the key corresponding to each label
00505 JetCorrectorParametersCollection::key_type 
00506 JetCorrectorParametersCollection::findKey( std::string const & label ) const {
00507 
00508   // First check L5 corrections
00509   std::vector<std::string>::const_iterator found1 =
00510     find( l5Flavors_.begin(), l5Flavors_.end(), label );
00511   if ( found1 != l5Flavors_.end() ) {
00512     return getL5Bin(label);
00513   } 
00514 
00515   // Next check L7 corrections
00516   std::vector<std::string>::const_iterator found2 =
00517     find( l7Partons_.begin(), l7Partons_.end(), label );
00518   if ( found2 != l7Partons_.end() ) {
00519     return getL7Bin(label);
00520   } 
00521 
00522   // Finally check the default corrections
00523   std::vector<std::string>::const_iterator found3 =
00524     find( labels_.begin(), labels_.end(), label );
00525   if ( found3 != labels_.end() ) {
00526     return static_cast<key_type>(found3 - labels_.begin());
00527   } 
00528 
00529   // Didn't find default corrections, throw exception
00530   throw cms::Exception("InvalidInput") << " Cannot find label " << label << std::endl;
00531 
00532 }
00533 
00534 
00535 //#include "FWCore/Framework/interface/EventSetup.h"
00536 //#include "FWCore/Framework/interface/ESHandle.h"
00537 //#include "FWCore/Framework/interface/ModuleFactory.h"
00538 #include "FWCore/Utilities/interface/typelookup.h"
00539  
00540 TYPELOOKUP_DATA_REG(JetCorrectorParameters);
00541 TYPELOOKUP_DATA_REG(JetCorrectorParametersCollection);