CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MillePedeAlignmentAlgorithm.cc
Go to the documentation of this file.
1 
12 
17 
19 
24 #include "Alignment/MillePedeAlignmentAlgorithm/src/Mille.h" // 'unpublished' interface located in src
29 
32 
36 
40 
41 // includes to make known that they inherit from Alignable:
45 
50 
56 
60 
62 
63 #include <fstream>
64 #include <sstream>
65 #include <algorithm>
66 #include <sys/stat.h>
67 
68 #include <TMath.h>
72 
73 // Includes for PXB survey
80 
82 
83 using namespace gbl;
84 
85 // Constructor ----------------------------------------------------------------
86 //____________________________________________________
88  : AlignmentAlgorithmBase(cfg, iC),
89  topoToken_(iC.esConsumes<TrackerTopology, TrackerTopologyRcd, edm::Transition::BeginRun>()),
90  aliThrToken_(iC.esConsumes<AlignPCLThresholds, AlignPCLThresholdsRcd, edm::Transition::BeginRun>()),
91  theConfig(cfg),
92  theMode(this->decodeMode(theConfig.getUntrackedParameter<std::string>("mode"))),
93  theDir(theConfig.getUntrackedParameter<std::string>("fileDir")),
94  theAlignmentParameterStore(nullptr),
95  theAlignables(),
96  theTrajectoryFactory(
97  TrajectoryFactoryPlugin::get()->create(theConfig.getParameter<edm::ParameterSet>("TrajectoryFactory")
98  .getParameter<std::string>("TrajectoryFactoryName"),
99  theConfig.getParameter<edm::ParameterSet>("TrajectoryFactory"),
100  iC)),
101  theMinNumHits(cfg.getParameter<unsigned int>("minNumHits")),
102  theMaximalCor2D(cfg.getParameter<double>("max2Dcorrelation")),
103  firstIOV_(cfg.getUntrackedParameter<AlignmentAlgorithmBase::RunNumber>("firstIOV")),
104  ignoreFirstIOVCheck_(cfg.getUntrackedParameter<bool>("ignoreFirstIOVCheck")),
105  enableAlignableUpdates_(cfg.getUntrackedParameter<bool>("enableAlignableUpdates")),
106  theLastWrittenIov(0),
107  theGblDoubleBinary(cfg.getParameter<bool>("doubleBinary")),
108  runAtPCL_(cfg.getParameter<bool>("runAtPCL")),
109  ignoreHitsWithoutGlobalDerivatives_(cfg.getParameter<bool>("ignoreHitsWithoutGlobalDerivatives")),
110  skipGlobalPositionRcdCheck_(cfg.getParameter<bool>("skipGlobalPositionRcdCheck")),
111  uniqueRunRanges_(align::makeUniqueRunRanges(cfg.getUntrackedParameter<edm::VParameterSet>("RunRangeSelection"),
112  cond::timeTypeSpecs[cond::runnumber].beginValue)),
113  enforceSingleIOVInput_(!(enableAlignableUpdates_ && areIOVsSpecified())),
114  lastProcessedRun_(cond::timeTypeSpecs[cond::runnumber].beginValue) {
115  if (!theDir.empty() && theDir.find_last_of('/') != theDir.size() - 1)
116  theDir += '/'; // may need '/'
117  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm"
118  << "Start in mode '" << theConfig.getUntrackedParameter<std::string>("mode")
119  << "' with output directory '" << theDir << "'.";
120  if (this->isMode(myMilleBit)) {
121  theMille = std::make_unique<Mille>(
122  (theDir + theConfig.getParameter<std::string>("binaryFile")).c_str()); // add ', false);' for text output);
123  // use same file for GBL
124  theBinary = std::make_unique<MilleBinary>((theDir + theConfig.getParameter<std::string>("binaryFile")).c_str(),
126  }
127 }
128 
129 // Destructor ----------------------------------------------------------------
130 //____________________________________________________
132 
133 // Call at beginning of job ---------------------------------------------------
134 //____________________________________________________
136  AlignableTracker *tracker,
138  AlignableExtras *extras,
139  AlignmentParameterStore *store) {
140  if (muon) {
141  edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::initialize"
142  << "Running with AlignabeMuon not yet tested.";
143  }
144 
148  const auto &iov_alignments = setup.get<TrackerAlignmentRcd>().validityInterval();
149  const auto &iov_surfaces = setup.get<TrackerSurfaceDeformationRcd>().validityInterval();
150  const auto &iov_errors = setup.get<TrackerAlignmentErrorExtendedRcd>().validityInterval();
151 
152  std::ostringstream message;
153  bool throwException{false};
154  if (iov_alignments.first().eventID().run() != MIN_VAL || iov_alignments.last().eventID().run() != MAX_VAL) {
155  message << "\nTrying to apply " << setup.get<TrackerAlignmentRcd>().key().name()
156  << " with multiple IOVs in tag without specifying 'RunRangeSelection'.\n"
157  << "Validity range is " << iov_alignments.first().eventID().run() << " - "
158  << iov_alignments.last().eventID().run() << "\n";
159  throwException = true;
160  }
161  if (iov_surfaces.first().eventID().run() != MIN_VAL || iov_surfaces.last().eventID().run() != MAX_VAL) {
162  message << "\nTrying to apply " << setup.get<TrackerSurfaceDeformationRcd>().key().name()
163  << " with multiple IOVs in tag without specifying 'RunRangeSelection'.\n"
164  << "Validity range is " << iov_surfaces.first().eventID().run() << " - "
165  << iov_surfaces.last().eventID().run() << "\n";
166  throwException = true;
167  }
168  if (iov_errors.first().eventID().run() != MIN_VAL || iov_errors.last().eventID().run() != MAX_VAL) {
169  message << "\nTrying to apply " << setup.get<TrackerAlignmentErrorExtendedRcd>().key().name()
170  << " with multiple IOVs in tag without specifying 'RunRangeSelection'.\n"
171  << "Validity range is " << iov_errors.first().eventID().run() << " - "
172  << iov_errors.last().eventID().run() << "\n";
173  throwException = true;
174  }
175  if (throwException) {
176  throw cms::Exception("DatabaseError") << "@SUB=MillePedeAlignmentAlgorithm::initialize" << message.str();
177  }
178  }
179 
180  //Retrieve tracker topology from geometry
181  const TrackerTopology *const tTopo = &setup.getData(topoToken_);
182 
183  //Retrieve the thresolds cuts from DB for the PCL
184  if (runAtPCL_) {
185  const auto &th = &setup.getData(aliThrToken_);
186  theThresholds = std::make_shared<AlignPCLThresholds>();
187  storeThresholds(th->getNrecords(), th->getThreshold_Map());
188  }
189 
190  theAlignableNavigator = std::make_unique<AlignableNavigator>(extras, tracker, muon);
192  theAlignables = theAlignmentParameterStore->alignables();
193 
194  edm::ParameterSet pedeLabelerCfg(theConfig.getParameter<edm::ParameterSet>("pedeLabeler"));
195  edm::VParameterSet RunRangeSelectionVPSet(theConfig.getUntrackedParameter<edm::VParameterSet>("RunRangeSelection"));
196  pedeLabelerCfg.addUntrackedParameter<edm::VParameterSet>("RunRangeSelection", RunRangeSelectionVPSet);
197 
198  std::string labelerPlugin = "PedeLabeler";
199  if (!RunRangeSelectionVPSet.empty()) {
200  labelerPlugin = "RunRangeDependentPedeLabeler";
201  if (pedeLabelerCfg.exists("plugin")) {
202  std::string labelerPluginCfg = pedeLabelerCfg.getParameter<std::string>("plugin");
203  if ((labelerPluginCfg != "PedeLabeler" && labelerPluginCfg != "RunRangeDependentPedeLabeler") ||
204  !pedeLabelerCfg.getUntrackedParameter<edm::VParameterSet>("parameterInstances").empty()) {
205  throw cms::Exception("BadConfig") << "MillePedeAlignmentAlgorithm::initialize"
206  << "both RunRangeSelection and generic labeler specified in config file. "
207  << "Please get rid of either one of them.\n";
208  }
209  }
210  } else {
211  if (pedeLabelerCfg.exists("plugin")) {
212  labelerPlugin = pedeLabelerCfg.getParameter<std::string>("plugin");
213  }
214  }
215 
216  if (!pedeLabelerCfg.exists("plugin")) {
217  pedeLabelerCfg.addUntrackedParameter<std::string>("plugin", labelerPlugin);
218  }
219 
220  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::initialize"
221  << "Using plugin '" << labelerPlugin << "' to generate labels.";
222 
223  thePedeLabels = std::shared_ptr<PedeLabelerBase>(PedeLabelerPluginFactory::get()->create(
224  labelerPlugin, PedeLabelerBase::TopLevelAlignables(tracker, muon, extras), pedeLabelerCfg));
225 
226  // 1) Create PedeSteerer: correct alignable positions for coordinate system selection
227  edm::ParameterSet pedeSteerCfg(theConfig.getParameter<edm::ParameterSet>("pedeSteerer"));
228  thePedeSteer = std::make_unique<PedeSteerer>(tracker,
229  muon,
230  extras,
232  thePedeLabels.get(),
233  pedeSteerCfg,
234  theDir,
235  !this->isMode(myPedeSteerBit));
236 
237  // 2) If requested, directly read in and apply result of previous pede run,
238  // assuming that correction from 1) was also applied to create the result:
239  const std::vector<edm::ParameterSet> mprespset(
240  theConfig.getParameter<std::vector<edm::ParameterSet> >("pedeReaderInputs"));
241  if (!mprespset.empty()) {
242  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::initialize"
243  << "Apply " << mprespset.end() - mprespset.begin()
244  << " previous MillePede constants from 'pedeReaderInputs'.";
245  }
246 
247  // FIXME: add selection of run range via 'pedeReaderInputs'
248  // Note: Results for parameters of IntegratedCalibration's cannot be treated...
250  for (std::vector<edm::ParameterSet>::const_iterator iSet = mprespset.begin(), iE = mprespset.end(); iSet != iE;
251  ++iSet) {
252  // This read will ignore calibrations as long as they are not yet passed to Millepede
253  // during/before initialize(..) - currently addCalibrations(..) is called later in AlignmentProducer
254  if (!this->readFromPede((*iSet), false, runrange)) { // false: do not erase SelectionUserVariables
255  throw cms::Exception("BadConfig")
256  << "MillePedeAlignmentAlgorithm::initialize: Problems reading input constants of "
257  << "pedeReaderInputs entry " << iSet - mprespset.begin() << '.';
258  }
259  theAlignmentParameterStore->applyParameters();
260  // Needed to shut up later warning from checkAliParams:
261  theAlignmentParameterStore->resetParameters();
262  }
263 
264  // 3) Now create steerings with 'final' start position:
265  thePedeSteer->buildSubSteer(tracker, muon, extras);
266 
267  // After (!) 1-3 of PedeSteerer which uses the SelectionUserVariables attached to the parameters:
268  this->buildUserVariables(theAlignables); // for hit statistics and/or pede result
269 
270  if (this->isMode(myMilleBit)) {
271  if (!theConfig.getParameter<std::vector<std::string> >("mergeBinaryFiles").empty() ||
272  !theConfig.getParameter<std::vector<std::string> >("mergeTreeFiles").empty()) {
273  throw cms::Exception("BadConfig") << "'vstring mergeTreeFiles' and 'vstring mergeBinaryFiles' must be empty for "
274  << "modes running mille.";
275  }
276  const std::string moniFile(theConfig.getUntrackedParameter<std::string>("monitorFile"));
277  if (!moniFile.empty())
278  theMonitor = std::make_unique<MillePedeMonitor>(tTopo, (theDir + moniFile).c_str());
279 
280  // Get trajectory factory. In case nothing found, FrameWork will throw...
281  }
282 
283  if (this->isMode(myPedeSteerBit)) {
284  // Get config for survey and set flag accordingly
285  const edm::ParameterSet pxbSurveyCfg(theConfig.getParameter<edm::ParameterSet>("surveyPixelBarrel"));
286  theDoSurveyPixelBarrel = pxbSurveyCfg.getParameter<bool>("doSurvey");
288  this->addPxbSurvey(pxbSurveyCfg);
289  }
290 }
291 
292 //____________________________________________________
294 
295 //____________________________________________________
296 bool MillePedeAlignmentAlgorithm::addCalibrations(const std::vector<IntegratedCalibrationBase *> &iCals) {
297  theCalibrations.insert(theCalibrations.end(), iCals.begin(), iCals.end());
298  thePedeLabels->addCalibrations(iCals);
299  return true;
300 }
301 
302 //____________________________________________________
304  const AlignPCLThresholds::threshold_map &thresholdMap) {
305  theThresholds->setAlignPCLThresholds(nRecords, thresholdMap);
306  return true;
307 }
308 
309 //_____________________________________________________________________________
311  if (isMode(myMilleBit)) {
312  return true;
313  } else {
314  return false;
315  }
316 }
317 
318 //_____________________________________________________________________________
320  if (isMode(myPedeReadBit)) {
321  if (runAtPCL_) {
322  MillePedeFileReader mpReader(
324  mpReader.read();
325  return mpReader.storeAlignments();
326  } else {
327  return true;
328  }
329  } else {
330  return false;
331  }
332 }
333 
334 //____________________________________________________
336  if (this->isMode(myPedeReadBit)) {
337  if (not theAlignmentParameterStore) {
338  return false;
339  }
340  // restore initial positions, rotations and deformations
342  theAlignmentParameterStore->restoreCachedTransformations(runrange.first);
343  } else {
344  theAlignmentParameterStore->restoreCachedTransformations();
345  }
346 
347  // Needed to shut up later warning from checkAliParams:
348  theAlignmentParameterStore->resetParameters();
349  // To avoid that they keep values from previous IOV if no new one in pede result
351 
352  if (!this->readFromPede(theConfig.getParameter<edm::ParameterSet>("pedeReader"), true, runrange)) {
353  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::setParametersForRunRange"
354  << "Problems reading pede result, but applying!";
355  }
356  theAlignmentParameterStore->applyParameters();
357 
358  this->doIO(++theLastWrittenIov); // pre-increment!
359  }
360 
361  return true;
362 }
363 
364 // Call at end of job ---------------------------------------------------------
365 //____________________________________________________
368  theMille.reset(); // delete to close binary before running pede below (flush would be enough...)
369  theBinary.reset();
370 
371  std::vector<std::string> files;
372  if (this->isMode(myMilleBit) || !theConfig.getParameter<std::string>("binaryFile").empty()) {
373  files.push_back(theDir + theConfig.getParameter<std::string>("binaryFile"));
374  } else {
375  const std::vector<std::string> plainFiles(theConfig.getParameter<std::vector<std::string> >("mergeBinaryFiles"));
376  files = getExistingFormattedFiles(plainFiles, theDir);
377  // Do some logging:
378  std::string filesForLogOutput;
379  for (const auto &file : files)
380  filesForLogOutput += " " + file + ",";
381  if (filesForLogOutput.length() != 0)
382  filesForLogOutput.pop_back();
383  edm::LogInfo("Alignment") << "Based on the config parameter mergeBinaryFiles, using the following "
384  << "files as input (assigned weights are indicated by ' -- <weight>'):"
385  << filesForLogOutput;
386  }
387 
389  return;
390 
391  // cache all positions, rotations and deformations
392  theAlignmentParameterStore->cacheTransformations();
394  if (lastProcessedRun_ < uniqueRunRanges_.back().first) {
395  throw cms::Exception("BadConfig") << "@SUB=MillePedeAlignmentAlgorithm::terminate\n"
396  << "Last IOV of 'RunRangeSelection' has not been processed. "
397  << "Please reconfigure your source to process the runs at least up to "
398  << uniqueRunRanges_.back().first << ".";
399  }
400  auto lastCachedRun = uniqueRunRanges_.front().first;
401  for (const auto &runRange : uniqueRunRanges_) {
402  const auto run = runRange.first;
403  if (std::find(cachedRuns_.begin(), cachedRuns_.end(), run) == cachedRuns_.end()) {
404  theAlignmentParameterStore->restoreCachedTransformations(lastCachedRun);
405  theAlignmentParameterStore->cacheTransformations(run);
406  } else {
407  lastCachedRun = run;
408  }
409  }
410  theAlignmentParameterStore->restoreCachedTransformations();
411  }
412 
413  const std::string masterSteer(thePedeSteer->buildMasterSteer(files)); // do only if myPedeSteerBit?
414  if (this->isMode(myPedeRunBit)) {
415  thePedeSteer->runPede(masterSteer);
416  }
417 
418  // parameters from pede are not yet applied,
419  // so we can still write start positions (but with hit statistics in case of mille):
420  this->doIO(0);
421  theLastWrittenIov = 0;
422 }
423 
425  const std::vector<std::string> &plainFiles, const std::string &theDir) {
426  std::vector<std::string> files;
427  for (const auto &plainFile : plainFiles) {
428  const std::string &theInputFileName = plainFile;
429  int theNumber = 0;
430  while (true) {
431  // Create a formatted version of the filename, with growing numbers
432  // If the parameter doesn't contain a formatting directive, it just stays unchanged
433  char theNumberedInputFileName[200];
434  sprintf(theNumberedInputFileName, theInputFileName.c_str(), theNumber);
435  std::string theCompleteInputFileName = theDir + theNumberedInputFileName;
436  const auto endOfStrippedFileName = theCompleteInputFileName.rfind(" --");
437  const auto strippedInputFileName = theCompleteInputFileName.substr(0, endOfStrippedFileName);
438  // Check if the file exists
439  struct stat buffer;
440  if (stat(strippedInputFileName.c_str(), &buffer) == 0) {
441  // If the file exists, add it to the list
442  files.push_back(theCompleteInputFileName);
443  if (theNumberedInputFileName == theInputFileName) {
444  // If the filename didn't contain a formatting directive, no reason to look any further, break out of the loop
445  break;
446  } else {
447  // Otherwise look for the next number
448  theNumber++;
449  }
450  } else {
451  // The file doesn't exist, break out of the loop
452  break;
453  }
454  }
455  // warning if unformatted (-> theNumber stays at 0) does not exist
456  if (theNumber == 0 && (files.empty() || files.back() != plainFile)) {
457  edm::LogWarning("Alignment") << "The input file '" << plainFile << "' does not exist.";
458  }
459  }
460  return files;
461 }
462 
463 // Run the algorithm on trajectories and tracks -------------------------------
464 //____________________________________________________
466  if (!this->isMode(myMilleBit))
467  return; // no theMille created...
468  const auto &tracks = eventInfo.trajTrackPairs();
469 
470  if (theMonitor) { // monitor input tracks
471  for (const auto &iTrajTrack : tracks) {
472  theMonitor->fillTrack(iTrajTrack.second);
473  }
474  }
475 
476  const RefTrajColl trajectories(theTrajectoryFactory->trajectories(setup, tracks, eventInfo.beamSpot()));
477 
478  // Now loop over ReferenceTrajectoryCollection
479  unsigned int refTrajCount = 0; // counter for track monitoring
480  const auto tracksPerTraj = theTrajectoryFactory->tracksPerTrajectory();
481  for (auto iRefTraj = trajectories.cbegin(), iRefTrajE = trajectories.cend(); iRefTraj != iRefTrajE;
482  ++iRefTraj, ++refTrajCount) {
483  RefTrajColl::value_type refTrajPtr = *iRefTraj;
484  if (theMonitor)
485  theMonitor->fillRefTrajectory(refTrajPtr);
486 
487  const auto nHitXy = this->addReferenceTrajectory(setup, eventInfo, refTrajPtr);
488 
489  if (theMonitor && (nHitXy.first || nHitXy.second)) {
490  // if trajectory used (i.e. some hits), fill monitoring
491  const auto offset = tracksPerTraj * refTrajCount;
492  for (unsigned int iTrack = 0; iTrack < tracksPerTraj; ++iTrack) {
493  theMonitor->fillUsedTrack(tracks[offset + iTrack].second, nHitXy.first, nHitXy.second);
494  }
495  }
496 
497  } // end of reference trajectory and track loop
498 }
499 
500 //____________________________________________________
501 std::pair<unsigned int, unsigned int> MillePedeAlignmentAlgorithm::addReferenceTrajectory(
502  const edm::EventSetup &setup, const EventInfo &eventInfo, const RefTrajColl::value_type &refTrajPtr) {
503  std::pair<unsigned int, unsigned int> hitResultXy(0, 0);
504  if (refTrajPtr->isValid()) {
505  // GblTrajectory?
506  if (!refTrajPtr->gblInput().empty()) {
507  // by construction: number of GblPoints == number of recHits or == zero !!!
508  unsigned int iHit = 0;
509  unsigned int numPointsWithMeas = 0;
510  std::vector<GblPoint>::iterator itPoint;
511  auto theGblInput = refTrajPtr->gblInput();
512  for (unsigned int iTraj = 0; iTraj < refTrajPtr->gblInput().size(); ++iTraj) {
513  for (itPoint = refTrajPtr->gblInput()[iTraj].first.begin(); itPoint < refTrajPtr->gblInput()[iTraj].first.end();
514  ++itPoint) {
515  if (this->addGlobalData(setup, eventInfo, refTrajPtr, iHit++, *itPoint) < 0)
516  return hitResultXy;
517  if (itPoint->numMeasurements() >= 1)
518  ++numPointsWithMeas;
519  }
520  }
521  hitResultXy.first = numPointsWithMeas;
522  // check #hits criterion
523  if (hitResultXy.first == 0 || hitResultXy.first < theMinNumHits)
524  return hitResultXy;
525  // construct GBL trajectory
526  if (refTrajPtr->gblInput().size() == 1) {
527  // from single track
528  GblTrajectory aGblTrajectory(refTrajPtr->gblInput()[0].first, refTrajPtr->nominalField() != 0);
529  // GBL fit trajectory
530  /*double Chi2;
531  int Ndf;
532  double lostWeight;
533  aGblTrajectory.fit(Chi2, Ndf, lostWeight);
534  std::cout << " GblFit: " << Chi2 << ", " << Ndf << ", " << lostWeight << std::endl; */
535  // write to MP binary file
536  if (aGblTrajectory.isValid() && aGblTrajectory.getNumPoints() >= theMinNumHits)
537  aGblTrajectory.milleOut(*theBinary);
538  }
539  if (refTrajPtr->gblInput().size() == 2) {
540  // from TwoBodyDecay
541  GblTrajectory aGblTrajectory(refTrajPtr->gblInput(),
542  refTrajPtr->gblExtDerivatives(),
543  refTrajPtr->gblExtMeasurements(),
544  refTrajPtr->gblExtPrecisions());
545  // write to MP binary file
546  if (aGblTrajectory.isValid() && aGblTrajectory.getNumPoints() >= theMinNumHits)
547  aGblTrajectory.milleOut(*theBinary);
548  }
549  } else {
550  // to add hits if all fine:
551  std::vector<AlignmentParameters *> parVec(refTrajPtr->recHits().size());
552  // collect hit statistics, assuming that there are no y-only hits
553  std::vector<bool> validHitVecY(refTrajPtr->recHits().size(), false);
554  // Use recHits from ReferenceTrajectory (since they have the right order!):
555  for (unsigned int iHit = 0; iHit < refTrajPtr->recHits().size(); ++iHit) {
556  const int flagXY = this->addMeasurementData(setup, eventInfo, refTrajPtr, iHit, parVec[iHit]);
557 
558  if (flagXY < 0) { // problem
559  hitResultXy.first = 0;
560  break;
561  } else { // hit is fine, increase x/y statistics
562  if (flagXY >= 1)
563  ++hitResultXy.first;
564  validHitVecY[iHit] = (flagXY >= 2);
565  }
566  } // end loop on hits
567 
568  // add virtual measurements
569  for (unsigned int iVirtualMeas = 0; iVirtualMeas < refTrajPtr->numberOfVirtualMeas(); ++iVirtualMeas) {
570  this->addVirtualMeas(refTrajPtr, iVirtualMeas);
571  }
572 
573  // kill or end 'track' for mille, depends on #hits criterion
574  if (hitResultXy.first == 0 || hitResultXy.first < theMinNumHits) {
575  theMille->kill();
576  hitResultXy.first = hitResultXy.second = 0; //reset
577  } else {
578  theMille->end();
579  // add x/y hit count to MillePedeVariables of parVec,
580  // returning number of y-hits of the reference trajectory
581  hitResultXy.second = this->addHitCount(parVec, validHitVecY);
582  //
583  }
584  }
585 
586  } // end if valid trajectory
587 
588  return hitResultXy;
589 }
590 
591 //____________________________________________________
592 unsigned int MillePedeAlignmentAlgorithm::addHitCount(const std::vector<AlignmentParameters *> &parVec,
593  const std::vector<bool> &validHitVecY) const {
594  // Loop on all hit information in the input arrays and count valid y-hits:
595  unsigned int nHitY = 0;
596  for (unsigned int iHit = 0; iHit < validHitVecY.size(); ++iHit) {
597  Alignable *ali = (parVec[iHit] ? parVec[iHit]->alignable() : nullptr);
598  // Loop upwards on hierarchy of alignables to add hits to all levels
599  // that are currently aligned. If only a non-selected alignable was hit,
600  // (i.e. flagXY == 0 in addReferenceTrajectory(..)), there is no loop at all...
601  while (ali) {
603  if (pars) { // otherwise hierarchy level not selected
604  // cast ensured by previous checks:
605  MillePedeVariables *mpVar = static_cast<MillePedeVariables *>(pars->userVariables());
606  // every hit has an x-measurement, cf. addReferenceTrajectory(..):
607  mpVar->increaseHitsX();
608  if (validHitVecY[iHit]) {
609  mpVar->increaseHitsY();
610  if (pars == parVec[iHit])
611  ++nHitY; // do not count hits twice
612  }
613  }
614  ali = ali->mother();
615  }
616  }
617 
618  return nHitY;
619 }
620 
622  if (run.run() < firstIOV_ && !ignoreFirstIOVCheck_) {
623  throw cms::Exception("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::beginRun\n"
624  << "Using data (run = " << run.run() << ") prior to the first defined IOV ("
625  << firstIOV_ << ").";
626  }
627 
628  lastProcessedRun_ = run.run();
629 
630  if (changed && enableAlignableUpdates_) {
631  const auto runNumber = run.run();
633  for (auto runRange = uniqueRunRanges_.crbegin(); runRange != uniqueRunRanges_.crend(); ++runRange) {
634  if (runNumber >= runRange->first) {
635  firstRun = runRange->first;
636  break;
637  }
638  }
639  if (std::find(cachedRuns_.begin(), cachedRuns_.end(), firstRun) != cachedRuns_.end()) {
640  const auto &geometryRcd = setup.get<IdealGeometryRecord>();
641  const auto &globalPosRcd = setup.get<GlobalPositionRcd>();
642  const auto &alignmentRcd = setup.get<TrackerAlignmentRcd>();
643  const auto &surfaceRcd = setup.get<TrackerSurfaceDeformationRcd>();
644  const auto &errorRcd = setup.get<TrackerAlignmentErrorExtendedRcd>();
645 
646  std::ostringstream message;
647  bool throwException{false};
648  message << "Trying to cache tracker alignment payloads for a run (" << runNumber << ") in an IOV (" << firstRun
649  << ") that was already cached.\n"
650  << "The following records in your input database tag have an IOV "
651  << "boundary that does not match your IOV definition:\n";
652  if (geometryRcd.validityInterval().first().eventID().run() > firstRun) {
653  message << " - IdealGeometryRecord '" << geometryRcd.key().name() << "' (since "
654  << geometryRcd.validityInterval().first().eventID().run() << ")\n";
655  throwException = true;
656  }
657  if (globalPosRcd.validityInterval().first().eventID().run() > firstRun) {
658  message << " - GlobalPositionRecord '" << globalPosRcd.key().name() << "' (since "
659  << globalPosRcd.validityInterval().first().eventID().run() << ")";
661  message << " --> ignored\n";
662  } else {
663  message << "\n";
664  throwException = true;
665  }
666  }
667  if (alignmentRcd.validityInterval().first().eventID().run() > firstRun) {
668  message << " - TrackerAlignmentRcd '" << alignmentRcd.key().name() << "' (since "
669  << alignmentRcd.validityInterval().first().eventID().run() << ")\n";
670  throwException = true;
671  }
672  if (surfaceRcd.validityInterval().first().eventID().run() > firstRun) {
673  message << " - TrackerSurfaceDeformationRcd '" << surfaceRcd.key().name() << "' (since "
674  << surfaceRcd.validityInterval().first().eventID().run() << ")\n";
675  throwException = true;
676  }
677  if (errorRcd.validityInterval().first().eventID().run() > firstRun) {
678  message << " - TrackerAlignmentErrorExtendedRcd '" << errorRcd.key().name() << "' (since "
679  << errorRcd.validityInterval().first().eventID().run() << ")\n";
680  throwException = true;
681  }
682 
683  if (throwException) {
684  throw cms::Exception("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::beginRun\n" << message.str();
685  }
686  } else {
687  cachedRuns_.push_back(firstRun);
689  }
690  }
691 }
692 
693 //____________________________________________________
695  const EndRunInfo &runInfo,
696  const edm::EventSetup &setup) {
697  if (runInfo.tkLasBeams() && runInfo.tkLasBeamTsoses()) {
698  // LAS beam treatment
699  this->addLaserData(eventInfo, *(runInfo.tkLasBeams()), *(runInfo.tkLasBeamTsoses()));
700  }
701  if (this->isMode(myMilleBit))
702  theMille->flushOutputFile();
703 }
704 
705 // Implementation of endRun that DOES get called. (Because we need it.)
707  if (this->isMode(myMilleBit))
708  theMille->flushOutputFile();
709 }
710 
711 //____________________________________________________
713  if (!runAtPCL_)
714  return;
715  if (this->isMode(myMilleBit))
716  theMille->resetOutputFile();
717 }
718 
719 //____________________________________________________
721  if (!runAtPCL_)
722  return;
723  if (this->isMode(myMilleBit))
724  theMille->flushOutputFile();
725 }
726 
727 //____________________________________________________
729  const EventInfo &eventInfo,
731  unsigned int iHit,
733  params = nullptr;
734  theFloatBufferX.clear();
735  theFloatBufferY.clear();
736  theIntBuffer.clear();
737 
738  const TrajectoryStateOnSurface &tsos = refTrajPtr->trajectoryStates()[iHit];
739  const ConstRecHitPointer &recHitPtr = refTrajPtr->recHits()[iHit];
740  // ignore invalid hits
741  if (!recHitPtr->isValid())
742  return 0;
743 
744  // First add the derivatives from IntegratedCalibration's,
745  // should even be OK if problems for "usual" derivatives from Alignables
746  this->globalDerivativesCalibration(recHitPtr,
747  tsos,
748  setup,
749  eventInfo, // input
752  theIntBuffer); // output
753 
754  // get AlignableDet/Unit for this hit
755  AlignableDetOrUnitPtr alidet(theAlignableNavigator->alignableFromDetId(recHitPtr->geographicalId()));
756 
757  if (!this->globalDerivativesHierarchy(eventInfo,
758  tsos,
759  alidet,
760  alidet,
761  theFloatBufferX, // 2x alidet, sic!
763  theIntBuffer,
764  params)) {
765  return -1; // problem
767  return 0; // empty for X: no alignable for hit, nor calibrations
768  } else {
769  // store measurement even if no alignable or calibrations
770  // -> measurement used for pede-internal track-fit
771  return this->callMille(refTrajPtr, iHit, theIntBuffer, theFloatBufferX, theFloatBufferY);
772  }
773 }
774 
775 //____________________________________________________
776 
778  const EventInfo &eventInfo,
780  unsigned int iHit,
781  GblPoint &gblPoint) {
782  AlignmentParameters *params = nullptr;
783  std::vector<double> theDoubleBufferX, theDoubleBufferY;
784  theDoubleBufferX.clear();
785  theDoubleBufferY.clear();
786  theIntBuffer.clear();
787  int iret = 0;
788 
789  const TrajectoryStateOnSurface &tsos = refTrajPtr->trajectoryStates()[iHit];
790  const ConstRecHitPointer &recHitPtr = refTrajPtr->recHits()[iHit];
791  // ignore invalid hits
792  if (!recHitPtr->isValid())
793  return 0;
794 
795  // get AlignableDet/Unit for this hit
796  AlignableDetOrUnitPtr alidet(theAlignableNavigator->alignableFromDetId(recHitPtr->geographicalId()));
797 
798  if (!this->globalDerivativesHierarchy(eventInfo,
799  tsos,
800  alidet,
801  alidet,
802  theDoubleBufferX, // 2x alidet, sic!
803  theDoubleBufferY,
804  theIntBuffer,
805  params)) {
806  return -1; // problem
807  }
808  //calibration parameters
809  int globalLabel;
810  std::vector<IntegratedCalibrationBase::ValuesIndexPair> derivs;
811  for (auto iCalib = theCalibrations.begin(); iCalib != theCalibrations.end(); ++iCalib) {
812  // get all derivatives of this calibration // const unsigned int num =
813  (*iCalib)->derivatives(derivs, *recHitPtr, tsos, setup, eventInfo);
814  for (auto iValuesInd = derivs.begin(); iValuesInd != derivs.end(); ++iValuesInd) {
815  // transfer label and x/y derivatives
816  globalLabel = thePedeLabels->calibrationLabel(*iCalib, iValuesInd->second);
817  if (globalLabel > 0 && globalLabel <= 2147483647) {
818  theIntBuffer.push_back(globalLabel);
819  theDoubleBufferX.push_back(iValuesInd->first.first);
820  theDoubleBufferY.push_back(iValuesInd->first.second);
821  } else {
822  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::addGlobalData"
823  << "Invalid label " << globalLabel << " <= 0 or > 2147483647";
824  }
825  }
826  }
827  unsigned int numGlobals = theIntBuffer.size();
828  if (numGlobals > 0) {
829  Eigen::Matrix<double, 2, Eigen::Dynamic> globalDer{2, numGlobals};
830  for (unsigned int i = 0; i < numGlobals; ++i) {
831  globalDer(0, i) = theDoubleBufferX[i];
832  globalDer(1, i) = theDoubleBufferY[i];
833  }
834  gblPoint.addGlobals(theIntBuffer, globalDer);
835  iret = 1;
836  }
837  return iret;
838 }
839 
840 //____________________________________________________
842  const TrajectoryStateOnSurface &tsos,
843  Alignable *ali,
844  const AlignableDetOrUnitPtr &alidet,
845  std::vector<float> &globalDerivativesX,
846  std::vector<float> &globalDerivativesY,
847  std::vector<int> &globalLabels,
848  AlignmentParameters *&lowestParams) const {
849  // derivatives and labels are recursively attached
850  if (!ali)
851  return true; // no mother might be OK
852 
853  if (false && theMonitor && alidet != ali)
854  theMonitor->fillFrameToFrame(alidet, ali);
855 
857 
858  if (params) {
859  if (!lowestParams)
860  lowestParams = params; // set parameters of lowest level
861 
862  bool hasSplitParameters = thePedeLabels->hasSplitParameters(ali);
863  const unsigned int alignableLabel = thePedeLabels->alignableLabel(ali);
864 
865  if (0 == alignableLabel) { // FIXME: what about regardAllHits in Markus' code?
866  edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::globalDerivativesHierarchy"
867  << "Label not found, skip Alignable.";
868  return false;
869  }
870 
871  const std::vector<bool> &selPars = params->selector();
872  const AlgebraicMatrix derivs(params->derivatives(tsos, alidet));
873 
874  // cols: 2, i.e. x&y, rows: parameters, usually RigidBodyAlignmentParameters::N_PARAM
875  for (unsigned int iSel = 0; iSel < selPars.size(); ++iSel) {
876  if (selPars[iSel]) {
877  globalDerivativesX.push_back(derivs[iSel][kLocalX] / thePedeSteer->cmsToPedeFactor(iSel));
878  if (hasSplitParameters == true) {
879  globalLabels.push_back(thePedeLabels->parameterLabel(ali, iSel, eventInfo, tsos));
880  } else {
881  globalLabels.push_back(thePedeLabels->parameterLabel(alignableLabel, iSel));
882  }
883  globalDerivativesY.push_back(derivs[iSel][kLocalY] / thePedeSteer->cmsToPedeFactor(iSel));
884  }
885  }
886  // Exclude mothers if Alignable selected to be no part of a hierarchy:
887  if (thePedeSteer->isNoHiera(ali))
888  return true;
889  }
890  // Call recursively for mother, will stop if mother == 0:
891  return this->globalDerivativesHierarchy(
892  eventInfo, tsos, ali->mother(), alidet, globalDerivativesX, globalDerivativesY, globalLabels, lowestParams);
893 }
894 
895 //____________________________________________________
897  const TrajectoryStateOnSurface &tsos,
898  Alignable *ali,
899  const AlignableDetOrUnitPtr &alidet,
900  std::vector<double> &globalDerivativesX,
901  std::vector<double> &globalDerivativesY,
902  std::vector<int> &globalLabels,
903  AlignmentParameters *&lowestParams) const {
904  // derivatives and labels are recursively attached
905  if (!ali)
906  return true; // no mother might be OK
907 
908  if (false && theMonitor && alidet != ali)
909  theMonitor->fillFrameToFrame(alidet, ali);
910 
912 
913  if (params) {
914  if (!lowestParams)
915  lowestParams = params; // set parameters of lowest level
916 
917  bool hasSplitParameters = thePedeLabels->hasSplitParameters(ali);
918  const unsigned int alignableLabel = thePedeLabels->alignableLabel(ali);
919 
920  if (0 == alignableLabel) { // FIXME: what about regardAllHits in Markus' code?
921  edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::globalDerivativesHierarchy"
922  << "Label not found, skip Alignable.";
923  return false;
924  }
925 
926  const std::vector<bool> &selPars = params->selector();
927  const AlgebraicMatrix derivs(params->derivatives(tsos, alidet));
928  int globalLabel;
929 
930  // cols: 2, i.e. x&y, rows: parameters, usually RigidBodyAlignmentParameters::N_PARAM
931  for (unsigned int iSel = 0; iSel < selPars.size(); ++iSel) {
932  if (selPars[iSel]) {
933  if (hasSplitParameters == true) {
934  globalLabel = thePedeLabels->parameterLabel(ali, iSel, eventInfo, tsos);
935  } else {
936  globalLabel = thePedeLabels->parameterLabel(alignableLabel, iSel);
937  }
938  if (globalLabel > 0 && globalLabel <= 2147483647) {
939  globalLabels.push_back(globalLabel);
940  globalDerivativesX.push_back(derivs[iSel][kLocalX] / thePedeSteer->cmsToPedeFactor(iSel));
941  globalDerivativesY.push_back(derivs[iSel][kLocalY] / thePedeSteer->cmsToPedeFactor(iSel));
942  } else {
943  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::globalDerivativesHierarchy"
944  << "Invalid label " << globalLabel << " <= 0 or > 2147483647";
945  }
946  }
947  }
948  // Exclude mothers if Alignable selected to be no part of a hierarchy:
949  if (thePedeSteer->isNoHiera(ali))
950  return true;
951  }
952  // Call recursively for mother, will stop if mother == 0:
953  return this->globalDerivativesHierarchy(
954  eventInfo, tsos, ali->mother(), alidet, globalDerivativesX, globalDerivativesY, globalLabels, lowestParams);
955 }
956 
957 //____________________________________________________
959  const TrajectoryStateOnSurface &tsos,
960  const edm::EventSetup &setup,
961  const EventInfo &eventInfo,
962  std::vector<float> &globalDerivativesX,
963  std::vector<float> &globalDerivativesY,
964  std::vector<int> &globalLabels) const {
965  std::vector<IntegratedCalibrationBase::ValuesIndexPair> derivs;
966  for (auto iCalib = theCalibrations.begin(); iCalib != theCalibrations.end(); ++iCalib) {
967  // get all derivatives of this calibration // const unsigned int num =
968  (*iCalib)->derivatives(derivs, *recHit, tsos, setup, eventInfo);
969  for (auto iValuesInd = derivs.begin(); iValuesInd != derivs.end(); ++iValuesInd) {
970  // transfer label and x/y derivatives
971  globalLabels.push_back(thePedeLabels->calibrationLabel(*iCalib, iValuesInd->second));
972  globalDerivativesX.push_back(iValuesInd->first.first);
973  globalDerivativesY.push_back(iValuesInd->first.second);
974  }
975  }
976 }
977 
978 // //____________________________________________________
979 // void MillePedeAlignmentAlgorithm
980 // ::callMille(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
981 // unsigned int iTrajHit, MeasurementDirection xOrY,
982 // const std::vector<float> &globalDerivatives, const std::vector<int> &globalLabels)
983 // {
984 // const unsigned int xyIndex = iTrajHit*2 + xOrY;
985 // // FIXME: here for residuum and sigma we could use KALMAN-Filter results
986 // const float residuum =
987 // refTrajPtr->measurements()[xyIndex] - refTrajPtr->trajectoryPositions()[xyIndex];
988 // const float covariance = refTrajPtr->measurementErrors()[xyIndex][xyIndex];
989 // const float sigma = (covariance > 0. ? TMath::Sqrt(covariance) : 0.);
990 
991 // const AlgebraicMatrix &locDerivMatrix = refTrajPtr->derivatives();
992 
993 // std::vector<float> localDerivs(locDerivMatrix.num_col());
994 // for (unsigned int i = 0; i < localDerivs.size(); ++i) {
995 // localDerivs[i] = locDerivMatrix[xyIndex][i];
996 // }
997 
998 // // &(vector[0]) is valid - as long as vector is not empty
999 // // cf. http://www.parashift.com/c++-faq-lite/containers.html#faq-34.3
1000 // theMille->mille(localDerivs.size(), &(localDerivs[0]),
1001 // globalDerivatives.size(), &(globalDerivatives[0]), &(globalLabels[0]),
1002 // residuum, sigma);
1003 // if (theMonitor) {
1004 // theMonitor->fillDerivatives(refTrajPtr->recHits()[iTrajHit],localDerivs, globalDerivatives,
1005 // (xOrY == kLocalY));
1006 // theMonitor->fillResiduals(refTrajPtr->recHits()[iTrajHit],
1007 // refTrajPtr->trajectoryStates()[iTrajHit],
1008 // iTrajHit, residuum, sigma, (xOrY == kLocalY));
1009 // }
1010 // }
1011 
1012 //____________________________________________________
1014  // FIXME: Check whether this is a reliable and recommended way to find out...
1015 
1016  if (recHit->dimension() < 2) {
1017  return false; // some muon and TIB/TOB stuff really has RecHit1D
1018  } else if (recHit->detUnit()) { // detunit in strip is 1D, in pixel 2D
1019  return recHit->detUnit()->type().isTrackerPixel();
1020  } else { // stereo strips (FIXME: endcap trouble due to non-parallel strips (wedge sensors)?)
1021  if (dynamic_cast<const ProjectedSiStripRecHit2D *>(recHit->hit())) { // check persistent hit
1022  // projected: 1D measurement on 'glued' module
1023  return false;
1024  } else {
1025  return true;
1026  }
1027  }
1028 }
1029 
1030 //__________________________________________________________________________________________________
1032  bool setUserVars,
1033  const RunRange &runrange) {
1034  bool allEmpty = this->areEmptyParams(theAlignables);
1035 
1036  PedeReader reader(mprespset, *thePedeSteer, *thePedeLabels, runrange);
1037  align::Alignables alis;
1038  bool okRead = reader.read(alis, setUserVars); // also may set params of IntegratedCalibration's
1039  bool numMatch = true;
1040 
1041  std::stringstream out;
1042  out << "Read " << alis.size() << " alignables";
1043  if (alis.size() != theAlignables.size()) {
1044  out << " while " << theAlignables.size() << " in store";
1045  numMatch = false; // FIXME: Should we check one by one? Or transfer 'alis' to the store?
1046  }
1047  if (!okRead)
1048  out << ", but problems in reading";
1049  if (!allEmpty)
1050  out << ", possibly overwriting previous settings";
1051  out << ".";
1052 
1053  if (okRead && allEmpty) {
1054  if (numMatch) { // as many alignables with result as trying to align
1055  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
1056  } else if (!alis.empty()) { // dead module do not get hits and no pede result
1057  edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
1058  } else { // serious problem: no result read - and not all modules can be dead...
1059  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
1060  return false;
1061  }
1062  return true;
1063  }
1064  // the rest is not OK:
1065  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
1066  return false;
1067 }
1068 
1069 //__________________________________________________________________________________________________
1071  for (const auto &iAli : alignables) {
1072  const AlignmentParameters *params = iAli->alignmentParameters();
1073  if (params) {
1074  const auto &parVec(params->parameters());
1075  const auto &parCov(params->covariance());
1076  for (int i = 0; i < parVec.num_row(); ++i) {
1077  if (parVec[i] != 0.)
1078  return false;
1079  for (int j = i; j < parCov.num_col(); ++j) {
1080  if (parCov[i][j] != 0.)
1081  return false;
1082  }
1083  }
1084  }
1085  }
1086 
1087  return true;
1088 }
1089 
1090 //__________________________________________________________________________________________________
1091 unsigned int MillePedeAlignmentAlgorithm::doIO(int loop) const {
1092  unsigned int result = 0;
1093 
1094  const std::string outFilePlain(theConfig.getParameter<std::string>("treeFile"));
1095  if (outFilePlain.empty()) {
1096  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
1097  << "treeFile parameter empty => skip writing for 'loop' " << loop;
1098  return result;
1099  }
1100 
1101  const std::string outFile(theDir + outFilePlain);
1102 
1103  AlignmentIORoot aliIO;
1104  int ioerr = 0;
1105  if (loop == 0) {
1106  aliIO.writeAlignableOriginalPositions(theAlignables, outFile.c_str(), loop, false, ioerr);
1107  if (ioerr) {
1108  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
1109  << "Problem " << ioerr << " in writeAlignableOriginalPositions";
1110  ++result;
1111  }
1112  } else if (loop == 1) {
1113  // only for first iov add hit counts, else 2x, 3x,... number of hits in IOV 2, 3,...
1114  const std::vector<std::string> inFiles(theConfig.getParameter<std::vector<std::string> >("mergeTreeFiles"));
1115  const std::vector<std::string> binFiles(theConfig.getParameter<std::vector<std::string> >("mergeBinaryFiles"));
1116  if (inFiles.size() != binFiles.size()) {
1117  edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
1118  << "'vstring mergeTreeFiles' and 'vstring mergeBinaryFiles' "
1119  << "differ in size";
1120  }
1121  this->addHitStatistics(0, outFile, inFiles); // add hit info from tree 0 in 'infiles'
1122  }
1123  MillePedeVariablesIORoot millePedeIO;
1124  millePedeIO.writeMillePedeVariables(theAlignables, outFile.c_str(), loop, false, ioerr);
1125  if (ioerr) {
1126  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
1127  << "Problem " << ioerr << " writing MillePedeVariables";
1128  ++result;
1129  }
1130 
1131  aliIO.writeOrigRigidBodyAlignmentParameters(theAlignables, outFile.c_str(), loop, false, ioerr);
1132  if (ioerr) {
1133  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
1134  << "Problem " << ioerr << " in writeOrigRigidBodyAlignmentParameters, " << loop;
1135  ++result;
1136  }
1137  aliIO.writeAlignableAbsolutePositions(theAlignables, outFile.c_str(), loop, false, ioerr);
1138  if (ioerr) {
1139  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
1140  << "Problem " << ioerr << " in writeAlignableAbsolutePositions, " << loop;
1141  ++result;
1142  }
1143 
1144  return result;
1145 }
1146 
1147 //__________________________________________________________________________________________________
1149  for (const auto &iAli : alis) {
1150  AlignmentParameters *params = iAli->alignmentParameters();
1151  if (!params) {
1152  throw cms::Exception("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::buildUserVariables"
1153  << "No parameters for alignable";
1154  }
1155  MillePedeVariables *userVars = dynamic_cast<MillePedeVariables *>(params->userVariables());
1156  if (userVars) { // Just re-use existing, keeping label and numHits:
1157  for (unsigned int iPar = 0; iPar < userVars->size(); ++iPar) {
1158  // if (params->hierarchyLevel() > 0) {
1159  //std::cout << params->hierarchyLevel() << "\nBefore: " << userVars->parameter()[iPar];
1160  // }
1161  userVars->setAllDefault(iPar);
1162  //std::cout << "\nAfter: " << userVars->parameter()[iPar] << std::endl;
1163  }
1164  } else { // Nothing yet or erase wrong type:
1165  userVars = new MillePedeVariables(
1166  params->size(),
1167  thePedeLabels->alignableLabel(iAli),
1168  thePedeLabels->alignableTracker()->objectIdProvider().typeToName(iAli->alignableObjectId()));
1169  params->setUserVariables(userVars);
1170  }
1171  }
1172 }
1173 
1174 //__________________________________________________________________________________________________
1176  if (mode == "full") {
1178  } else if (mode == "mille") {
1179  return myMilleBit; // + myPedeSteerBit; // sic! Including production of steerig file. NO!
1180  } else if (mode == "pede") {
1182  } else if (mode == "pedeSteer") {
1183  return myPedeSteerBit;
1184  } else if (mode == "pedeRun") {
1185  return myPedeSteerBit + myPedeRunBit + myPedeReadBit; // sic! Including steering and reading of result.
1186  } else if (mode == "pedeRead") {
1187  return myPedeReadBit;
1188  }
1189 
1190  throw cms::Exception("BadConfig") << "Unknown mode '" << mode
1191  << "', use 'full', 'mille', 'pede', 'pedeRun', 'pedeSteer' or 'pedeRead'.";
1192 
1193  return 0;
1194 }
1195 
1196 //__________________________________________________________________________________________________
1198  const std::string &outFile,
1199  const std::vector<std::string> &inFiles) const {
1200  bool allOk = true;
1201  int ierr = 0;
1202  MillePedeVariablesIORoot millePedeIO;
1203  for (std::vector<std::string>::const_iterator iFile = inFiles.begin(); iFile != inFiles.end(); ++iFile) {
1204  const std::string inFile(theDir + *iFile);
1205  const std::vector<AlignmentUserVariables *> mpVars =
1206  millePedeIO.readMillePedeVariables(theAlignables, inFile.c_str(), fromIov, ierr);
1207  if (ierr || !this->addHits(theAlignables, mpVars)) {
1208  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::addHitStatistics"
1209  << "Error " << ierr << " reading from " << inFile << ", tree " << fromIov
1210  << ", or problems in addHits";
1211  allOk = false;
1212  }
1213  for (std::vector<AlignmentUserVariables *>::const_iterator i = mpVars.begin(); i != mpVars.end(); ++i) {
1214  delete *i; // clean created objects
1215  }
1216  }
1217 
1218  return allOk;
1219 }
1220 
1221 //__________________________________________________________________________________________________
1223  const std::vector<AlignmentUserVariables *> &mpVars) const {
1224  bool allOk = (mpVars.size() == alis.size());
1225  std::vector<AlignmentUserVariables *>::const_iterator iUser = mpVars.begin();
1226  for (auto iAli = alis.cbegin(); iAli != alis.cend() && iUser != mpVars.end(); ++iAli, ++iUser) {
1227  MillePedeVariables *mpVarNew = dynamic_cast<MillePedeVariables *>(*iUser);
1228  AlignmentParameters *ps = (*iAli)->alignmentParameters();
1229  MillePedeVariables *mpVarOld = (ps ? dynamic_cast<MillePedeVariables *>(ps->userVariables()) : nullptr);
1230  if (!mpVarNew || !mpVarOld || mpVarOld->size() != mpVarNew->size()) {
1231  allOk = false;
1232  continue; // FIXME error etc.?
1233  }
1234 
1235  mpVarOld->increaseHitsX(mpVarNew->hitsX());
1236  mpVarOld->increaseHitsY(mpVarNew->hitsY());
1237  }
1238 
1239  return allOk;
1240 }
1241 
1242 //__________________________________________________________________________________________________
1243 template <typename GlobalDerivativeMatrix>
1244 void MillePedeAlignmentAlgorithm::makeGlobDerivMatrix(const std::vector<float> &globalDerivativesx,
1245  const std::vector<float> &globalDerivativesy,
1246  Eigen::MatrixBase<GlobalDerivativeMatrix> &aGlobalDerivativesM) {
1247  static_assert(GlobalDerivativeMatrix::RowsAtCompileTime == 2, "global derivative matrix must have two rows");
1248 
1249  for (size_t i = 0; i < globalDerivativesx.size(); ++i) {
1250  aGlobalDerivativesM(0, i) = globalDerivativesx[i];
1251  aGlobalDerivativesM(1, i) = globalDerivativesy[i];
1252  }
1253 }
1254 
1255 //__________________________________________________________________________________________________
1256 template <typename CovarianceMatrix,
1257  typename LocalDerivativeMatrix,
1258  typename ResidualMatrix,
1259  typename GlobalDerivativeMatrix>
1260 void MillePedeAlignmentAlgorithm::diagonalize(Eigen::MatrixBase<CovarianceMatrix> &aHitCovarianceM,
1261  Eigen::MatrixBase<LocalDerivativeMatrix> &aLocalDerivativesM,
1262  Eigen::MatrixBase<ResidualMatrix> &aHitResidualsM,
1263  Eigen::MatrixBase<GlobalDerivativeMatrix> &aGlobalDerivativesM) const {
1265  "'aLocalDerivativesM' and 'aHitResidualsM' must have the "
1266  "same underlying scalar type");
1268  "'aLocalDerivativesM' and 'aGlobalDerivativesM' must have the "
1269  "same underlying scalar type");
1270 
1271  Eigen::SelfAdjointEigenSolver<typename CovarianceMatrix::PlainObject> myDiag{aHitCovarianceM};
1272  // eigenvectors of real symmetric matrices are orthogonal, i.e. invert == transpose
1273  auto aTranfoToDiagonalSystemInv =
1274  myDiag.eigenvectors().transpose().template cast<typename LocalDerivativeMatrix::Scalar>();
1275 
1276  aHitCovarianceM = myDiag.eigenvalues().asDiagonal();
1277  aLocalDerivativesM = aTranfoToDiagonalSystemInv * aLocalDerivativesM;
1278  aHitResidualsM = aTranfoToDiagonalSystemInv * aHitResidualsM;
1279  if (aGlobalDerivativesM.size() > 0) {
1280  // diagonalize only if measurement depends on alignables or calibrations
1281  aGlobalDerivativesM = aTranfoToDiagonalSystemInv * aGlobalDerivativesM;
1282  }
1283 }
1284 
1285 //__________________________________________________________________________________________________
1286 template <typename CovarianceMatrix, typename ResidualMatrix, typename LocalDerivativeMatrix>
1289  unsigned int iVirtualMeas,
1290  Eigen::MatrixBase<CovarianceMatrix> &aHitCovarianceM,
1291  Eigen::MatrixBase<ResidualMatrix> &aHitResidualsM,
1292  Eigen::MatrixBase<LocalDerivativeMatrix> &aLocalDerivativesM) {
1293  // This Method is valid for 1D measurements only
1294 
1295  const unsigned int xIndex = iVirtualMeas + refTrajPtr->numberOfHitMeas();
1296 
1297  aHitCovarianceM(0, 0) = refTrajPtr->measurementErrors()[xIndex][xIndex];
1298  aHitResidualsM(0, 0) = refTrajPtr->measurements()[xIndex];
1299 
1300  const auto &locDerivMatrix = refTrajPtr->derivatives();
1301  for (int i = 0; i < locDerivMatrix.num_col(); ++i) {
1302  aLocalDerivativesM(0, i) = locDerivMatrix[xIndex][i];
1303  }
1304 }
1305 
1306 //__________________________________________________________________________________________________
1307 template <typename CovarianceMatrix, typename ResidualMatrix, typename LocalDerivativeMatrix>
1309  unsigned int iTrajHit,
1310  Eigen::MatrixBase<CovarianceMatrix> &aHitCovarianceM,
1311  Eigen::MatrixBase<ResidualMatrix> &aHitResidualsM,
1312  Eigen::MatrixBase<LocalDerivativeMatrix> &aLocalDerivativesM) {
1313  // This Method is valid for 2D measurements only
1314 
1315  const unsigned int xIndex = iTrajHit * 2;
1316  const unsigned int yIndex = iTrajHit * 2 + 1;
1317 
1318  aHitCovarianceM(0, 0) = refTrajPtr->measurementErrors()[xIndex][xIndex];
1319  aHitCovarianceM(0, 1) = refTrajPtr->measurementErrors()[xIndex][yIndex];
1320  aHitCovarianceM(1, 0) = refTrajPtr->measurementErrors()[yIndex][xIndex];
1321  aHitCovarianceM(1, 1) = refTrajPtr->measurementErrors()[yIndex][yIndex];
1322 
1323  aHitResidualsM(0, 0) = refTrajPtr->measurements()[xIndex] - refTrajPtr->trajectoryPositions()[xIndex];
1324  aHitResidualsM(1, 0) = refTrajPtr->measurements()[yIndex] - refTrajPtr->trajectoryPositions()[yIndex];
1325 
1326  const auto &locDerivMatrix = refTrajPtr->derivatives();
1327  for (int i = 0; i < locDerivMatrix.num_col(); ++i) {
1328  aLocalDerivativesM(0, i) = locDerivMatrix[xIndex][i];
1329  aLocalDerivativesM(1, i) = locDerivMatrix[yIndex][i];
1330  }
1331 }
1332 
1333 //__________________________________________________________________________________________________
1335  unsigned int iTrajHit,
1336  const std::vector<int> &globalLabels,
1337  const std::vector<float> &globalDerivativesX,
1338  const std::vector<float> &globalDerivativesY) {
1339  const ConstRecHitPointer aRecHit(refTrajPtr->recHits()[iTrajHit]);
1340 
1341  if ((aRecHit)->dimension() == 1) {
1342  return this->callMille1D(refTrajPtr, iTrajHit, globalLabels, globalDerivativesX);
1343  } else {
1344  return this->callMille2D(refTrajPtr, iTrajHit, globalLabels, globalDerivativesX, globalDerivativesY);
1345  }
1346 }
1347 
1348 //__________________________________________________________________________________________________
1350  unsigned int iTrajHit,
1351  const std::vector<int> &globalLabels,
1352  const std::vector<float> &globalDerivativesX) {
1353  const ConstRecHitPointer aRecHit(refTrajPtr->recHits()[iTrajHit]);
1354  const unsigned int xIndex = iTrajHit * 2; // the even ones are local x
1355 
1356  // local derivatives
1357  const AlgebraicMatrix &locDerivMatrix = refTrajPtr->derivatives();
1358  const int nLocal = locDerivMatrix.num_col();
1359  std::vector<float> localDerivatives(nLocal);
1360  for (unsigned int i = 0; i < localDerivatives.size(); ++i) {
1361  localDerivatives[i] = locDerivMatrix[xIndex][i];
1362  }
1363 
1364  // residuum and error
1365  float residX = refTrajPtr->measurements()[xIndex] - refTrajPtr->trajectoryPositions()[xIndex];
1366  float hitErrX = TMath::Sqrt(refTrajPtr->measurementErrors()[xIndex][xIndex]);
1367 
1368  // number of global derivatives
1369  const int nGlobal = globalDerivativesX.size();
1370 
1371  // &(localDerivatives[0]) etc. are valid - as long as vector is not empty
1372  // cf. http://www.parashift.com/c++-faq-lite/containers.html#faq-34.3
1373  theMille->mille(
1374  nLocal, &(localDerivatives[0]), nGlobal, &(globalDerivativesX[0]), &(globalLabels[0]), residX, hitErrX);
1375 
1376  if (theMonitor) {
1377  theMonitor->fillDerivatives(
1378  aRecHit, &(localDerivatives[0]), nLocal, &(globalDerivativesX[0]), nGlobal, &(globalLabels[0]));
1379  theMonitor->fillResiduals(aRecHit, refTrajPtr->trajectoryStates()[iTrajHit], iTrajHit, residX, hitErrX, false);
1380  }
1381 
1382  return 1;
1383 }
1384 
1385 //__________________________________________________________________________________________________
1387  unsigned int iTrajHit,
1388  const std::vector<int> &globalLabels,
1389  const std::vector<float> &globalDerivativesx,
1390  const std::vector<float> &globalDerivativesy) {
1391  const ConstRecHitPointer aRecHit(refTrajPtr->recHits()[iTrajHit]);
1392 
1393  if ((aRecHit)->dimension() != 2) {
1394  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::callMille2D"
1395  << "You try to call method for 2D hits for a " << (aRecHit)->dimension()
1396  << "D Hit. Hit gets ignored!";
1397  return -1;
1398  }
1399 
1400  Eigen::Matrix<double, 2, 2> aHitCovarianceM;
1401  Eigen::Matrix<float, 2, 1> aHitResidualsM;
1402  Eigen::Matrix<float, 2, Eigen::Dynamic> aLocalDerivativesM{2, refTrajPtr->derivatives().num_col()};
1403  // below method fills above 3 matrices
1404  this->addRefTrackData2D(refTrajPtr, iTrajHit, aHitCovarianceM, aHitResidualsM, aLocalDerivativesM);
1405  Eigen::Matrix<float, 2, Eigen::Dynamic> aGlobalDerivativesM{2, globalDerivativesx.size()};
1406  this->makeGlobDerivMatrix(globalDerivativesx, globalDerivativesy, aGlobalDerivativesM);
1407 
1408  // calculates correlation between Hit measurements
1409  // FIXME: Should take correlation (and resulting transformation) from original hit,
1410  // not 2x2 matrix from ReferenceTrajectory: That can come from error propagation etc.!
1411  const double corr = aHitCovarianceM(0, 1) / sqrt(aHitCovarianceM(0, 0) * aHitCovarianceM(1, 1));
1412  if (theMonitor)
1413  theMonitor->fillCorrelations2D(corr, aRecHit);
1414  bool diag = false; // diagonalise only tracker TID, TEC
1415  switch (aRecHit->geographicalId().subdetId()) {
1416  case SiStripDetId::TID:
1417  case SiStripDetId::TEC:
1418  if (aRecHit->geographicalId().det() == DetId::Tracker && TMath::Abs(corr) > theMaximalCor2D) {
1419  this->diagonalize(aHitCovarianceM, aLocalDerivativesM, aHitResidualsM, aGlobalDerivativesM);
1420  diag = true;
1421  }
1422  break;
1423  default:;
1424  }
1425 
1426  float newResidX = aHitResidualsM(0, 0);
1427  float newResidY = aHitResidualsM(1, 0);
1428  float newHitErrX = TMath::Sqrt(aHitCovarianceM(0, 0));
1429  float newHitErrY = TMath::Sqrt(aHitCovarianceM(1, 1));
1430 
1431  // change from column major (Eigen default) to row major to have row entries
1432  // in continuous memory
1433  std::vector<float> newLocalDerivs(aLocalDerivativesM.size());
1434  Eigen::Map<Eigen::Matrix<float, 2, Eigen::Dynamic, Eigen::RowMajor> >(
1435  newLocalDerivs.data(), aLocalDerivativesM.rows(), aLocalDerivativesM.cols()) = aLocalDerivativesM;
1436  float *newLocalDerivsX = &(newLocalDerivs[0]);
1437  float *newLocalDerivsY = &(newLocalDerivs[aLocalDerivativesM.cols()]);
1438 
1439  // change from column major (Eigen default) to row major to have row entries
1440  // in continuous memory
1441  std::vector<float> newGlobDerivs(aGlobalDerivativesM.size());
1442  Eigen::Map<Eigen::Matrix<float, 2, Eigen::Dynamic, Eigen::RowMajor> >(
1443  newGlobDerivs.data(), aGlobalDerivativesM.rows(), aGlobalDerivativesM.cols()) = aGlobalDerivativesM;
1444  float *newGlobDerivsX = &(newGlobDerivs[0]);
1445  float *newGlobDerivsY = &(newGlobDerivs[aGlobalDerivativesM.cols()]);
1446 
1447  const int nLocal = aLocalDerivativesM.cols();
1448  const int nGlobal = aGlobalDerivativesM.cols();
1449 
1450  if (diag && (newHitErrX > newHitErrY)) { // also for 2D hits?
1451  // measurement with smaller error is x-measurement (for !is2D do not fill y-measurement):
1452  std::swap(newResidX, newResidY);
1453  std::swap(newHitErrX, newHitErrY);
1454  std::swap(newLocalDerivsX, newLocalDerivsY);
1455  std::swap(newGlobDerivsX, newGlobDerivsY);
1456  }
1457 
1458  // &(globalLabels[0]) is valid - as long as vector is not empty
1459  // cf. http://www.parashift.com/c++-faq-lite/containers.html#faq-34.3
1460  theMille->mille(nLocal, newLocalDerivsX, nGlobal, newGlobDerivsX, &(globalLabels[0]), newResidX, newHitErrX);
1461 
1462  if (theMonitor) {
1463  theMonitor->fillDerivatives(aRecHit, newLocalDerivsX, nLocal, newGlobDerivsX, nGlobal, &(globalLabels[0]));
1464  theMonitor->fillResiduals(
1465  aRecHit, refTrajPtr->trajectoryStates()[iTrajHit], iTrajHit, newResidX, newHitErrX, false);
1466  }
1467  const bool isReal2DHit = this->is2D(aRecHit); // strip is 1D (except matched hits)
1468  if (isReal2DHit) {
1469  theMille->mille(nLocal, newLocalDerivsY, nGlobal, newGlobDerivsY, &(globalLabels[0]), newResidY, newHitErrY);
1470  if (theMonitor) {
1471  theMonitor->fillDerivatives(aRecHit, newLocalDerivsY, nLocal, newGlobDerivsY, nGlobal, &(globalLabels[0]));
1472  theMonitor->fillResiduals(
1473  aRecHit, refTrajPtr->trajectoryStates()[iTrajHit], iTrajHit, newResidY, newHitErrY, true); // true: y
1474  }
1475  }
1476 
1477  return (isReal2DHit ? 2 : 1);
1478 }
1479 
1480 //__________________________________________________________________________________________________
1482  unsigned int iVirtualMeas) {
1483  Eigen::Matrix<double, 1, 1> aHitCovarianceM;
1484  Eigen::Matrix<float, 1, 1> aHitResidualsM;
1485  Eigen::Matrix<float, 1, Eigen::Dynamic> aLocalDerivativesM{1, refTrajPtr->derivatives().num_col()};
1486  // below method fills above 3 'matrices'
1487  this->addRefTrackVirtualMeas1D(refTrajPtr, iVirtualMeas, aHitCovarianceM, aHitResidualsM, aLocalDerivativesM);
1488 
1489  // no global parameters (use dummy 0)
1490  auto aGlobalDerivativesM = Eigen::Matrix<float, 1, 1>::Zero();
1491 
1492  float newResidX = aHitResidualsM(0, 0);
1493  float newHitErrX = TMath::Sqrt(aHitCovarianceM(0, 0));
1494  std::vector<float> newLocalDerivsX(aLocalDerivativesM.size());
1495  Eigen::Map<Eigen::Matrix<float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> >(
1496  newLocalDerivsX.data(), aLocalDerivativesM.rows(), aLocalDerivativesM.cols()) = aLocalDerivativesM;
1497 
1498  std::vector<float> newGlobDerivsX(aGlobalDerivativesM.size());
1499  Eigen::Map<Eigen::Matrix<float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> >(
1500  newGlobDerivsX.data(), aGlobalDerivativesM.rows(), aGlobalDerivativesM.cols()) = aGlobalDerivativesM;
1501 
1502  const int nLocal = aLocalDerivativesM.cols();
1503  const int nGlobal = 0;
1504 
1505  theMille->mille(nLocal, newLocalDerivsX.data(), nGlobal, newGlobDerivsX.data(), &nGlobal, newResidX, newHitErrX);
1506 }
1507 
1508 //____________________________________________________
1510  const TkFittedLasBeamCollection &lasBeams,
1511  const TsosVectorCollection &lasBeamTsoses) {
1512  TsosVectorCollection::const_iterator iTsoses = lasBeamTsoses.begin();
1513  for (TkFittedLasBeamCollection::const_iterator iBeam = lasBeams.begin(), iEnd = lasBeams.end(); iBeam != iEnd;
1514  ++iBeam, ++iTsoses) { // beam/tsoses parallel!
1515 
1516  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::addLaserData"
1517  << "Beam " << iBeam->getBeamId() << " with " << iBeam->parameters().size()
1518  << " parameters and " << iBeam->getData().size() << " hits.\n There are "
1519  << iTsoses->size() << " TSOSes.";
1520 
1521  this->addLasBeam(eventInfo, *iBeam, *iTsoses);
1522  }
1523 }
1524 
1525 //____________________________________________________
1527  const TkFittedLasBeam &lasBeam,
1528  const std::vector<TrajectoryStateOnSurface> &tsoses) {
1529  AlignmentParameters *dummyPtr = nullptr; // for globalDerivativesHierarchy()
1530  std::vector<float> lasLocalDerivsX; // buffer for local derivatives
1531  const unsigned int beamLabel = thePedeLabels->lasBeamLabel(lasBeam.getBeamId()); // for global par
1532 
1533  for (unsigned int iHit = 0; iHit < tsoses.size(); ++iHit) {
1534  if (!tsoses[iHit].isValid())
1535  continue;
1536  // clear buffer
1537  theFloatBufferX.clear();
1538  theFloatBufferY.clear();
1539  theIntBuffer.clear();
1540  lasLocalDerivsX.clear();
1541  // get alignables and global parameters
1542  const SiStripLaserRecHit2D &hit = lasBeam.getData()[iHit];
1543  AlignableDetOrUnitPtr lasAli(theAlignableNavigator->alignableFromDetId(hit.getDetId()));
1545  eventInfo, tsoses[iHit], lasAli, lasAli, theFloatBufferX, theFloatBufferY, theIntBuffer, dummyPtr);
1546  // fill derivatives vector from derivatives matrix
1547  for (unsigned int nFitParams = 0; nFitParams < static_cast<unsigned int>(lasBeam.parameters().size());
1548  ++nFitParams) {
1549  const float derivative = lasBeam.derivatives()[iHit][nFitParams];
1550  if (nFitParams < lasBeam.firstFixedParameter()) { // first local beam parameters
1551  lasLocalDerivsX.push_back(derivative);
1552  } else { // now global ones
1553  const unsigned int numPar = nFitParams - lasBeam.firstFixedParameter();
1554  theIntBuffer.push_back(thePedeLabels->parameterLabel(beamLabel, numPar));
1555  theFloatBufferX.push_back(derivative);
1556  }
1557  } // end loop over parameters
1558 
1559  const float residual = hit.localPosition().x() - tsoses[iHit].localPosition().x();
1560  // error from file or assume 0.003
1561  const float error = 0.003; // hit.localPositionError().xx(); sqrt???
1562 
1563  theMille->mille(lasLocalDerivsX.size(),
1564  &(lasLocalDerivsX[0]),
1565  theFloatBufferX.size(),
1566  &(theFloatBufferX[0]),
1567  &(theIntBuffer[0]),
1568  residual,
1569  error);
1570  } // end of loop over hits
1571 
1572  theMille->end();
1573 }
1574 
1576  // do some printing, if requested
1577  const bool doOutputOnStdout(pxbSurveyCfg.getParameter<bool>("doOutputOnStdout"));
1578  if (doOutputOnStdout) {
1579  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::addPxbSurvey"
1580  << "# Output from addPxbSurvey follows below because "
1581  << "doOutputOnStdout is set to True";
1582  }
1583 
1584  // instantiate a dicer object
1585  SurveyPxbDicer dicer(pxbSurveyCfg.getParameter<std::vector<edm::ParameterSet> >("toySurveyParameters"),
1586  pxbSurveyCfg.getParameter<unsigned int>("toySurveySeed"));
1587  std::ofstream outfile(pxbSurveyCfg.getUntrackedParameter<std::string>("toySurveyFile").c_str());
1588 
1589  // read data from file
1590  std::vector<SurveyPxbImageLocalFit> measurements;
1591  std::string filename(pxbSurveyCfg.getParameter<edm::FileInPath>("infile").fullPath());
1593 
1594  // loop over photographs (=measurements) and perform the fit
1595  for (std::vector<SurveyPxbImageLocalFit>::size_type i = 0; i != measurements.size(); i++) {
1596  if (doOutputOnStdout) {
1597  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::addPxbSurvey"
1598  << "Module " << i << ": ";
1599  }
1600 
1601  // get the Alignables and their surfaces
1602  AlignableDetOrUnitPtr mod1(theAlignableNavigator->alignableFromDetId(measurements[i].getIdFirst()));
1603  AlignableDetOrUnitPtr mod2(theAlignableNavigator->alignableFromDetId(measurements[i].getIdSecond()));
1604  const AlignableSurface &surf1 = mod1->surface();
1605  const AlignableSurface &surf2 = mod2->surface();
1606 
1607  // the position of the fiducial points in local frame of a PXB module
1608  const LocalPoint fidpoint0(-0.91, +3.30);
1609  const LocalPoint fidpoint1(+0.91, +3.30);
1610  const LocalPoint fidpoint2(+0.91, -3.30);
1611  const LocalPoint fidpoint3(-0.91, -3.30);
1612 
1613  // We choose the local frame of the first module as reference,
1614  // so take the fidpoints of the second module and calculate their
1615  // positions in the reference frame
1616  const GlobalPoint surf2point0(surf2.toGlobal(fidpoint0));
1617  const GlobalPoint surf2point1(surf2.toGlobal(fidpoint1));
1618  const LocalPoint fidpoint0inSurf1frame(surf1.toLocal(surf2point0));
1619  const LocalPoint fidpoint1inSurf1frame(surf1.toLocal(surf2point1));
1620 
1621  // Create the vector for the fit
1623  fidpointvec.push_back(fidpoint0inSurf1frame);
1624  fidpointvec.push_back(fidpoint1inSurf1frame);
1625  fidpointvec.push_back(fidpoint2);
1626  fidpointvec.push_back(fidpoint3);
1627 
1628  // if toy survey is requested, dice the values now
1629  if (pxbSurveyCfg.getParameter<bool>("doToySurvey")) {
1630  dicer.doDice(fidpointvec, measurements[i].getIdPair(), outfile);
1631  }
1632 
1633  // do the fit
1634  measurements[i].doFit(fidpointvec, thePedeLabels->alignableLabel(mod1), thePedeLabels->alignableLabel(mod2));
1635  SurveyPxbImageLocalFit::localpars_t a; // local pars from fit
1636  a = measurements[i].getLocalParameters();
1637  const SurveyPxbImageLocalFit::value_t chi2 = measurements[i].getChi2();
1638 
1639  // do some reporting, if requested
1640  if (doOutputOnStdout) {
1641  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::addPxbSurvey"
1642  << "a: " << a[0] << ", " << a[1] << ", " << a[2] << ", " << a[3]
1643  << " S= " << sqrt(a[2] * a[2] + a[3] * a[3]) << " phi= " << atan(a[3] / a[2])
1644  << " chi2= " << chi2 << std::endl;
1645  }
1646  if (theMonitor) {
1647  theMonitor->fillPxbSurveyHistsChi2(chi2);
1648  theMonitor->fillPxbSurveyHistsLocalPars(a[0], a[1], sqrt(a[2] * a[2] + a[3] * a[3]), atan(a[3] / a[2]));
1649  }
1650 
1651  // pass the results from the local fit to mille
1653  theMille->mille((int)measurements[i].getLocalDerivsSize(),
1654  measurements[i].getLocalDerivsPtr(j),
1655  (int)measurements[i].getGlobalDerivsSize(),
1656  measurements[i].getGlobalDerivsPtr(j),
1657  measurements[i].getGlobalDerivsLabelPtr(j),
1658  measurements[i].getResiduum(j),
1659  measurements[i].getSigma(j));
1660  }
1661  theMille->end();
1662  }
1663  outfile.close();
1664 }
1665 
1667  const auto runRangeSelection = theConfig.getUntrackedParameter<edm::VParameterSet>("RunRangeSelection");
1668 
1669  if (runRangeSelection.empty())
1670  return false;
1671 
1672  const auto runRanges =
1674 
1675  return !(runRanges.empty());
1676 }
RunNumber_t run() const
Definition: EventID.h:38
bool processesEvents() override
Returns whether MP should process events in the current configuration.
unsigned int hitsX() const
get number of hits for x-measurement
const TimeTypeSpecs timeTypeSpecs[]
Definition: Time.cc:16
unsigned int firstFixedParameter() const
void globalDerivativesCalibration(const TransientTrackingRecHit::ConstRecHitPointer &recHit, const TrajectoryStateOnSurface &tsos, const edm::EventSetup &setup, const EventInfo &eventInfo, std::vector< float > &globalDerivativesX, std::vector< float > &globalDerivativesY, std::vector< int > &globalLabels) const
adding derivatives from integrated calibrations
T getUntrackedParameter(std::string const &, T const &) const
unsigned int size() const
number of parameters
LocalPoint localPosition() const override
bool supportsCalibrations() override
Returns whether MP supports calibrations.
std::unique_ptr< MillePedeMonitor > theMonitor
tuple cfg
Definition: looper.py:296
void increaseHitsX(unsigned int add=1)
increase hits for x-measurement
bool setParametersForRunRange(const RunRange &runrange) override
RunNumber_t run() const
Definition: RunBase.h:40
Derivative< X, A >::type derivative(const A &_)
Definition: Derivative.h:18
virtual void endRun(const EventInfo &, const EndRunInfo &, const edm::EventSetup &)
Run on run products, e.g. TkLAS.
const bool isValid(const Frame &aFrame, const FrameQuality &aQuality, const uint16_t aExpectedPos)
std::vector< coord_t > fidpoint_t
std::shared_ptr< AlignPCLThresholds > theThresholds
void diagonalize(Eigen::MatrixBase< CovarianceMatrix > &aHitCovarianceM, Eigen::MatrixBase< LocalDerivativeMatrix > &aLocalDerivativesM, Eigen::MatrixBase< ResidualMatrix > &aHitResidualsM, Eigen::MatrixBase< GlobalDerivativeMatrix > &aGlobalDerivativesM) const
std::map< std::string, AlignPCLThreshold > threshold_map
const EventID & eventID() const
Definition: IOVSyncValue.h:40
bool areEmptyParams(const align::Alignables &alignables) const
void increaseHitsY(unsigned int add=1)
increase hits for y-measurement
Time_t beginValue
Definition: Time.h:41
void endLuminosityBlock(const edm::EventSetup &) override
called at end of luminosity block
static constexpr auto TID
Definition: SiStripDetId.h:38
void buildUserVariables(const align::Alignables &alignables) const
add MillePedeVariables for each AlignmentParameters (exception if no parameters...)
TrajectoryFactoryBase::ReferenceTrajectoryCollection RefTrajColl
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:34
std::unique_ptr< TrajectoryFactoryBase > theTrajectoryFactory
std::unique_ptr< gbl::MilleBinary > theBinary
const Time_t MIN_VAL(0)
int addGlobalData(const edm::EventSetup &setup, const EventInfo &eventInfo, const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iHit, gbl::GblPoint &gblPoint)
bool exists(std::string const &parameterName) const
checks if a parameter exists
bool globalDerivativesHierarchy(const EventInfo &eventInfo, const TrajectoryStateOnSurface &tsos, Alignable *ali, const AlignableDetOrUnitPtr &alidet, std::vector< float > &globalDerivativesX, std::vector< float > &globalDerivativesY, std::vector< int > &globalLabels, AlignmentParameters *&lowestParams) const
recursively adding derivatives and labels, false if problems
std::shared_ptr< PedeLabelerBase > thePedeLabels
string firstRun
Definition: dataset.py:936
std::vector< std::string > getExistingFormattedFiles(const std::vector< std::string > &plainFiles, const std::string &theDir)
auto const & tracks
cannot be loose
unsigned int addHitCount(const std::vector< AlignmentParameters * > &parVec, const std::vector< bool > &validHitVecY) const
bool addCalibrations(const std::vector< IntegratedCalibrationBase * > &iCals) override
Pass integrated calibrations to Millepede (they are not owned by Millepede!)
const AlgebraicMatrix & derivatives() const
matrix of local derivatives: columns are parameters, rows are hits
const std::vector< bool > & selector(void) const
Get alignment parameter selector vector.
int callMille2D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iTrajHit, const std::vector< int > &globalLabels, const std::vector< float > &globalDerivativesx, const std::vector< float > &globalDerivativesy)
math::Error< 5 >::type CovarianceMatrix
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
void writeMillePedeVariables(const align::Alignables &alivec, const char *filename, int iter, bool validCheck, int &ierr)
void writeAlignableOriginalPositions(const align::Alignables &alivec, const char *filename, int iter, bool validCheck, int &ierr) override
write Alignable original (before misalignment) absolute positions
uint16_t size_type
void addLasBeam(const EventInfo &eventInfo, const TkFittedLasBeam &lasBeam, const std::vector< TrajectoryStateOnSurface > &tsoses)
AlignmentParameters * alignmentParameters() const
Get the AlignmentParameters.
Definition: Alignable.h:58
const ConstTrajTrackPairCollection & trajTrackPairs() const
define event information passed to algorithms
tuple result
Definition: mps_fire.py:311
unsigned int doIO(int loop) const
int callMille(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iTrajHit, const std::vector< int > &globalLabels, const std::vector< float > &globalDerivativesX, const std::vector< float > &globalDerivativesY)
calls callMille1D or callMille2D
TransientTrackingRecHit::ConstRecHitPointer ConstRecHitPointer
bool getData(T &iHolder) const
Definition: EventSetup.h:122
U second(std::pair< T, U > const &p)
const AlgebraicVector & parameters(void) const
Get alignment parameters.
void beginLuminosityBlock(const edm::EventSetup &) override
called at begin of luminosity block (resets Mille binary in mille mode)
std::vector< align::RunNumber > cachedRuns_
std::string doDice(const fidpoint_t &fidpointvec, const idPair_t &id, const bool rotate=false)
AlignmentUserVariables * userVariables(void) const
Get pointer to user variables.
int callMille1D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iTrajHit, const std::vector< int > &globalLabels, const std::vector< float > &globalDerivativesX)
calls Mille for 1D hits
void addLaserData(const EventInfo &eventInfo, const TkFittedLasBeamCollection &tkLasBeams, const TsosVectorCollection &tkLasBeamTsoses)
int addMeasurementData(const edm::EventSetup &setup, const EventInfo &eventInfo, const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iHit, AlignmentParameters *&params)
bool isMode(unsigned int testMode) const
void beginRun(const edm::Run &run, const edm::EventSetup &setup, bool changed) override
called at begin of run
align::RotationType toLocal(const align::RotationType &) const
Return in local frame a rotation given in global frame.
virtual AlgebraicMatrix derivatives(const TrajectoryStateOnSurface &tsos, const AlignableDetOrUnitPtr &alidet) const =0
Get derivatives of selected parameters.
bool addHits(const align::Alignables &alis, const std::vector< AlignmentUserVariables * > &mpVars) const
CLHEP::HepMatrix AlgebraicMatrix
const std::vector< Scalar > & parameters() const
parallel to derivatives()
RunRanges makeNonOverlappingRunRanges(const edm::VParameterSet &runRanges, const RunNumber &defaultRun)
Definition: Utilities.cc:202
void makeGlobDerivMatrix(const std::vector< float > &globalDerivativesx, const std::vector< float > &globalDerivativesy, Eigen::MatrixBase< GlobalDerivativeMatrix > &aGlobalDerivativesM)
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
void writeOrigRigidBodyAlignmentParameters(const align::Alignables &alivec, const char *filename, int iter, bool validCheck, int &ierr) override
write RigidBodyAlignmentParameters as applied on top of original positions
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoToken_
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
T sqrt(T t)
Definition: SSEVec.h:19
RunRanges makeUniqueRunRanges(const edm::VParameterSet &runRanges, const RunNumber &defaultRun)
Definition: Utilities.cc:241
unsigned int getBeamId(void) const
return the full beam identifier
Definition: TkLasBeam.h:23
AlignmentParameterStore * theAlignmentParameterStore
directory for all kind of files
unsigned int hitsY() const
get number of hits for y-measurement
bool is2D(const TransientTrackingRecHit::ConstRecHitPointer &recHit) const
true if hit belongs to 2D detector (currently tracker specific)
tuple key
prepare the HTCondor submission files and eventually submit them
std::vector< AlignmentUserVariables * > readMillePedeVariables(const align::Alignables &alivec, const char *filename, int iter, int &ierr)
Transition
Definition: Transition.h:12
void addRefTrackVirtualMeas1D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iVirtualMeas, Eigen::MatrixBase< CovarianceMatrix > &aHitCovarianceM, Eigen::MatrixBase< ResidualMatrix > &aHitResidualsM, Eigen::MatrixBase< LocalDerivativeMatrix > &aLocalDerivativesM)
adds data for a specific virtual measurement from reference trajectory
void addRefTrackData2D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iTrajHit, Eigen::MatrixBase< CovarianceMatrix > &aHitCovarianceM, Eigen::MatrixBase< ResidualMatrix > &aHitResidualsM, Eigen::MatrixBase< LocalDerivativeMatrix > &aLocalDerivativesM)
adds data from reference trajectory from a specific Hit
tuple reader
Definition: DQM.py:105
const TkFittedLasBeamCollection * tkLasBeams() const
std::unique_ptr< PedeSteerer > thePedeSteer
bool storeAlignments() override
Returns whether MP produced results to be stored.
void terminate() override
Called at end of job.
void initialize(const edm::EventSetup &setup, AlignableTracker *tracker, AlignableMuon *muon, AlignableExtras *extras, AlignmentParameterStore *store) override
Called at beginning of job.
const AlignableSurface & surface() const
Return the Surface (global position and orientation) of the object.
Definition: Alignable.h:132
std::vector< ConstRecHitPointer > ConstRecHitContainer
bool addHitStatistics(int fromLoop, const std::string &outFile, const std::vector< std::string > &inFiles) const
Log< level::Info, false > LogInfo
void addVirtualMeas(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iVirtualMeas)
adds data for virtual measurements from reference trajectory
bool readFromPede(const edm::ParameterSet &mprespset, bool setUserVars, const RunRange &runrange)
read pede input defined by &#39;psetName&#39;, flag to create/not create MillePedeVariables ...
std::pair< unsigned int, unsigned int > addReferenceTrajectory(const edm::EventSetup &setup, const EventInfo &eventInfo, const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr)
fill mille for a trajectory, returning number of x/y hits ([0,0] if &#39;bad&#39; trajectory) ...
const SiStripDetId & getDetId(void) const
virtual bool storeThresholds(const int &nRecords, const AlignPCLThresholds::threshold_map &thresholdMap)
void setUserVariables(AlignmentUserVariables *auv)
Set pointer to user variables.
int size(void) const
Get number of parameters.
void addUntrackedParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:192
std::vector< IntegratedCalibrationBase * > theCalibrations
Class to hold one picture of the BPix survey.
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< value_t > localpars_t
ValidityInterval validityInterval() const
void writeAlignableAbsolutePositions(const align::Alignables &alivec, const char *filename, int iter, bool validCheck, int &ierr) override
write Alignable current absolute positions
std::vector< Alignable * > Alignables
Definition: Utilities.h:31
void run(const edm::EventSetup &setup, const EventInfo &eventInfo) override
Run the algorithm on trajectories and tracks.
~MillePedeAlignmentAlgorithm() override
Destructor.
MillePedeAlignmentAlgorithm(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
Constructor.
std::vector< TkFittedLasBeam > TkFittedLasBeamCollection
double a
Definition: hdecay.h:119
const reco::BeamSpot & beamSpot() const
T get() const
Definition: EventSetup.h:82
align::GlobalPoints toGlobal(const align::LocalPoints &) const
Return in global coord given a set of local points.
void cacheTransformations(void)
cache the current position, rotation and other parameters
tuple filename
Definition: lut2db_cfg.py:20
const edm::ESGetToken< AlignPCLThresholds, AlignPCLThresholdsRcd > aliThrToken_
TransientTrackingRecHit::ConstRecHitContainer ConstRecHitContainer
std::string fullPath() const
Definition: FileInPath.cc:161
static const count_t nMsrmts
void addPxbSurvey(const edm::ParameterSet &pxbSurveyCfg)
add measurement data from PXB survey
const IOVSyncValue & first() const
unsigned int decodeMode(const std::string &mode) const
#define get
std::unique_ptr< AlignableNavigator > theAlignableNavigator
Log< level::Warning, false > LogWarning
Time_t endValue
Definition: Time.h:42
const std::vector< SiStripLaserRecHit2D > & getData(void) const
access the collection of hits
Definition: TkLasBeam.h:26
std::vector< std::vector< TrajectoryStateOnSurface > > TsosVectorCollection
define run information passed to algorithms (in endRun)
const TsosVectorCollection * tkLasBeamTsoses() const
might be null!
const AlgebraicSymMatrix & covariance(void) const
Get parameter covariance matrix.
Constructor of the full muon geometry.
Definition: AlignableMuon.h:38
uint32_t dimension(pat::CandKinResolution::Parametrization parametrization)
Returns the number of free parameters in a parametrization (3 or 4)
bool read(align::Alignables &alignables, bool setUserVars)
Definition: PedeReader.cc:56
Alignable * mother() const
Return pointer to container alignable (if any)
Definition: Alignable.h:91
const Time_t MAX_VAL(std::numeric_limits< Time_t >::max())
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
std::vector< ReferenceTrajectoryPtr > ReferenceTrajectoryCollection
void throwException(const std::string &message, const std::string &methodName)
Definition: Exception.cc:12
bool setAllDefault(unsigned int nParam)
set default values for all data concerning nParam (false if nParam out of range)
static constexpr auto TEC
Definition: SiStripDetId.h:40
Definition: Run.h:45