CMS 3D CMS Logo

MESetUtils.cc
Go to the documentation of this file.
2 
12 
15 
16 namespace ecaldqm
17 {
18  MESet*
19  createMESet(edm::ParameterSet const& _MEParam)
20  {
25 
26  binning::AxisSpecs xaxis, yaxis, zaxis;
27  bool hasXaxis(_MEParam.existsAs<edm::ParameterSet>("xaxis", false));
28  if(hasXaxis)
29  xaxis = binning::formAxis(_MEParam.getUntrackedParameterSet("xaxis"));
30  bool hasYaxis(_MEParam.existsAs<edm::ParameterSet>("yaxis", false));
31  if(hasYaxis)
32  yaxis = binning::formAxis(_MEParam.getUntrackedParameterSet("yaxis"));
33  bool hasZaxis(_MEParam.existsAs<edm::ParameterSet>("zaxis", false));
34  if(hasZaxis)
35  zaxis = binning::formAxis(_MEParam.getUntrackedParameterSet("zaxis"));
36 
37  MESet* set(nullptr);
38 
39  if(btype == binning::kTrend){
40  MESetTrend* setTrend(new MESetTrend(path, otype, btype, kind, hasYaxis ? &yaxis : nullptr));
41  if(_MEParam.existsAs<bool>("minutely", false) && _MEParam.getUntrackedParameter<bool>("minutely")) setTrend->setMinutely();
42  if(_MEParam.existsAs<bool>("cumulative", false) && _MEParam.getUntrackedParameter<bool>("cumulative")) setTrend->setCumulative();
43  if(_MEParam.existsAs<bool>("shiftAxis", false) && _MEParam.getUntrackedParameter<bool>("shiftAxis")) setTrend->setShiftAxis();
44  set = setTrend;
45  }
46  else if(otype == binning::nObjType)
47  set = new MESetNonObject(path, otype, btype, kind, hasXaxis ? &xaxis : nullptr, hasYaxis ? &yaxis : nullptr, hasZaxis ? &zaxis : nullptr);
48  else if(otype == binning::kChannel)
49 // Class removed until concurrency issue is finalized
50 #if 0
51  set = new MESetChannel(path, otype, btype, kind);
52 #else
53  set = nullptr;
54 #endif
55  else if(btype == binning::kProjEta || btype == binning::kProjPhi)
56  set = new MESetProjection(path, otype, btype, kind, hasYaxis ? &yaxis : nullptr);
57  else{
58  unsigned logicalDimensions(-1);
59  switch(kind){
61  logicalDimensions = 0;
62  break;
65  logicalDimensions = 1;
66  break;
69  logicalDimensions = 2;
70  break;
71  default:
72  break;
73  }
74 
75  // example case: Ecal/TriggerPrimitives/EmulMatching/TrigPrimTask matching index
76  if(logicalDimensions == 2 && hasYaxis && btype != binning::kUser) logicalDimensions = 1;
77 
78  if(logicalDimensions > 2 || (btype == binning::kReport && logicalDimensions != 0))
79  throw cms::Exception("InvalidConfiguration") << "Cannot create MESet at " << path;
80 
81  if(btype == binning::kUser)
82  set = new MESetEcal(path, otype, btype, kind, logicalDimensions, hasXaxis ? &xaxis : nullptr, hasYaxis ? &yaxis : nullptr, hasZaxis ? &zaxis : nullptr);
83  else if(logicalDimensions == 0)
84  set = new MESetDet0D(path, otype, btype, kind);
85  else if(logicalDimensions == 1)
86  set = new MESetDet1D(path, otype, btype, kind, hasYaxis ? &yaxis : nullptr);
87  else if(logicalDimensions == 2)
88  set = new MESetDet2D(path, otype, btype, kind, hasZaxis ? &zaxis : nullptr);
89  }
90 
91  if(_MEParam.existsAs<edm::ParameterSet>("multi", false)){
92  typedef std::vector<std::string> VString;
93 
94  edm::ParameterSet const& multiParams(_MEParam.getUntrackedParameterSet("multi"));
95  VString replacementNames(multiParams.getParameterNames());
96  if(replacementNames.empty())
97  throw cms::Exception("InvalidConfiguration") << "0 multiplicity for MESet at " << path;
98 
100  for(unsigned iD(0); iD != replacementNames.size(); ++iD){
101  VString reps;
102  if(multiParams.existsAs<VString>(replacementNames[iD], false))
103  reps = multiParams.getUntrackedParameter<VString>(replacementNames[iD]);
104  else if(multiParams.existsAs<std::vector<int> >(replacementNames[iD], false)){
105  std::vector<int> repInts(multiParams.getUntrackedParameter<std::vector<int> >(replacementNames[iD]));
106  for(unsigned iR(0); iR != repInts.size(); ++iR)
107  reps.push_back(std::to_string(repInts[iR]));
108  }
109 
110  if(reps.empty())
111  throw cms::Exception("InvalidConfiguration") << "0 multiplicity for MESet at " << path;
112 
113  candidates[replacementNames[iD]] = reps;
114  }
115  MESetMulti* multi(new MESetMulti(*set, candidates));
116  delete set;
117  set = multi;
118  }
119 
120  if(!set) throw cms::Exception("InvalidConfiguration") << "MESet " << path << " could not be initialized";
121 
122  if(_MEParam.getUntrackedParameter<bool>("perLumi"))
123  set->setLumiFlag();
124 
125  return set;
126  }
127 
128  void
130  {
131  _desc.addUntracked<std::string>("path");
132  _desc.addUntracked<std::string>("kind");
133  _desc.addUntracked<std::string>("otype");
134  _desc.addUntracked<std::string>("btype");
135  _desc.addUntracked<std::string>("description");
136  _desc.addUntracked<bool>("online", false);
137  _desc.addUntracked<bool>("perLumi", false);
138  _desc.addOptionalUntracked<bool>("minutely");
139  _desc.addOptionalUntracked<bool>("cumulative");
140  _desc.addOptionalUntracked<bool>("shiftAxis");
141 
142  edm::ParameterSetDescription axisParameters;
143  binning::fillAxisDescriptions(axisParameters);
144  _desc.addOptionalUntracked("xaxis", axisParameters);
145  _desc.addOptionalUntracked("yaxis", axisParameters);
146  _desc.addOptionalUntracked("zaxis", axisParameters);
147 
148  edm::ParameterSetDescription multiParameters;
149  multiParameters.addWildcardUntracked<std::vector<std::string> >("*");
150  multiParameters.addWildcardUntracked<std::vector<int> >("*");
151  _desc.addOptionalUntracked("multi", multiParameters);
152  }
153 }
T getUntrackedParameter(std::string const &, T const &) const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
ParameterSet getUntrackedParameterSet(std::string const &name, ParameterSet const &defaultValue) const
ParameterWildcardBase * addWildcardUntracked(U const &pattern)
ObjectType translateObjectType(std::string const &)
std::map< std::string, std::vector< std::string > > ReplCandidates
Definition: MESetMulti.h:16
void fillAxisDescriptions(edm::ParameterSetDescription &)
MonitorElement::Kind translateKind(std::string const &)
ParameterDescriptionBase * addOptionalUntracked(U const &iLabel, T const &value)
BinningType translateBinningType(std::string const &)
AxisSpecs formAxis(edm::ParameterSet const &)
MESet * createMESet(edm::ParameterSet const &)
Definition: MESetUtils.cc:19
void fillMESetDescriptions(edm::ParameterSetDescription &)
Definition: MESetUtils.cc:129