CMS 3D CMS Logo

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