CMS 3D CMS Logo

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