CMS 3D CMS Logo

Classes | Public Member Functions | Private Attributes

JetCorrectorParameters Class Reference

#include <JetCorrectorParameters.h>

List of all members.

Classes

class  Definitions
class  Record

Public Member Functions

std::vector< float > binCenters (unsigned fVar) const
int binIndex (const std::vector< float > &fX) const
const Definitionsdefinitions () const
bool isValid () const
 JetCorrectorParameters ()
 JetCorrectorParameters (const std::string &fFile, const std::string &fSection="")
 JetCorrectorParameters (const JetCorrectorParameters::Definitions &fDefinitions, const std::vector< JetCorrectorParameters::Record > &fRecords)
int neighbourBin (unsigned fIndex, unsigned fVar, bool fNext) const
void printFile (const std::string &fFileName) const
void printScreen () const
const Recordrecord (unsigned fBin) const
unsigned size (unsigned fVar) const
unsigned size () const

Private Attributes

JetCorrectorParameters::Definitions mDefinitions
std::vector
< JetCorrectorParameters::Record
mRecords
bool valid_

Detailed Description

Definition at line 16 of file JetCorrectorParameters.h.


Constructor & Destructor Documentation

JetCorrectorParameters::JetCorrectorParameters ( ) [inline]

Definition at line 74 of file JetCorrectorParameters.h.

References valid_.

{ valid_ = false;}
JetCorrectorParameters::JetCorrectorParameters ( const std::string &  fFile,
const std::string &  fSection = "" 
)

Definition at line 109 of file JetCorrectorParameters.cc.

References CastorDataFrameFilter_impl::check(), definitions(), JetCorrectorParameters::Definitions::formula(), insertMaterial::getSection(), i, LaserDQM_cfg::input, geometryCSVtoXML::line, mDefinitions, mRecords, JetCorrectorParameters::Definitions::nBinVar(), JetCorrectorParameters::Record::nParameters(), record(), python::multivaluedict::sort(), AlCaHLTBitMon_QueryRunRegistry::string, tmp, valid_, JetCorrectorParameters::Record::xMax(), and JetCorrectorParameters::Record::xMin().

{
  std::ifstream input(fFile.c_str());
  std::string currentSection = "";
  std::string line;
  std::string currentDefinitions = "";
  while (std::getline(input,line)) 
    {
      std::string section = getSection(line);
      std::string tmp = getDefinitions(line);
      if (!section.empty() && tmp.empty()) 
        {
          currentSection = section;
          continue;
        }
      if (currentSection == fSection) 
        {
          if (!tmp.empty()) 
            {
              currentDefinitions = tmp;
              continue; 
            }
          Definitions definitions(currentDefinitions);
          if (!(definitions.nBinVar()==0 && definitions.formula()==""))
            mDefinitions = definitions;
          Record record(line,mDefinitions.nBinVar());
          bool check(true);
          for(unsigned i=0;i<mDefinitions.nBinVar();++i)
            if (record.xMin(i)==0 && record.xMax(i)==0)
              check = false;
          if (record.nParameters() == 0)
            check = false;  
          if (check)
            mRecords.push_back(record);
        } 
    }
  if (currentDefinitions=="")
    handleError("JetCorrectorParameters","No definitions found!!!");
  if (mRecords.empty() && currentSection == "") mRecords.push_back(Record());
  if (mRecords.empty() && currentSection != "") 
    {
      std::stringstream sserr; 
      sserr<<"the requested section "<<fSection<<" doesn't exist!";
      handleError("JetCorrectorParameters",sserr.str()); 
    }
  std::sort(mRecords.begin(), mRecords.end());
  valid_ = true;
}
JetCorrectorParameters::JetCorrectorParameters ( const JetCorrectorParameters::Definitions fDefinitions,
const std::vector< JetCorrectorParameters::Record > &  fRecords 
) [inline]

Definition at line 76 of file JetCorrectorParameters.h.

References valid_.

      : mDefinitions(fDefinitions),mRecords(fRecords) { valid_ = true;}

Member Function Documentation

std::vector< float > JetCorrectorParameters::binCenters ( unsigned  fVar) const

Definition at line 252 of file JetCorrectorParameters.cc.

References i, record(), query::result, and size().

{
  std::vector<float> result;
  for (unsigned i = 0; i < size(); ++i)
    result.push_back(record(i).xMiddle(fVar));
  return result;
}
int JetCorrectorParameters::binIndex ( const std::vector< float > &  fX) const

Definition at line 160 of file JetCorrectorParameters.cc.

References i, j, mDefinitions, N, JetCorrectorParameters::Definitions::nBinVar(), record(), query::result, size(), tmp, JetCorrectorParameters::Record::xMax(), and JetCorrectorParameters::Record::xMin().

Referenced by SimpleJetCorrector::correction(), and SimpleJetCorrectionUncertainty::uncertainty().

{
  int result = -1;
  unsigned N = mDefinitions.nBinVar();
  if (N != fX.size()) 
    {
      std::stringstream sserr; 
      sserr<<"# bin variables "<<N<<" doesn't correspont to requested #: "<<fX.size();
      handleError("JetCorrectorParameters",sserr.str());
    }
  unsigned tmp;
  for (unsigned i = 0; i < size(); ++i) 
    {
      tmp = 0;
      for (unsigned j=0;j<N;j++)
        if (fX[j] >= record(i).xMin(j) && fX[j] < record(i).xMax(j))
          tmp+=1;
      if (tmp==N)
        { 
          result = i;
          break;
        }
    } 
  return result;
}
const Definitions& JetCorrectorParameters::definitions ( ) const [inline]
bool JetCorrectorParameters::isValid ( void  ) const [inline]

Definition at line 89 of file JetCorrectorParameters.h.

References valid_.

{ return valid_; }
int JetCorrectorParameters::neighbourBin ( unsigned  fIndex,
unsigned  fVar,
bool  fNext 
) const

Definition at line 188 of file JetCorrectorParameters.cc.

References i, j, mDefinitions, N, JetCorrectorParameters::Definitions::nBinVar(), record(), query::result, size(), and tmp.

Referenced by SimpleJetCorrector::correction().

{
  int result = -1;
  unsigned N = mDefinitions.nBinVar();
  if (fVar >= N) 
    {
      std::stringstream sserr; 
      sserr<<"# of bin variables "<<N<<" doesn't correspond to requested #: "<<fVar;
      handleError("JetCorrectorParameters",sserr.str()); 
    }
  unsigned tmp;
  for (unsigned i = 0; i < size(); ++i) 
    {
      tmp = 0;
      for (unsigned j=0;j<fVar;j++)
        if (fabs(record(i).xMin(j)-record(fIndex).xMin(j))<0.0001)
          tmp+=1;
      for (unsigned j=fVar+1;j<N;j++)
        if (fabs(record(i).xMin(j)-record(fIndex).xMin(j))<0.0001)
          tmp+=1;
      if (tmp<N-1)
        continue; 
      if (tmp==N-1)
        {
          if (fNext)
            if (fabs(record(i).xMin(fVar)-record(fIndex).xMax(fVar))<0.0001)
              tmp+=1;
          if (!fNext)
            if (fabs(record(i).xMax(fVar)-record(fIndex).xMin(fVar))<0.0001)
              tmp+=1;
        } 
      if (tmp==N)
        { 
          result = i;
          break;
        }
    } 
  return result;
}
void JetCorrectorParameters::printFile ( const std::string &  fFileName) const

Definition at line 300 of file JetCorrectorParameters.cc.

References JetCorrectorParameters::Definitions::binVar(), definitions(), JetCorrectorParameters::Definitions::formula(), EcalElecEmulExample_cfg::formula, i, j, JetCorrectorParameters::Definitions::level(), JetCorrectorParameters::Definitions::nBinVar(), JetCorrectorParameters::Record::nParameters(), JetCorrectorParameters::Definitions::nParVar(), JetCorrectorParameters::Record::parameter(), JetCorrectorParameters::Definitions::parVar(), record(), size(), JetCorrectorParameters::Record::xMax(), and JetCorrectorParameters::Record::xMin().

Referenced by JetCorrectorDBReader::analyze().

{
  std::ofstream txtFile;
  txtFile.open(fFileName.c_str());
  txtFile.setf(std::ios::right);
  txtFile<<"{"<<definitions().nBinVar()<<std::setw(15);
  for(unsigned i=0;i<definitions().nBinVar();i++)
    txtFile<<definitions().binVar(i)<<std::setw(15);
  txtFile<<definitions().nParVar()<<std::setw(15);
  for(unsigned i=0;i<definitions().nParVar();i++)
    txtFile<<definitions().parVar(i)<<std::setw(15);
  txtFile<<std::setw(definitions().formula().size()+15)<<definitions().formula()<<std::setw(15);
  if (definitions().isResponse())
    txtFile<<"Response"<<std::setw(15);
  else
    txtFile<<"Correction"<<std::setw(15);
  txtFile<<definitions().level()<<"}"<<"\n";
  for(unsigned i=0;i<size();i++)
    {
      for(unsigned j=0;j<definitions().nBinVar();j++)
        txtFile<<record(i).xMin(j)<<std::setw(15)<<record(i).xMax(j)<<std::setw(15);
      txtFile<<record(i).nParameters()<<std::setw(15);
      for(unsigned j=0;j<record(i).nParameters();j++)
        txtFile<<record(i).parameter(j)<<std::setw(15);
      txtFile<<"\n";
    }
  txtFile.close();
}
void JetCorrectorParameters::printScreen ( ) const

Definition at line 262 of file JetCorrectorParameters.cc.

References JetCorrectorParameters::Definitions::binVar(), gather_cfg::cout, definitions(), JetCorrectorParameters::Definitions::formula(), i, j, JetCorrectorParameters::Definitions::level(), JetCorrectorParameters::Definitions::nBinVar(), JetCorrectorParameters::Record::nParameters(), JetCorrectorParameters::Definitions::nParVar(), JetCorrectorParameters::Record::parameter(), JetCorrectorParameters::Definitions::parVar(), record(), size(), JetCorrectorParameters::Record::xMax(), and JetCorrectorParameters::Record::xMin().

Referenced by JetCorrectorDBReader::analyze().

{
  std::cout<<"--------------------------------------------"<<std::endl;
  std::cout<<"////////  PARAMETERS: //////////////////////"<<std::endl;
  std::cout<<"--------------------------------------------"<<std::endl;
  std::cout<<"Number of binning variables:   "<<definitions().nBinVar()<<std::endl;
  std::cout<<"Names of binning variables:    ";
  for(unsigned i=0;i<definitions().nBinVar();i++)
    std::cout<<definitions().binVar(i)<<" ";
  std::cout<<std::endl;
  std::cout<<"--------------------------------------------"<<std::endl;
  std::cout<<"Number of parameter variables: "<<definitions().nParVar()<<std::endl;
  std::cout<<"Names of parameter variables:  ";
  for(unsigned i=0;i<definitions().nParVar();i++)
    std::cout<<definitions().parVar(i)<<" ";
  std::cout<<std::endl;
  std::cout<<"--------------------------------------------"<<std::endl;
  std::cout<<"Parametrization Formula:       "<<definitions().formula()<<std::endl;
  if (definitions().isResponse())
    std::cout<<"Type (Response or Correction): "<<"Response"<<std::endl;
  else
    std::cout<<"Type (Response or Correction): "<<"Correction"<<std::endl;
  std::cout<<"Correction Level:              "<<definitions().level()<<std::endl;
  std::cout<<"--------------------------------------------"<<std::endl;
  std::cout<<"------- Bin contents -----------------------"<<std::endl;
  for(unsigned i=0;i<size();i++)
    {
      for(unsigned j=0;j<definitions().nBinVar();j++)
        std::cout<<record(i).xMin(j)<<" "<<record(i).xMax(j)<<" ";
      std::cout<<record(i).nParameters()<<" ";
      for(unsigned j=0;j<record(i).nParameters();j++)
        std::cout<<record(i).parameter(j)<<" ";
      std::cout<<std::endl;
    }  
}
const Record& JetCorrectorParameters::record ( unsigned  fBin) const [inline]
unsigned JetCorrectorParameters::size ( unsigned  fVar) const

Definition at line 230 of file JetCorrectorParameters.cc.

References i, mDefinitions, JetCorrectorParameters::Definitions::nBinVar(), record(), query::result, size(), JetCorrectorParameters::Record::xMax(), and JetCorrectorParameters::Record::xMin().

{
  if (fVar >= mDefinitions.nBinVar()) 
    { 
      std::stringstream sserr; 
      sserr<<"requested bin variable index "<<fVar<<" is greater than number of variables "<<mDefinitions.nBinVar();
      handleError("JetCorrectorParameters",sserr.str()); 
    }    
  unsigned result = 0;
  float tmpMin(-9999),tmpMax(-9999);
  for (unsigned i = 0; i < size(); ++i)
    if (record(i).xMin(fVar) > tmpMin && record(i).xMax(fVar) > tmpMax)
      { 
        result++;
        tmpMin = record(i).xMin(fVar);
        tmpMax = record(i).xMax(fVar);
      }
  return result; 
}
unsigned JetCorrectorParameters::size ( void  ) const [inline]

Member Data Documentation

Definition at line 94 of file JetCorrectorParameters.h.

Referenced by JetCorrectorParameters(), record(), and size().

Definition at line 95 of file JetCorrectorParameters.h.

Referenced by isValid(), and JetCorrectorParameters().