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