CMS 3D CMS Logo

JetCorrectorParameters.h
Go to the documentation of this file.
1 //
2 // Original Author: Fedor Ratnikov Nov 9, 2007
3 // $Id: JetCorrectorParameters.h,v 1.15 2012/03/01 18:24:52 srappocc Exp $
4 //
5 // Generic parameters for Jet corrections
6 //
7 #ifndef JetCorrectorParameters_h
8 #define JetCorrectorParameters_h
9 
12 
13 #include <string>
14 #include <vector>
15 #include <tuple>
16 #include <algorithm>
17 #include <functional>
18 #include <iostream>
19 #include <memory>
22 
24 
26 {
27  //---------------- JetCorrectorParameters class ----------------
28  //-- Encapsulates all the information of the parametrization ---
29  public:
30  //---------------- Definitions class ---------------------------
31  //-- Global iformation about the parametrization is kept here --
32  class Definitions
33  {
34  public:
35  //-------- Constructors --------------
37  Definitions(const std::vector<std::string>& fBinVar, const std::vector<std::string>& fParVar, const std::string& fFormula, bool fIsResponse);
38  Definitions(const std::string& fLine);
39  //-------- Member functions ----------
40  unsigned nBinVar() const {return mBinVar.size(); }
41  unsigned nParVar() const {return mParVar.size(); }
42  std::vector<std::string> parVar() const {return mParVar; }
43  std::vector<std::string> binVar() const {return mBinVar; }
44  std::string parVar(unsigned fIndex) const {return mParVar[fIndex];}
45  std::string binVar(unsigned fIndex) const {return mBinVar[fIndex];}
46  std::string formula() const {return mFormula; }
47  std::string level() const {return mLevel; }
48  bool isResponse() const {return mIsResponse; }
49  private:
50  //-------- Member variables ----------
51  bool mIsResponse;
54  std::vector<std::string> mParVar;
55  std::vector<std::string> mBinVar;
56 
58  };
59  //---------------- Record class --------------------------------
60  //-- Each Record holds the properties of a bin -----------------
61  class Record
62  {
63  public:
64  //-------- Constructors --------------
65  Record() : mNvar(0),mMin(0),mMax(0) {}
66  Record(unsigned fNvar, const std::vector<float>& fXMin, const std::vector<float>& fXMax, const std::vector<float>& fParameters) : mNvar(fNvar),mMin(fXMin),mMax(fXMax),mParameters(fParameters) {}
67  Record(const std::string& fLine, unsigned fNvar);
68  //-------- Member functions ----------
69  unsigned nVar() const {return mNvar; }
70  float xMin(unsigned fVar) const {return mMin[fVar]; }
71  float xMax(unsigned fVar) const {return mMax[fVar]; }
72  float xMiddle(unsigned fVar) const {return 0.5*(xMin(fVar)+xMax(fVar));}
73  float parameter(unsigned fIndex) const {return mParameters[fIndex]; }
74  std::vector<float> parameters() const {return mParameters; }
75  unsigned nParameters() const {return mParameters.size(); }
76  bool operator< (const Record& other) const
77  {
78  if (xMin(0) < other.xMin(0)) return true;
79  if (xMin(0) > other.xMin(0)) return false;
80  if (xMin(1) < other.xMin(1)) return true;
81  if (xMin(1) > other.xMin(1)) return false;
82  return (xMin(2) < other.xMin(2));
83  }
84  private:
85  //-------- Member variables ----------
86  unsigned mNvar;
87  std::vector<float> mMin;
88  std::vector<float> mMax;
89  std::vector<float> mParameters;
90 
92  };
93 
94  //-------- Constructors --------------
96  JetCorrectorParameters(const std::string& fFile, const std::string& fSection = "");
98  const std::vector<JetCorrectorParameters::Record>& fRecords)
99  : mDefinitions(fDefinitions),mRecords(fRecords) { valid_ = true;}
100  //-------- Member functions ----------
101  const Record& record(unsigned fBin) const {return mRecords[fBin]; }
102  const Definitions& definitions() const {return mDefinitions; }
103  unsigned size() const {return mRecords.size();}
104  unsigned size(unsigned fVar) const;
105  int binIndex(const std::vector<float>& fX) const;
106  int binIndexN(const std::vector<float>& fX) const;
107  int neighbourBin(unsigned fIndex, unsigned fVar, bool fNext) const;
108  std::vector<float> binCenters(unsigned fVar) const;
109  void printScreen() const;
110  void printFile(const std::string& fFileName) const;
111  bool isValid() const { return valid_; }
112  void init();
113 
115 
116  private:
117  //-------- Member variables ----------
119  std::vector<JetCorrectorParameters::Record> mRecords;
120  bool valid_;
121 
122  std::shared_ptr<JetCorrectorParametersHelper> helper COND_TRANSIENT;
123 
125 };
126 std::ostream& operator<<(std::ostream& out, const JetCorrectorParameters::Record& fBin);
127 
128 
130  //---------------- JetCorrectorParametersCollection class ----------------
131  //-- Adds several JetCorrectorParameters together by algorithm type ---
132  //-- to reduce the number of payloads in the Database ---
133  //-- NB: The enum is listed in "human-logical" order, but the actual
134  // enum value is in order of appearance when people thought of them.
135  public:
136  enum Level_t { L1Offset=0,
137  L1JPTOffset=7,
138  L1FastJet = 10,
139  L2Relative=1,
140  L3Absolute=2,
141  L2L3Residual=8,
142  L4EMF=3,
143  L5Flavor=4,
144  L6UE=5,
145  L7Parton=6,
146  Uncertainty=9,
147  UncertaintyAbsolute=11,
148  UncertaintyHighPtExtra=12,
149  UncertaintySinglePionECAL=13,
150  UncertaintySinglePionHCAL=27,
151  UncertaintyFlavor=14,
152  UncertaintyTime=15,
153  UncertaintyRelativeJEREC1=16,
154  UncertaintyRelativeJEREC2=17,
155  UncertaintyRelativeJERHF=18,
156  UncertaintyRelativePtEC1=28,
157  UncertaintyRelativePtEC2=29,
158  UncertaintyRelativePtHF=30,
159  UncertaintyRelativeStatEC2=19,
160  UncertaintyRelativeStatHF=20,
161  UncertaintyRelativeFSR=21,
162  UncertaintyRelativeSample=31,
163  UncertaintyPileUpDataMC=22,
164  UncertaintyPileUpOOT=23,
165  UncertaintyPileUpPtBB=24,
166  UncertaintyPileUpPtEC=32,
167  UncertaintyPileUpPtHF=33,
168  UncertaintyPileUpBias=25,
169  UncertaintyPileUpJetRate=26,
170  L1RC=34,
171  L1Residual=35,
172  UncertaintyAux3=36,
173  UncertaintyAux4=37,
174  N_LEVELS=38
175  };
176 
177 
178  enum L5_Species_t {L5_bJ=0,L5_cJ,L5_qJ,L5_gJ,L5_bT,L5_cT,L5_qT,L5_gT,N_L5_SPECIES};
179  enum L7_Species_t {L7_gJ=0,L7_qJ,L7_cJ,L7_bJ,L7_jJ,L7_qT,L7_cT,L7_bT,L7_jT,N_L7_SPECIES};
180  typedef int key_type;
183  typedef std::pair<key_type,value_type> pair_type;
184  typedef std::vector<pair_type> collection_type;
185 
186 
187  // Constructor... initialize all three vectors to zero
188  JetCorrectorParametersCollection() { corrections_.clear(); correctionsL5_.clear(); correctionsL7_.clear(); }
189 
190  // Add a JetCorrectorParameter object, possibly with flavor.
191  void push_back( key_type i, value_type const & j, label_type const & flav = "" );
192 
193  // Access the JetCorrectorParameter via the key k.
194  // key_type is hashed to deal with the three collections
195  JetCorrectorParameters const & operator[]( key_type k ) const;
196 
197  // Access the JetCorrectorParameter via a string.
198  // Will find the hashed value for the label, and call via that
199  // operator.
201  return operator[]( findKey(label) );
202  }
203 
204  // Get a list of valid keys. These will contain hashed keys
205  // that are aware of all three collections.
206  void validKeys(std::vector<key_type> & keys ) const;
207 
208 
209 
210  // Helper method to find all of the sections in a given
211  // parameters file
212  static void getSections(std::string inputFile,
213  std::vector<std::string> & outputs );
214 
215  // Find the L5 bin for hashing
216  static key_type getL5Bin( std::string const & flav );
217  // Find the L7 bin for hashing
218  static key_type getL7Bin( std::string const & flav );
219  // Check if this is an L5 hashed value
220  static bool isL5( key_type k );
221  // Check if this is an L7 hashed value
222  static bool isL7( key_type k );
223 
224  static std::string findLabel( key_type k );
225  static std::string findL5Flavor( key_type k );
226  static std::string findL7Parton( key_type k );
227 
228  protected:
229 
230  // Find the key corresponding to each label
231  key_type findKey( std::string const & label ) const;
232 
233  collection_type corrections_;
234  collection_type correctionsL5_;
235  collection_type correctionsL7_;
236 
237  collection_type& getCorrections() {return corrections_;}
238  collection_type& getCorrectionsL5() {return correctionsL5_;}
239  collection_type& getCorrectionsL7() {return correctionsL7_;}
240 
242 
244 
245 };
246 
249  for (auto & ptype : jcpc.getCorrections()) {ptype.second.init();}
250  for (auto & ptype : jcpc.getCorrectionsL5()) {ptype.second.init();}
251  for (auto & ptype : jcpc.getCorrectionsL7()) {ptype.second.init();}
252  }
253 };
254 
255 #endif
float xMin(unsigned fVar) const
Definition: helper.py:1
JetCorrectorParameters const & operator[](std::string const &label) const
std::pair< key_type, value_type > pair_type
std::string parVar(unsigned fIndex) const
std::vector< float > parameters() const
std::string binVar(unsigned fIndex) const
const Definitions & definitions() const
std::vector< pair_type > collection_type
std::vector< JetCorrectorParameters::Record > mRecords
const Record & record(unsigned fBin) const
std::vector< std::string > parVar() const
float parameter(unsigned fIndex) const
char const * label
void operator()(JetCorrectorParametersCollection &jcpc)
void printFile(const std::string &fFileName) const
int binIndex(const std::vector< float > &fX) const
std::vector< float > binCenters(unsigned fVar) const
Record(unsigned fNvar, const std::vector< float > &fXMin, const std::vector< float > &fXMax, const std::vector< float > &fParameters)
T operator[](int i) const
static const int MAX_SIZE_DIMENSIONALITY
float xMiddle(unsigned fVar) const
float xMax(unsigned fVar) const
int k[5][pyjets_maxn]
#define COND_TRANSIENT
Definition: Serializable.h:61
#define COND_SERIALIZABLE
Definition: Serializable.h:38
bool operator<(DTCELinkId const &lhs, DTCELinkId const &rhs)
Definition: DTCELinkId.h:73
int binIndexN(const std::vector< float > &fX) const
std::ostream & operator<<(std::ostream &out, const JetCorrectorParameters::Record &fBin)
std::vector< std::string > binVar() const
JetCorrectorParameters::Definitions mDefinitions
int neighbourBin(unsigned fIndex, unsigned fVar, bool fNext) const
JetCorrectorParameters(const JetCorrectorParameters::Definitions &fDefinitions, const std::vector< JetCorrectorParameters::Record > &fRecords)