CMS 3D CMS Logo

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