CMS 3D CMS Logo

TkModuleGroupSelector.cc
Go to the documentation of this file.
1 
16 
17 #include <vector>
18 #include <map>
19 #include <set>
20 
21 //============================================================================
23  const edm::ParameterSet &cfg,
24  const std::vector<int> &sdets)
25  : nparameters_(0), subdetids_(sdets) {
26  //verify that all provided options are known
27  std::vector<std::string> parameterNames = cfg.getParameterNames();
28  for (std::vector<std::string>::const_iterator iParam = parameterNames.begin(); iParam != parameterNames.end();
29  ++iParam) {
30  const std::string name = (*iParam);
31  if (name != "RunRange" && name != "ReferenceRun" && name != "Granularity") {
32  throw cms::Exception("BadConfig") << "@SUB=TkModuleGroupSelector::TkModuleGroupSelector:"
33  << " Unknown parameter name '" << name << "' in PSet. Maybe a typo?";
34  }
35  }
36 
37  //extract the reference run range if defined
38  const edm::RunNumber_t defaultReferenceRun =
39  (cfg.exists("ReferenceRun") ? cfg.getParameter<edm::RunNumber_t>("ReferenceRun") : 0);
40 
41  //extract run range to be used for all module groups (if not locally overwritten)
42  const std::vector<edm::RunNumber_t> defaultRunRange =
43  (cfg.exists("RunRange") ? cfg.getParameter<std::vector<edm::RunNumber_t> >("RunRange")
44  : std::vector<edm::RunNumber_t>());
45 
46  // finally create everything from configuration
47  this->createModuleGroups(
48  aliTracker, cfg.getParameter<edm::VParameterSet>("Granularity"), defaultRunRange, defaultReferenceRun);
49 }
50 
51 //============================================================================
52 void TkModuleGroupSelector::fillDetIdMap(const unsigned int detid, const unsigned int groupid) {
53  //only add new entries
54  if (mapDetIdGroupId_.find(detid) == mapDetIdGroupId_.end()) {
55  mapDetIdGroupId_.insert(std::pair<unsigned int, unsigned int>(detid, groupid));
56  } else {
57  throw cms::Exception("BadConfig") << "@SUB=TkModuleGroupSelector:fillDetIdMap:"
58  << " Module with det ID " << detid << " configured in group " << groupid
59  << " but it was already selected"
60  << " in group " << mapDetIdGroupId_[detid] << ".";
61  }
62 }
63 
64 //============================================================================
66  bool split = false;
67  if (pset.exists("split")) {
68  split = pset.getParameter<bool>("split");
69  }
70  return split;
71 }
72 
73 //============================================================================
74 bool TkModuleGroupSelector::createGroup(unsigned int &Id,
75  const std::vector<edm::RunNumber_t> &range,
76  const std::list<Alignable *> &selected_alis,
77  const edm::RunNumber_t refrun) {
78  bool modules_selected = false;
79 
80  referenceRun_.push_back(refrun);
81  firstId_.push_back(Id);
82  runRange_.push_back(range);
83  for (std::list<Alignable *>::const_iterator it = selected_alis.begin(); it != selected_alis.end(); ++it) {
84  this->fillDetIdMap((*it)->id(), firstId_.size() - 1);
85  modules_selected = true;
86  }
87  if (refrun > 0 && !range.empty()) { //FIXME: last condition not really needed?
88  Id += range.size() - 1;
89  nparameters_ += range.size() - 1;
90  } else {
91  Id += range.size();
92  nparameters_ += range.size();
93  }
94 
95  if (refrun > 0 && range.front() > refrun) { //range.size() > 0 checked before
96  throw cms::Exception("BadConfig") << "@SUB=TkModuleGroupSelector::createGroup:\n"
97  << "Invalid combination of reference run number and specified run dependence"
98  << "\n in module group " << firstId_.size() << "."
99  << "\n Reference run number (" << refrun << ") is smaller than starting run "
100  << "\n number (" << range.front() << ") of first IOV.";
101  }
102  return modules_selected;
103 }
104 
105 //============================================================================
106 void TkModuleGroupSelector::verifyParameterNames(const edm::ParameterSet &pset, unsigned int psetnr) const {
107  std::vector<std::string> parameterNames = pset.getParameterNames();
108  for (std::vector<std::string>::const_iterator iParam = parameterNames.begin(); iParam != parameterNames.end();
109  ++iParam) {
110  const std::string name = (*iParam);
111  if (name != "levels" && name != "RunRange" && name != "split" && name != "ReferenceRun") {
112  throw cms::Exception("BadConfig") << "@SUB=TkModuleGroupSelector::verifyParameterNames:"
113  << " Unknown parameter name '" << name << "' in PSet number " << psetnr
114  << ". Maybe a typo?";
115  }
116  }
117 }
118 
119 //============================================================================
121  const edm::VParameterSet &granularityConfig,
122  const std::vector<edm::RunNumber_t> &defaultRunRange,
123  edm::RunNumber_t defaultReferenceRun) {
124  std::set<edm::RunNumber_t> localRunRange;
125  nparameters_ = 0;
126  unsigned int Id = 0;
127  unsigned int psetnr = 0;
128  //loop over all LA groups
129  for (edm::VParameterSet::const_iterator pset = granularityConfig.begin(); pset != granularityConfig.end(); ++pset) {
130  //test for unknown parameters
131  this->verifyParameterNames((*pset), psetnr);
132  psetnr++;
133 
134  bool modules_selected = false; //track whether at all a module has been selected in this group
135  const std::vector<edm::RunNumber_t> range =
136  ((*pset).exists("RunRange") ? pset->getParameter<std::vector<edm::RunNumber_t> >("RunRange") : defaultRunRange);
137  if (range.empty()) {
138  throw cms::Exception("BadConfig") << "@SUB=TkModuleGroupSelector::createModuleGroups:\n"
139  << "Run range array empty!";
140  }
141  const bool split = this->testSplitOption((*pset));
142 
143  edm::RunNumber_t refrun = 0;
144  if ((*pset).exists("ReferenceRun")) {
145  refrun = (*pset).getParameter<edm::RunNumber_t>("ReferenceRun");
146  } else {
147  refrun = defaultReferenceRun;
148  }
149 
150  AlignmentParameterSelector selector(aliTracker);
151  selector.clear();
152  selector.addSelections((*pset).getParameter<edm::ParameterSet>("levels"));
153 
154  const auto &alis = selector.selectedAlignables();
155  std::list<Alignable *> selected_alis;
156  for (const auto &it : alis) {
157  const auto &aliDaughts = it->deepComponents();
158  for (const auto &iD : aliDaughts) {
159  if (iD->alignableObjectId() == align::AlignableDetUnit || iD->alignableObjectId() == align::AlignableDet) {
160  if (split) {
161  modules_selected = this->createGroup(Id, range, std::list<Alignable *>(1, iD), refrun);
162  } else {
163  selected_alis.push_back(iD);
164  }
165  }
166  }
167  }
168 
169  if (!split) {
170  modules_selected = this->createGroup(Id, range, selected_alis, refrun);
171  }
172 
174  for (std::vector<edm::RunNumber_t>::const_iterator iRun = range.begin(); iRun != range.end(); ++iRun) {
175  localRunRange.insert((*iRun));
176  if ((*iRun) > firstRun) {
177  firstRun = (*iRun);
178  } else {
179  throw cms::Exception("BadConfig") << "@SUB=TkModuleGroupSelector::createModuleGroups:"
180  << " Run range not sorted.";
181  }
182  }
183 
184  if (!modules_selected) {
185  throw cms::Exception("BadConfig") << "@SUB=TkModuleGroupSelector:createModuleGroups:"
186  << " No module was selected in the module group selector in group "
187  << (firstId_.size() - 1) << ".";
188  }
189  }
190 
191  //copy local set into the global vector of run boundaries
192  for (std::set<edm::RunNumber_t>::const_iterator itRun = localRunRange.begin(); itRun != localRunRange.end();
193  ++itRun) {
194  globalRunRange_.push_back((*itRun));
195  }
196 }
197 
198 //============================================================================
200 
201 //============================================================================
202 unsigned int TkModuleGroupSelector::numIovs() const { return globalRunRange_.size(); }
203 
204 //============================================================================
206  return iovNum < this->numIovs() ? globalRunRange_.at(iovNum) : 0;
207 }
208 
209 //======================================================================
211  // Return the index of the parameter that is used for this DetId.
212  // If this DetId is not treated, return values < 0.
213 
214  const DetId temp_id(detId);
215 
216  int index = -1;
217 
218  bool sel = false;
219  for (std::vector<int>::const_iterator itSubDets = subdetids_.begin(); itSubDets != subdetids_.end(); ++itSubDets) {
220  if (temp_id.det() == DetId::Tracker && temp_id.subdetId() == (*itSubDets)) {
221  sel = true;
222  break;
223  }
224  }
225 
226  if (temp_id.det() != DetId::Tracker || !sel)
227  return -1;
228 
229  std::map<unsigned int, unsigned int>::const_iterator it = mapDetIdGroupId_.find(detId);
230  if (it != mapDetIdGroupId_.end()) {
231  const unsigned int iAlignableGroup = (*it).second;
232  const std::vector<edm::RunNumber_t> &runs = runRange_.at(iAlignableGroup);
233  const unsigned int id0 = firstId_.at(iAlignableGroup);
234  const edm::RunNumber_t refrun = referenceRun_.at(iAlignableGroup);
235 
236  unsigned int iovNum = 0;
237  for (; iovNum < runs.size(); ++iovNum) {
238  if (runs[iovNum] > run)
239  break;
240  }
241  if (iovNum == 0) {
242  throw cms::Exception("BadConfig") << "@SUB=TkModuleGroupSelector::getParameterIndexFromDetId:\n"
243  << "Run " << run << " not foreseen for detid '" << detId << "'"
244  << " in module group " << iAlignableGroup << ".";
245  } else {
246  --iovNum;
247  }
248 
249  //test whether the iov contains the reference run
250  if (refrun > 0) { //if > 0 a reference run number has been provided
251  if (iovNum + 1 == runs.size()) {
252  if (refrun >= runs[iovNum])
253  return -1;
254  } else if ((iovNum + 1) < runs.size()) {
255  if (refrun >= runs[iovNum] && refrun < runs[iovNum + 1]) {
256  return -1;
257  }
258  }
259  if (run > refrun) {
260  //iovNum > 0 due to checks in createGroup(...) and createModuleGroups(...)
261  //remove IOV in which the reference run can be found
262  iovNum -= 1;
263  }
264  }
265 
266  index = id0 + iovNum;
267  }
268  return index;
269 }
T getParameter(std::string const &) const
std::vector< int > subdetids_
std::vector< edm::RunNumber_t > globalRunRange_
edm::RunNumber_t firstRunOfIOV(unsigned int iovNum) const
First run of iov (0 if iovNum not treated).
std::vector< edm::RunNumber_t > referenceRun_
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:33
void fillDetIdMap(const unsigned int detid, const unsigned int groupid)
bool exists(std::string const &parameterName) const
checks if a parameter exists
int getParameterIndexFromDetId(unsigned int detId, edm::RunNumber_t run) const
std::map< unsigned int, unsigned int > mapDetIdGroupId_
void clear()
remove all selected Alignables and geometrical restrictions
void verifyParameterNames(const edm::ParameterSet &pset, unsigned int psetnr) const
void createModuleGroups(AlignableTracker *aliTracker, const edm::VParameterSet &granularityConfig, const std::vector< edm::RunNumber_t > &defaultRunRange, edm::RunNumber_t defaultReferenceRun)
std::vector< std::vector< edm::RunNumber_t > > runRange_
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
std::vector< std::string > getParameterNames() const
const align::Alignables & selectedAlignables() const
vector of alignables selected so far
Definition: DetId.h:18
TkModuleGroupSelector(AlignableTracker *aliTracker, const edm::ParameterSet &cfg, const std::vector< int > &sdets)
Constructor.
const bool testSplitOption(const edm::ParameterSet &pset) const
bool createGroup(unsigned int &Id, const std::vector< edm::RunNumber_t > &range, const std::list< Alignable * > &selected_alis, const edm::RunNumber_t refrun)
unsigned int getNumberOfParameters() const
unsigned int numIovs() const
Total number of IOVs.
unsigned int addSelections(const edm::ParameterSet &pSet)
firstRun
Definition: dataset.py:940
unsigned int RunNumber_t
std::vector< unsigned int > firstId_
double split
Definition: MVATrainer.cc:139
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:39