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 //#include "MillePedeAlignmentAlgorithm.h"
13 
17 
19 // in header, too
20 // end in header, too
21 
25 #include "Mille.h" // 'unpublished' interface located in src
26 #include "PedeSteerer.h" // dito
27 #include "PedeReader.h" // dito
30 
33 
37 
41 
42 // includes to make known that they inherit from Alignable:
46 
52 
58 
60 
61 #include <fstream>
62 #include <sstream>
63 #include <algorithm>
64 #include <sys/stat.h>
65 
66 #include <TMath.h>
67 #include <TMatrixDSymEigen.h>
71 
72 // Includes for PXB survey
73 #include <iostream>
80 
82 
83 using namespace gbl;
84 
85 
86 // Constructor ----------------------------------------------------------------
87 //____________________________________________________
90  theConfig(cfg),
91  theMode(this->decodeMode(theConfig.getUntrackedParameter<std::string>("mode"))),
92  theDir(theConfig.getUntrackedParameter<std::string>("fileDir")),
93  theAlignmentParameterStore(nullptr),
94  theAlignables(),
95  theMinNumHits(cfg.getParameter<unsigned int>("minNumHits")),
96  theMaximalCor2D(cfg.getParameter<double>("max2Dcorrelation")),
97  theLastWrittenIov(0),
98  theGblDoubleBinary(cfg.getParameter<bool>("doubleBinary")),
99  runAtPCL_(cfg.getParameter<bool>("runAtPCL"))
100 {
101  if (!theDir.empty() && theDir.find_last_of('/') != theDir.size()-1) theDir += '/';// may need '/'
102  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm" << "Start in mode '"
104  << "' with output directory '" << theDir << "'.";
105  if (this->isMode(myMilleBit)) {
106  theMille = std::make_unique<Mille>((theDir + theConfig.getParameter<std::string>("binaryFile")).c_str());// add ', false);' for text output);
107  // use same file for GBL
108  theBinary = std::make_unique<MilleBinary>((theDir + theConfig.getParameter<std::string>("binaryFile")).c_str(), theGblDoubleBinary);
109  }
110 }
111 
112 // Destructor ----------------------------------------------------------------
113 //____________________________________________________
115 {
116 }
117 
118 // Call at beginning of job ---------------------------------------------------
119 //____________________________________________________
123 {
124  if (muon) {
125  edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::initialize"
126  << "Running with AlignabeMuon not yet tested.";
127  }
128 
129  //Retrieve tracker topology from geometry
130  edm::ESHandle<TrackerTopology> tTopoHandle;
131  setup.get<TrackerTopologyRcd>().get(tTopoHandle);
132  const TrackerTopology* const tTopo = tTopoHandle.product();
133 
134  theAlignableNavigator = std::make_unique<AlignableNavigator>(extras, tracker, muon);
136  theAlignables = theAlignmentParameterStore->alignables();
137 
138  edm::ParameterSet pedeLabelerCfg(theConfig.getParameter<edm::ParameterSet>("pedeLabeler"));
139  edm::VParameterSet RunRangeSelectionVPSet(theConfig.getUntrackedParameter<edm::VParameterSet>("RunRangeSelection"));
140  pedeLabelerCfg.addUntrackedParameter<edm::VParameterSet>("RunRangeSelection",
141  RunRangeSelectionVPSet);
142 
143  std::string labelerPlugin = "PedeLabeler";
144  if (RunRangeSelectionVPSet.size()>0) {
145  labelerPlugin = "RunRangeDependentPedeLabeler";
146  if (pedeLabelerCfg.exists("plugin")) {
147  std::string labelerPluginCfg = pedeLabelerCfg.getParameter<std::string>("plugin");
148  if ((labelerPluginCfg!="PedeLabeler" && labelerPluginCfg!="RunRangeDependentPedeLabeler") ||
149  pedeLabelerCfg.getUntrackedParameter<edm::VParameterSet>("parameterInstances").size()>0) {
150  throw cms::Exception("BadConfig")
151  << "MillePedeAlignmentAlgorithm::initialize"
152  << "both RunRangeSelection and generic labeler specified in config file. "
153  << "Please get rid of either one of them.\n";
154  }
155  }
156  } else {
157  if (pedeLabelerCfg.exists("plugin")) {
158  labelerPlugin = pedeLabelerCfg.getParameter<std::string>("plugin");
159  }
160  }
161 
162  if (!pedeLabelerCfg.exists("plugin")) {
163  pedeLabelerCfg.addUntrackedParameter<std::string>("plugin", labelerPlugin);
164  }
165 
166  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::initialize"
167  << "Using plugin '" << labelerPlugin << "' to generate labels.";
168 
169  thePedeLabels = std::unique_ptr<PedeLabelerBase>(PedeLabelerPluginFactory::get()
170  ->create(labelerPlugin,
171  PedeLabelerBase::TopLevelAlignables(tracker, muon, extras),
172  pedeLabelerCfg));
173 
174  // 1) Create PedeSteerer: correct alignable positions for coordinate system selection
175  edm::ParameterSet pedeSteerCfg(theConfig.getParameter<edm::ParameterSet>("pedeSteerer"));
176  thePedeSteer = std::make_unique<PedeSteerer>(tracker, muon, extras,
178  pedeSteerCfg, theDir, !this->isMode(myPedeSteerBit));
179 
180  // 2) If requested, directly read in and apply result of previous pede run,
181  // assuming that correction from 1) was also applied to create the result:
182  const std::vector<edm::ParameterSet> mprespset
183  (theConfig.getParameter<std::vector<edm::ParameterSet> >("pedeReaderInputs"));
184  if (!mprespset.empty()) {
185  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::initialize"
186  << "Apply " << mprespset.end() - mprespset.begin()
187  << " previous MillePede constants from 'pedeReaderInputs'.";
188  }
189 
190  // FIXME: add selection of run range via 'pedeReaderInputs'
191  // Note: Results for parameters of IntegratedCalibration's cannot be treated...
192  RunRange runrange(cond::timeTypeSpecs[cond::runnumber].beginValue,
194  for (std::vector<edm::ParameterSet>::const_iterator iSet = mprespset.begin(), iE = mprespset.end();
195  iSet != iE; ++iSet) {
196  // This read will ignore calibrations as long as they are not yet passed to Millepede
197  // during/before initialize(..) - currently addCalibrations(..) is called later in AlignmentProducer
198  if (!this->readFromPede((*iSet), false, runrange)) { // false: do not erase SelectionUserVariables
199  throw cms::Exception("BadConfig")
200  << "MillePedeAlignmentAlgorithm::initialize: Problems reading input constants of "
201  << "pedeReaderInputs entry " << iSet - mprespset.begin() << '.';
202  }
203  theAlignmentParameterStore->applyParameters();
204  // Needed to shut up later warning from checkAliParams:
205  theAlignmentParameterStore->resetParameters();
206  }
207 
208  // 3) Now create steerings with 'final' start position:
209  thePedeSteer->buildSubSteer(tracker, muon, extras);
210 
211  // After (!) 1-3 of PedeSteerer which uses the SelectionUserVariables attached to the parameters:
212  this->buildUserVariables(theAlignables); // for hit statistics and/or pede result
213 
214  if (this->isMode(myMilleBit)) {
215  if (!theConfig.getParameter<std::vector<std::string> >("mergeBinaryFiles").empty() ||
216  !theConfig.getParameter<std::vector<std::string> >("mergeTreeFiles").empty()) {
217  throw cms::Exception("BadConfig")
218  << "'vstring mergeTreeFiles' and 'vstring mergeBinaryFiles' must be empty for "
219  << "modes running mille.";
220  }
221  const std::string moniFile(theConfig.getUntrackedParameter<std::string>("monitorFile"));
222  if (moniFile.size()) theMonitor = std::make_unique<MillePedeMonitor>(tTopo, (theDir + moniFile).c_str());
223 
224  // Get trajectory factory. In case nothing found, FrameWork will throw...
225  const edm::ParameterSet fctCfg(theConfig.getParameter<edm::ParameterSet>("TrajectoryFactory"));
226  const std::string fctName(fctCfg.getParameter<std::string>("TrajectoryFactoryName"));
227  theTrajectoryFactory = std::unique_ptr<TrajectoryFactoryBase>(TrajectoryFactoryPlugin::get()
228  ->create(fctName, fctCfg));
229  }
230 
231  if (this->isMode(myPedeSteerBit)) {
232  // Get config for survey and set flag accordingly
233  const edm::ParameterSet pxbSurveyCfg(theConfig.getParameter<edm::ParameterSet>("surveyPixelBarrel"));
234  theDoSurveyPixelBarrel = pxbSurveyCfg.getParameter<bool>("doSurvey");
235  if (theDoSurveyPixelBarrel) this->addPxbSurvey(pxbSurveyCfg);
236  }
237 }
238 
239 //____________________________________________________
241  return true;
242 }
243 
244 //____________________________________________________
245 bool MillePedeAlignmentAlgorithm::addCalibrations(const std::vector<IntegratedCalibrationBase*> &iCals)
246 {
247  theCalibrations.insert(theCalibrations.end(), iCals.begin(), iCals.end());
248  thePedeLabels->addCalibrations(iCals);
249  return true;
250 }
251 
252 //_____________________________________________________________________________
254 {
255  if (isMode(myMilleBit)) {
256  return true;
257  } else {
258  return false;
259  }
260 }
261 
262 //____________________________________________________
264 {
265  if (this->isMode(myPedeReadBit)) {
266  // restore initial positions, rotations and deformations
268 
269  // Needed to shut up later warning from checkAliParams:
271  // To avoid that they keep values from previous IOV if no new one in pede result
273 
274  if (!this->readFromPede(theConfig.getParameter<edm::ParameterSet>("pedeReader"), true, runrange)) {
275  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::setParametersForRunRange"
276  << "Problems reading pede result, but applying!";
277  }
279 
280  this->doIO(++theLastWrittenIov); // pre-increment!
281  }
282 
283  return true;
284 }
285 
286 // Call at end of job ---------------------------------------------------------
287 //____________________________________________________
289 {
290  terminate();
291 }
293 {
294  theMille.reset(); // delete to close binary before running pede below (flush would be enough...)
295  theBinary.reset();
296 
297  std::vector<std::string> files;
298  if (this->isMode(myMilleBit) || !theConfig.getParameter<std::string>("binaryFile").empty()) {
299  files.push_back(theDir + theConfig.getParameter<std::string>("binaryFile"));
300  } else {
301  const std::vector<std::string> plainFiles(theConfig.getParameter<std::vector<std::string> >("mergeBinaryFiles"));
302  files = getExistingFormattedFiles(plainFiles, theDir);
303  // Do some logging:
304  std::string filesForLogOutput;
305  for (const auto& file: files) filesForLogOutput += " " + file + ",";
306  if (filesForLogOutput.length() != 0) filesForLogOutput.pop_back();
307  edm::LogInfo("Alignment")
308  << "Based on the config parameter mergeBinaryFiles, using the following "
309  << "files as input (assigned weights are indicated by ' -- <weight>'):"
310  << filesForLogOutput;
311  }
312 
313  // cache all positions, rotations and deformations
315 
316  const std::string masterSteer(thePedeSteer->buildMasterSteer(files));// do only if myPedeSteerBit?
317  if (this->isMode(myPedeRunBit)) {
318  thePedeSteer->runPede(masterSteer);
319  }
320 
321  // parameters from pede are not yet applied,
322  // so we can still write start positions (but with hit statistics in case of mille):
323  this->doIO(0);
324  theLastWrittenIov = 0;
325 }
326 
327 std::vector<std::string> MillePedeAlignmentAlgorithm::getExistingFormattedFiles(const std::vector<std::string>& plainFiles, const std::string& theDir) {
328  std::vector<std::string> files;
329  for (const auto& plainFile: plainFiles) {
330  std::string theInputFileName = plainFile;
331  int theNumber = 0;
332  while (true) {
333  // Create a formatted version of the filename, with growing numbers
334  // If the parameter doesn't contain a formatting directive, it just stays unchanged
335  char theNumberedInputFileName[200];
336  sprintf(theNumberedInputFileName, theInputFileName.c_str(), theNumber);
337  std::string theCompleteInputFileName = theDir + theNumberedInputFileName;
338  const auto endOfStrippedFileName = theCompleteInputFileName.rfind(" --");
339  const auto strippedInputFileName = theCompleteInputFileName.substr(0, endOfStrippedFileName);
340  // Check if the file exists
341  struct stat buffer;
342  if (stat (strippedInputFileName.c_str(), &buffer) == 0) {
343  // If the file exists, add it to the list
344  files.push_back(theCompleteInputFileName);
345  if (theNumberedInputFileName == theInputFileName) {
346  // If the filename didn't contain a formatting directive, no reason to look any further, break out of the loop
347  break;
348  } else {
349  // Otherwise look for the next number
350  theNumber++;
351  }
352  } else {
353  // The file doesn't exist, break out of the loop
354  break;
355  }
356  }
357  // warning if unformatted (-> theNumber stays at 0) does not exist
358  if (theNumber == 0 && (files.size() == 0 || files.back() != plainFile)) {
359  edm::LogWarning("Alignment")
360  << "The input file '" << plainFile << "' does not exist.";
361  }
362  }
363  return files;
364 }
365 
366 // Run the algorithm on trajectories and tracks -------------------------------
367 //____________________________________________________
369 {
370  if (!this->isMode(myMilleBit)) return; // no theMille created...
372 
373  if (theMonitor) { // monitor input tracks
374  for (ConstTrajTrackPairCollection::const_iterator iTrajTrack = tracks.begin();
375  iTrajTrack != tracks.end(); ++iTrajTrack) {
376  theMonitor->fillTrack((*iTrajTrack).second);
377  }
378  }
379 
380  const RefTrajColl trajectories(theTrajectoryFactory->trajectories(setup, tracks, eventInfo.beamSpot()));
381 
382  // Now loop over ReferenceTrajectoryCollection
383  unsigned int refTrajCount = 0; // counter for track monitoring if 1 track per trajectory
384  for (RefTrajColl::const_iterator iRefTraj = trajectories.begin(), iRefTrajE = trajectories.end();
385  iRefTraj != iRefTrajE; ++iRefTraj, ++refTrajCount) {
386 
387  RefTrajColl::value_type refTrajPtr = *iRefTraj;
388  if (theMonitor) theMonitor->fillRefTrajectory(refTrajPtr);
389 
390  const std::pair<unsigned int, unsigned int> nHitXy
391  = this->addReferenceTrajectory(setup, eventInfo, refTrajPtr);
392 
393  if (theMonitor && (nHitXy.first || nHitXy.second)) {
394  // if track used (i.e. some hits), fill monitoring
395  // track NULL ptr if trajectories and tracks do not match
396  const reco::Track *trackPtr =
397  (trajectories.size() == tracks.size() ? tracks[refTrajCount].second : 0);
398  theMonitor->fillUsedTrack(trackPtr, nHitXy.first, nHitXy.second);
399  }
400 
401  } // end of reference trajectory and track loop
402 }
403 
404 //____________________________________________________
405 std::pair<unsigned int, unsigned int>
407  const EventInfo &eventInfo,
408  const RefTrajColl::value_type &refTrajPtr)
409 {
410  std::pair<unsigned int, unsigned int> hitResultXy(0,0);
411  if (refTrajPtr->isValid()) {
412 
413 
414  // GblTrajectory?
415  if (refTrajPtr->gblInput().size() > 0) {
416  // by construction: number of GblPoints == number of recHits or == zero !!!
417  unsigned int iHit = 0;
418  unsigned int numPointsWithMeas = 0;
419  std::vector<GblPoint>::iterator itPoint;
420  std::vector<std::pair<std::vector<GblPoint>, TMatrixD> > theGblInput = refTrajPtr->gblInput();
421  for (unsigned int iTraj = 0; iTraj < refTrajPtr->gblInput().size(); ++iTraj) {
422  for (itPoint = refTrajPtr->gblInput()[iTraj].first.begin(); itPoint < refTrajPtr->gblInput()[iTraj].first.end(); ++itPoint) {
423  if (this->addGlobalData(setup, eventInfo, refTrajPtr, iHit++, *itPoint) < 0) return hitResultXy;
424  if (itPoint->hasMeasurement() >= 1) ++numPointsWithMeas;
425  }
426  }
427  hitResultXy.first = numPointsWithMeas;
428  // check #hits criterion
429  if (hitResultXy.first == 0 || hitResultXy.first < theMinNumHits) return hitResultXy;
430  // construct GBL trajectory
431  if (refTrajPtr->gblInput().size() == 1) {
432  // from single track
433  GblTrajectory aGblTrajectory( refTrajPtr->gblInput()[0].first, refTrajPtr->nominalField() != 0 );
434  // GBL fit trajectory
435  /*double Chi2;
436  int Ndf;
437  double lostWeight;
438  aGblTrajectory.fit(Chi2, Ndf, lostWeight);
439  std::cout << " GblFit: " << Chi2 << ", " << Ndf << ", " << lostWeight << std::endl; */
440  // write to MP binary file
441  if (aGblTrajectory.isValid() && aGblTrajectory.getNumPoints() >= theMinNumHits) aGblTrajectory.milleOut(*theBinary);
442  }
443  if (refTrajPtr->gblInput().size() == 2) {
444  // from TwoBodyDecay
445  GblTrajectory aGblTrajectory( refTrajPtr->gblInput(), refTrajPtr->gblExtDerivatives(), refTrajPtr->gblExtMeasurements(), refTrajPtr->gblExtPrecisions() );
446  // write to MP binary file
447  if (aGblTrajectory.isValid() && aGblTrajectory.getNumPoints() >= theMinNumHits) aGblTrajectory.milleOut(*theBinary);
448  }
449  } else {
450  // to add hits if all fine:
451  std::vector<AlignmentParameters*> parVec(refTrajPtr->recHits().size());
452  // collect hit statistics, assuming that there are no y-only hits
453  std::vector<bool> validHitVecY(refTrajPtr->recHits().size(), false);
454  // Use recHits from ReferenceTrajectory (since they have the right order!):
455  for (unsigned int iHit = 0; iHit < refTrajPtr->recHits().size(); ++iHit) {
456  const int flagXY = this->addMeasurementData(setup, eventInfo, refTrajPtr, iHit, parVec[iHit]);
457 
458  if (flagXY < 0) { // problem
459  hitResultXy.first = 0;
460  break;
461  } else { // hit is fine, increase x/y statistics
462  if (flagXY >= 1) ++hitResultXy.first;
463  validHitVecY[iHit] = (flagXY >= 2);
464  }
465  } // end loop on hits
466 
467  // add virtual measurements
468  for (unsigned int iVirtualMeas = 0; iVirtualMeas < refTrajPtr->numberOfVirtualMeas(); ++iVirtualMeas) {
469  this->addVirtualMeas(refTrajPtr, iVirtualMeas);
470  }
471 
472  // kill or end 'track' for mille, depends on #hits criterion
473  if (hitResultXy.first == 0 || hitResultXy.first < theMinNumHits) {
474  theMille->kill();
475  hitResultXy.first = hitResultXy.second = 0; //reset
476  } else {
477  theMille->end();
478  // add x/y hit count to MillePedeVariables of parVec,
479  // returning number of y-hits of the reference trajectory
480  hitResultXy.second = this->addHitCount(parVec, validHitVecY);
481  //
482  }
483  }
484 
485  } // end if valid trajectory
486 
487  return hitResultXy;
488 }
489 
490 //____________________________________________________
491 unsigned int
492 MillePedeAlignmentAlgorithm::addHitCount(const std::vector<AlignmentParameters*> &parVec,
493  const std::vector<bool> &validHitVecY) const
494 {
495  // Loop on all hit information in the input arrays and count valid y-hits:
496  unsigned int nHitY = 0;
497  for (unsigned int iHit = 0; iHit < validHitVecY.size(); ++iHit) {
498  Alignable *ali = (parVec[iHit] ? parVec[iHit]->alignable() : 0);
499  // Loop upwards on hierarchy of alignables to add hits to all levels
500  // that are currently aligned. If only a non-selected alignable was hit,
501  // (i.e. flagXY == 0 in addReferenceTrajectory(..)), there is no loop at all...
502  while (ali) {
504  if (pars) { // otherwise hierarchy level not selected
505  // cast ensured by previous checks:
506  MillePedeVariables *mpVar = static_cast<MillePedeVariables*>(pars->userVariables());
507  // every hit has an x-measurement, cf. addReferenceTrajectory(..):
508  mpVar->increaseHitsX();
509  if (validHitVecY[iHit]) {
510  mpVar->increaseHitsY();
511  if (pars == parVec[iHit]) ++nHitY; // do not count hits twice
512  }
513  }
514  ali = ali->mother();
515  }
516  }
517 
518  return nHitY;
519 }
520 
521 //____________________________________________________
523  const edm::EventSetup &setup)
524 {
525  if(runInfo.tkLasBeams() && runInfo.tkLasBeamTsoses()){
526  // LAS beam treatment
527  this->addLaserData(eventInfo, *(runInfo.tkLasBeams()), *(runInfo.tkLasBeamTsoses()));
528  }
529  if(this->isMode(myMilleBit)) theMille->flushOutputFile();
530 }
531 
532 // Implementation of endRun that DOES get called. (Because we need it.)
534  if(this->isMode(myMilleBit)) theMille->flushOutputFile();
535 }
536 
537 //____________________________________________________
539 {
540  if (!runAtPCL_) return;
541  if(this->isMode(myMilleBit)) theMille->resetOutputFile();
542 }
543 
544 //____________________________________________________
546 {
547  if (!runAtPCL_) return;
548  if(this->isMode(myMilleBit)) theMille->flushOutputFile();
549 }
550 
551 
552 //____________________________________________________
554  const EventInfo &eventInfo,
556  unsigned int iHit,
557  AlignmentParameters *&params)
558 {
559  params = 0;
560  theFloatBufferX.clear();
561  theFloatBufferY.clear();
562  theIntBuffer.clear();
563 
564  const TrajectoryStateOnSurface &tsos = refTrajPtr->trajectoryStates()[iHit];
565  const ConstRecHitPointer &recHitPtr = refTrajPtr->recHits()[iHit];
566  // ignore invalid hits
567  if (!recHitPtr->isValid()) return 0;
568 
569  // First add the derivatives from IntegratedCalibration's,
570  // should even be OK if problems for "usual" derivatives from Alignables
571  this->globalDerivativesCalibration(recHitPtr, tsos, setup, eventInfo, // input
573 
574  // get AlignableDet/Unit for this hit
575  AlignableDetOrUnitPtr alidet(theAlignableNavigator->alignableFromDetId(recHitPtr->geographicalId()));
576 
577  if (!this->globalDerivativesHierarchy(eventInfo,
578  tsos, alidet, alidet, theFloatBufferX, // 2x alidet, sic!
579  theFloatBufferY, theIntBuffer, params)) {
580  return -1; // problem
581  } else if (theFloatBufferX.empty()) {
582  return 0; // empty for X: no alignable for hit, nor calibrations
583  } else { // now even if no alignable, but calibrations!
584  return this->callMille(refTrajPtr, iHit, theIntBuffer, theFloatBufferX, theFloatBufferY);
585  }
586 }
587 
588 //____________________________________________________
589 
592  unsigned int iHit, GblPoint &gblPoint)
593 {
594  AlignmentParameters* params = 0;
595  std::vector<double> theDoubleBufferX, theDoubleBufferY;
596  theDoubleBufferX.clear();
597  theDoubleBufferY.clear();
598  theIntBuffer.clear();
599  int iret = 0;
600 
601  const TrajectoryStateOnSurface &tsos = refTrajPtr->trajectoryStates()[iHit];
602  const ConstRecHitPointer &recHitPtr = refTrajPtr->recHits()[iHit];
603  // ignore invalid hits
604  if (!recHitPtr->isValid()) return 0;
605 
606  // get AlignableDet/Unit for this hit
607  AlignableDetOrUnitPtr alidet(theAlignableNavigator->alignableFromDetId(recHitPtr->geographicalId()));
608 
609  if (!this->globalDerivativesHierarchy(eventInfo,
610  tsos, alidet, alidet, theDoubleBufferX, // 2x alidet, sic!
611  theDoubleBufferY, theIntBuffer, params)) {
612  return -1; // problem
613  }
614  //calibration parameters
615  int globalLabel;
616  std::vector<IntegratedCalibrationBase::ValuesIndexPair> derivs;
617  for (auto iCalib = theCalibrations.begin(); iCalib != theCalibrations.end(); ++iCalib) {
618  // get all derivatives of this calibration // const unsigned int num =
619  (*iCalib)->derivatives(derivs, *recHitPtr, tsos, setup, eventInfo);
620  for (auto iValuesInd = derivs.begin(); iValuesInd != derivs.end(); ++iValuesInd) {
621  // transfer label and x/y derivatives
622  globalLabel = thePedeLabels->calibrationLabel(*iCalib, iValuesInd->second);
623  if (globalLabel > 0 && globalLabel <= 2147483647) {
624  theIntBuffer.push_back(globalLabel);
625  theDoubleBufferX.push_back(iValuesInd->first.first);
626  theDoubleBufferY.push_back(iValuesInd->first.second);
627  } else {
628  std::cerr << "MillePedeAlignmentAlgorithm::addGlobalData: Invalid label " << globalLabel << " <= 0 or > 2147483647" << std::endl;
629  }
630  }
631  }
632  unsigned int numGlobals = theIntBuffer.size();
633  if (numGlobals > 0)
634  {
635  TMatrixD globalDer(2,numGlobals);
636  for (unsigned int i = 0; i < numGlobals; ++i) {
637  globalDer(0,i) = theDoubleBufferX[i];
638  globalDer(1,i) = theDoubleBufferY[i];
639  }
640  gblPoint.addGlobals( theIntBuffer, globalDer );
641  iret = 1;
642  }
643  return iret;
644 }
645 
646 //____________________________________________________
649  const TrajectoryStateOnSurface &tsos,
650  Alignable *ali, const AlignableDetOrUnitPtr &alidet,
651  std::vector<float> &globalDerivativesX,
652  std::vector<float> &globalDerivativesY,
653  std::vector<int> &globalLabels,
654  AlignmentParameters *&lowestParams) const
655 {
656  // derivatives and labels are recursively attached
657  if (!ali) return true; // no mother might be OK
658 
659  if (false && theMonitor && alidet != ali) theMonitor->fillFrameToFrame(alidet, ali);
660 
661  AlignmentParameters *params = ali->alignmentParameters();
662 
663  if (params) {
664  if (!lowestParams) lowestParams = params; // set parameters of lowest level
665 
666  bool hasSplitParameters = thePedeLabels->hasSplitParameters(ali);
667  const unsigned int alignableLabel = thePedeLabels->alignableLabel(ali);
668 
669  if (0 == alignableLabel) { // FIXME: what about regardAllHits in Markus' code?
670  edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::globalDerivativesHierarchy"
671  << "Label not found, skip Alignable.";
672  return false;
673  }
674 
675  const std::vector<bool> &selPars = params->selector();
676  const AlgebraicMatrix derivs(params->derivatives(tsos, alidet));
677 
678  // cols: 2, i.e. x&y, rows: parameters, usually RigidBodyAlignmentParameters::N_PARAM
679  for (unsigned int iSel = 0; iSel < selPars.size(); ++iSel) {
680  if (selPars[iSel]) {
681  globalDerivativesX.push_back(derivs[iSel][kLocalX]
682  /thePedeSteer->cmsToPedeFactor(iSel));
683  if (hasSplitParameters==true) {
684  globalLabels.push_back(thePedeLabels->parameterLabel(ali, iSel, eventInfo, tsos));
685  } else {
686  globalLabels.push_back(thePedeLabels->parameterLabel(alignableLabel, iSel));
687  }
688  globalDerivativesY.push_back(derivs[iSel][kLocalY]
689  /thePedeSteer->cmsToPedeFactor(iSel));
690  }
691  }
692  // Exclude mothers if Alignable selected to be no part of a hierarchy:
693  if (thePedeSteer->isNoHiera(ali)) return true;
694  }
695  // Call recursively for mother, will stop if mother == 0:
696  return this->globalDerivativesHierarchy(eventInfo,
697  tsos, ali->mother(), alidet,
698  globalDerivativesX, globalDerivativesY,
699  globalLabels, lowestParams);
700 }
701 
702 //____________________________________________________
705  const TrajectoryStateOnSurface &tsos,
706  Alignable *ali, const AlignableDetOrUnitPtr &alidet,
707  std::vector<double> &globalDerivativesX,
708  std::vector<double> &globalDerivativesY,
709  std::vector<int> &globalLabels,
710  AlignmentParameters *&lowestParams) const
711 {
712  // derivatives and labels are recursively attached
713  if (!ali) return true; // no mother might be OK
714 
715  if (false && theMonitor && alidet != ali) theMonitor->fillFrameToFrame(alidet, ali);
716 
717  AlignmentParameters *params = ali->alignmentParameters();
718 
719  if (params) {
720  if (!lowestParams) lowestParams = params; // set parameters of lowest level
721 
722  bool hasSplitParameters = thePedeLabels->hasSplitParameters(ali);
723  const unsigned int alignableLabel = thePedeLabels->alignableLabel(ali);
724 
725  if (0 == alignableLabel) { // FIXME: what about regardAllHits in Markus' code?
726  edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::globalDerivativesHierarchy"
727  << "Label not found, skip Alignable.";
728  return false;
729  }
730 
731  const std::vector<bool> &selPars = params->selector();
732  const AlgebraicMatrix derivs(params->derivatives(tsos, alidet));
733  int globalLabel;
734 
735  // cols: 2, i.e. x&y, rows: parameters, usually RigidBodyAlignmentParameters::N_PARAM
736  for (unsigned int iSel = 0; iSel < selPars.size(); ++iSel) {
737  if (selPars[iSel]) {
738  if (hasSplitParameters==true) {
739  globalLabel = thePedeLabels->parameterLabel(ali, iSel, eventInfo, tsos);
740  } else {
741  globalLabel = thePedeLabels->parameterLabel(alignableLabel, iSel);
742  }
743  if (globalLabel > 0 && globalLabel <= 2147483647) {
744  globalLabels.push_back(globalLabel);
745  globalDerivativesX.push_back(derivs[iSel][kLocalX] / thePedeSteer->cmsToPedeFactor(iSel));
746  globalDerivativesY.push_back(derivs[iSel][kLocalY] / thePedeSteer->cmsToPedeFactor(iSel));
747  } else {
748  std::cerr << "MillePedeAlignmentAlgorithm::globalDerivativesHierarchy: Invalid label " << globalLabel << " <= 0 or > 2147483647" << std::endl;
749  }
750  }
751  }
752  // Exclude mothers if Alignable selected to be no part of a hierarchy:
753  if (thePedeSteer->isNoHiera(ali)) return true;
754  }
755  // Call recursively for mother, will stop if mother == 0:
756  return this->globalDerivativesHierarchy(eventInfo,
757  tsos, ali->mother(), alidet,
758  globalDerivativesX, globalDerivativesY,
759  globalLabels, lowestParams);
760 }
761 
762 //____________________________________________________
765  const TrajectoryStateOnSurface &tsos,
766  const edm::EventSetup &setup, const EventInfo &eventInfo,
767  std::vector<float> &globalDerivativesX,
768  std::vector<float> &globalDerivativesY,
769  std::vector<int> &globalLabels) const
770 {
771  std::vector<IntegratedCalibrationBase::ValuesIndexPair> derivs;
772  for (auto iCalib = theCalibrations.begin(); iCalib != theCalibrations.end(); ++iCalib) {
773  // get all derivatives of this calibration // const unsigned int num =
774  (*iCalib)->derivatives(derivs, *recHit, tsos, setup, eventInfo);
775  for (auto iValuesInd = derivs.begin(); iValuesInd != derivs.end(); ++iValuesInd) {
776  // transfer label and x/y derivatives
777  globalLabels.push_back(thePedeLabels->calibrationLabel(*iCalib, iValuesInd->second));
778  globalDerivativesX.push_back(iValuesInd->first.first);
779  globalDerivativesY.push_back(iValuesInd->first.second);
780  }
781  }
782 }
783 
784 // //____________________________________________________
785 // void MillePedeAlignmentAlgorithm
786 // ::callMille(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
787 // unsigned int iTrajHit, MeasurementDirection xOrY,
788 // const std::vector<float> &globalDerivatives, const std::vector<int> &globalLabels)
789 // {
790 // const unsigned int xyIndex = iTrajHit*2 + xOrY;
791 // // FIXME: here for residuum and sigma we could use KALMAN-Filter results
792 // const float residuum =
793 // refTrajPtr->measurements()[xyIndex] - refTrajPtr->trajectoryPositions()[xyIndex];
794 // const float covariance = refTrajPtr->measurementErrors()[xyIndex][xyIndex];
795 // const float sigma = (covariance > 0. ? TMath::Sqrt(covariance) : 0.);
796 
797 // const AlgebraicMatrix &locDerivMatrix = refTrajPtr->derivatives();
798 
799 // std::vector<float> localDerivs(locDerivMatrix.num_col());
800 // for (unsigned int i = 0; i < localDerivs.size(); ++i) {
801 // localDerivs[i] = locDerivMatrix[xyIndex][i];
802 // }
803 
804 // // &(vector[0]) is valid - as long as vector is not empty
805 // // cf. http://www.parashift.com/c++-faq-lite/containers.html#faq-34.3
806 // theMille->mille(localDerivs.size(), &(localDerivs[0]),
807 // globalDerivatives.size(), &(globalDerivatives[0]), &(globalLabels[0]),
808 // residuum, sigma);
809 // if (theMonitor) {
810 // theMonitor->fillDerivatives(refTrajPtr->recHits()[iTrajHit],localDerivs, globalDerivatives,
811 // (xOrY == kLocalY));
812 // theMonitor->fillResiduals(refTrajPtr->recHits()[iTrajHit],
813 // refTrajPtr->trajectoryStates()[iTrajHit],
814 // iTrajHit, residuum, sigma, (xOrY == kLocalY));
815 // }
816 // }
817 
818 //____________________________________________________
820 {
821  // FIXME: Check whether this is a reliable and recommended way to find out...
822 
823  if (recHit->dimension() < 2) {
824  return false; // some muon and TIB/TOB stuff really has RecHit1D
825  } else if (recHit->detUnit()) { // detunit in strip is 1D, in pixel 2D
826  return recHit->detUnit()->type().isTrackerPixel();
827  } else { // stereo strips (FIXME: endcap trouble due to non-parallel strips (wedge sensors)?)
828  if (dynamic_cast<const ProjectedSiStripRecHit2D*>(recHit->hit())) { // check persistent hit
829  // projected: 1D measurement on 'glued' module
830  return false;
831  } else {
832  return true;
833  }
834  }
835 }
836 
837 //__________________________________________________________________________________________________
838 bool MillePedeAlignmentAlgorithm::readFromPede(const edm::ParameterSet &mprespset, bool setUserVars,
839  const RunRange &runrange)
840 {
841  bool allEmpty = this->areEmptyParams(theAlignables);
842 
843  PedeReader reader(mprespset, *thePedeSteer, *thePedeLabels, runrange);
844  std::vector<Alignable*> alis;
845  bool okRead = reader.read(alis, setUserVars); // also may set params of IntegratedCalibration's
846  bool numMatch = true;
847 
848  std::stringstream out;
849  out << "Read " << alis.size() << " alignables";
850  if (alis.size() != theAlignables.size()) {
851  out << " while " << theAlignables.size() << " in store";
852  numMatch = false; // FIXME: Should we check one by one? Or transfer 'alis' to the store?
853  }
854  if (!okRead) out << ", but problems in reading";
855  if (!allEmpty) out << ", possibly overwriting previous settings";
856  out << ".";
857 
858  if (okRead && allEmpty) {
859  if (numMatch) { // as many alignables with result as trying to align
860  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
861  } else if (alis.size()) { // dead module do not get hits and no pede result
862  edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
863  } else { // serious problem: no result read - and not all modules can be dead...
864  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
865  return false;
866  }
867  return true;
868  }
869  // the rest is not OK:
870  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
871  return false;
872 }
873 
874 //__________________________________________________________________________________________________
875 bool MillePedeAlignmentAlgorithm::areEmptyParams(const std::vector<Alignable*> &alignables) const
876 {
877 
878  for (std::vector<Alignable*>::const_iterator iAli = alignables.begin();
879  iAli != alignables.end(); ++iAli) {
880  const AlignmentParameters *params = (*iAli)->alignmentParameters();
881  if (params) {
882  const AlgebraicVector &parVec(params->parameters());
883  const AlgebraicMatrix &parCov(params->covariance());
884  for (int i = 0; i < parVec.num_row(); ++i) {
885  if (parVec[i] != 0.) return false;
886  for (int j = i; j < parCov.num_col(); ++j) {
887  if (parCov[i][j] != 0.) return false;
888  }
889  }
890  }
891  }
892 
893  return true;
894 }
895 
896 //__________________________________________________________________________________________________
898 {
899  unsigned int result = 0;
900 
901  const std::string outFilePlain(theConfig.getParameter<std::string>("treeFile"));
902  if (outFilePlain.empty()) {
903  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
904  << "treeFile parameter empty => skip writing for 'loop' " << loop;
905  return result;
906  }
907 
908  const std::string outFile(theDir + outFilePlain);
909 
910  AlignmentIORoot aliIO;
911  int ioerr = 0;
912  if (loop == 0) {
913  aliIO.writeAlignableOriginalPositions(theAlignables, outFile.c_str(), loop, false, ioerr);
914  if (ioerr) {
915  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
916  << "Problem " << ioerr << " in writeAlignableOriginalPositions";
917  ++result;
918  }
919  } else if (loop == 1) {
920  // only for first iov add hit counts, else 2x, 3x,... number of hits in IOV 2, 3,...
921  const std::vector<std::string> inFiles
922  (theConfig.getParameter<std::vector<std::string> >("mergeTreeFiles"));
923  const std::vector<std::string> binFiles
924  (theConfig.getParameter<std::vector<std::string> >("mergeBinaryFiles"));
925  if (inFiles.size() != binFiles.size()) {
926  edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
927  << "'vstring mergeTreeFiles' and 'vstring mergeBinaryFiles' "
928  << "differ in size";
929  }
930  this->addHitStatistics(0, outFile, inFiles); // add hit info from tree 0 in 'infiles'
931  }
932  MillePedeVariablesIORoot millePedeIO;
933  millePedeIO.writeMillePedeVariables(theAlignables, outFile.c_str(), loop, false, ioerr);
934  if (ioerr) {
935  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
936  << "Problem " << ioerr << " writing MillePedeVariables";
937  ++result;
938  }
939 
940  aliIO.writeOrigRigidBodyAlignmentParameters(theAlignables, outFile.c_str(), loop, false, ioerr);
941  if (ioerr) {
942  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO" << "Problem " << ioerr
943  << " in writeOrigRigidBodyAlignmentParameters, " << loop;
944  ++result;
945  }
946  aliIO.writeAlignableAbsolutePositions(theAlignables, outFile.c_str(), loop, false, ioerr);
947  if (ioerr) {
948  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO" << "Problem " << ioerr
949  << " in writeAlignableAbsolutePositions, " << loop;
950  ++result;
951  }
952 
953  return result;
954 }
955 
956 //__________________________________________________________________________________________________
957 void MillePedeAlignmentAlgorithm::buildUserVariables(const std::vector<Alignable*> &alis) const
958 {
959  for (std::vector<Alignable*>::const_iterator iAli = alis.begin(); iAli != alis.end(); ++iAli) {
960  AlignmentParameters *params = (*iAli)->alignmentParameters();
961  if (!params) {
962  throw cms::Exception("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::buildUserVariables"
963  << "No parameters for alignable";
964  }
965  MillePedeVariables *userVars = dynamic_cast<MillePedeVariables*>(params->userVariables());
966  if (userVars) { // Just re-use existing, keeping label and numHits:
967  for (unsigned int iPar = 0; iPar < userVars->size(); ++iPar) {
968  // if (params->hierarchyLevel() > 0) {
969  //std::cout << params->hierarchyLevel() << "\nBefore: " << userVars->parameter()[iPar];
970  // }
971  userVars->setAllDefault(iPar);
972  //std::cout << "\nAfter: " << userVars->parameter()[iPar] << std::endl;
973  }
974  } else { // Nothing yet or erase wrong type:
975  userVars = new MillePedeVariables(params->size(), thePedeLabels->alignableLabel(*iAli));
976  params->setUserVariables(userVars);
977  }
978  }
979 }
980 
981 //__________________________________________________________________________________________________
983 {
984  if (mode == "full") {
986  } else if (mode == "mille") {
987  return myMilleBit; // + myPedeSteerBit; // sic! Including production of steerig file. NO!
988  } else if (mode == "pede") {
990  } else if (mode == "pedeSteer") {
991  return myPedeSteerBit;
992  } else if (mode == "pedeRun") {
993  return myPedeSteerBit + myPedeRunBit + myPedeReadBit; // sic! Including steering and reading of result.
994  } else if (mode == "pedeRead") {
995  return myPedeReadBit;
996  }
997 
998  throw cms::Exception("BadConfig")
999  << "Unknown mode '" << mode
1000  << "', use 'full', 'mille', 'pede', 'pedeRun', 'pedeSteer' or 'pedeRead'.";
1001 
1002  return 0;
1003 }
1004 
1005 //__________________________________________________________________________________________________
1007  const std::vector<std::string> &inFiles) const
1008 {
1009  bool allOk = true;
1010  int ierr = 0;
1011  MillePedeVariablesIORoot millePedeIO;
1012  for (std::vector<std::string>::const_iterator iFile = inFiles.begin();
1013  iFile != inFiles.end(); ++iFile) {
1014  const std::string inFile(theDir + *iFile);
1015  const std::vector<AlignmentUserVariables*> mpVars =
1016  millePedeIO.readMillePedeVariables(theAlignables, inFile.c_str(), fromIov, ierr);
1017  if (ierr || !this->addHits(theAlignables, mpVars)) {
1018  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::addHitStatistics"
1019  << "Error " << ierr << " reading from " << inFile
1020  << ", tree " << fromIov << ", or problems in addHits";
1021  allOk = false;
1022  }
1023  for (std::vector<AlignmentUserVariables*>::const_iterator i = mpVars.begin();
1024  i != mpVars.end(); ++i) {
1025  delete *i; // clean created objects
1026  }
1027  }
1028 
1029  return allOk;
1030 }
1031 
1032 //__________________________________________________________________________________________________
1033 bool MillePedeAlignmentAlgorithm::addHits(const std::vector<Alignable*> &alis,
1034  const std::vector<AlignmentUserVariables*> &mpVars) const
1035 {
1036  bool allOk = (mpVars.size() == alis.size());
1037  std::vector<AlignmentUserVariables*>::const_iterator iUser = mpVars.begin();
1038  for (std::vector<Alignable*>::const_iterator iAli = alis.begin();
1039  iAli != alis.end() && iUser != mpVars.end(); ++iAli, ++iUser) {
1040  MillePedeVariables *mpVarNew = dynamic_cast<MillePedeVariables*>(*iUser);
1041  AlignmentParameters *ps = (*iAli)->alignmentParameters();
1042  MillePedeVariables *mpVarOld = (ps ? dynamic_cast<MillePedeVariables*>(ps->userVariables()) : 0);
1043  if (!mpVarNew || !mpVarOld || mpVarOld->size() != mpVarNew->size()) {
1044  allOk = false;
1045  continue; // FIXME error etc.?
1046  }
1047 
1048  mpVarOld->increaseHitsX(mpVarNew->hitsX());
1049  mpVarOld->increaseHitsY(mpVarNew->hitsY());
1050  }
1051 
1052  return allOk;
1053 }
1054 
1055 //__________________________________________________________________________________________________
1056 void MillePedeAlignmentAlgorithm::makeGlobDerivMatrix(const std::vector<float> &globalDerivativesx,
1057  const std::vector<float> &globalDerivativesy,
1058  TMatrixF &aGlobalDerivativesM)
1059 {
1060 
1061  for (unsigned int i = 0; i < globalDerivativesx.size(); ++i) {
1062  aGlobalDerivativesM(0,i) = globalDerivativesx[i];
1063  aGlobalDerivativesM(1,i) = globalDerivativesy[i];
1064  }
1065 }
1066 
1067 //__________________________________________________________________________________________________
1069 (TMatrixDSym &aHitCovarianceM, TMatrixF &aLocalDerivativesM, TMatrixF &aHitResidualsM,
1070  TMatrixF &aGlobalDerivativesM) const
1071 {
1072  TMatrixDSymEigen myDiag(aHitCovarianceM);
1073  TMatrixD aTranfoToDiagonalSystem = myDiag.GetEigenVectors();
1074  TMatrixD aTranfoToDiagonalSystemInv = myDiag.GetEigenVectors( );
1075  TMatrixF aTranfoToDiagonalSystemInvF = myDiag.GetEigenVectors( );
1076  TMatrixD aMatrix = aTranfoToDiagonalSystemInv.Invert() * aHitCovarianceM * aTranfoToDiagonalSystem;
1077  // Tranformation of matrix M is done by A^T*M*A, not A^{-1}*M*A.
1078  // But here A^T == A^{-1}, so we would only save CPU by Transpose()...
1079  // FIXME this - I guess simply use T(), not Transpose()...
1080  // TMatrixD aMatrix = aTranfoToDiagonalSystemInv.Transpose() * aHitCovarianceM
1081  // * aTranfoToDiagonalSystem;
1082  aHitCovarianceM = TMatrixDSym(2, aMatrix.GetMatrixArray());
1083  aTranfoToDiagonalSystemInvF.Invert();
1084  //edm::LogInfo("Alignment") << "NEW HIT loca in matrix"<<aLocalDerivativesM(0,0);
1085  aLocalDerivativesM = aTranfoToDiagonalSystemInvF * aLocalDerivativesM;
1086 
1087  //edm::LogInfo("Alignment") << "NEW HIT loca in matrix after diag:"<<aLocalDerivativesM(0,0);
1088  aHitResidualsM = aTranfoToDiagonalSystemInvF * aHitResidualsM;
1089  aGlobalDerivativesM = aTranfoToDiagonalSystemInvF * aGlobalDerivativesM;
1090 }
1091 
1092 //__________________________________________________________________________________________________
1095  unsigned int iVirtualMeas, TMatrixDSym &aHitCovarianceM,
1096  TMatrixF &aHitResidualsM, TMatrixF &aLocalDerivativesM)
1097 {
1098  // This Method is valid for 1D measurements only
1099 
1100  const unsigned int xIndex = iVirtualMeas + refTrajPtr->numberOfHitMeas();
1101  // Covariance into a TMatrixDSym
1102 
1103  //aHitCovarianceM = new TMatrixDSym(1);
1104  aHitCovarianceM(0,0)=refTrajPtr->measurementErrors()[xIndex][xIndex];
1105 
1106  //theHitResidualsM= new TMatrixF(1,1);
1107  aHitResidualsM(0,0)= refTrajPtr->measurements()[xIndex];
1108 
1109  // Local Derivatives into a TMatrixDSym (to use matrix operations)
1110  const AlgebraicMatrix &locDerivMatrix = refTrajPtr->derivatives();
1111  // theLocalDerivativeNumber = locDerivMatrix.num_col();
1112 
1113  //theLocalDerivativesM = new TMatrixF(1,locDerivMatrix.num_col());
1114  for (int i = 0; i < locDerivMatrix.num_col(); ++i) {
1115  aLocalDerivativesM(0,i) = locDerivMatrix[xIndex][i];
1116  }
1117 }
1118 
1119 //__________________________________________________________________________________________________
1122  unsigned int iTrajHit, TMatrixDSym &aHitCovarianceM,
1123  TMatrixF &aHitResidualsM, TMatrixF &aLocalDerivativesM)
1124 {
1125  // This Method is valid for 2D measurements only
1126 
1127  const unsigned int xIndex = iTrajHit*2;
1128  const unsigned int yIndex = iTrajHit*2+1;
1129  // Covariance into a TMatrixDSym
1130 
1131  //aHitCovarianceM = new TMatrixDSym(2);
1132  aHitCovarianceM(0,0)=refTrajPtr->measurementErrors()[xIndex][xIndex];
1133  aHitCovarianceM(0,1)=refTrajPtr->measurementErrors()[xIndex][yIndex];
1134  aHitCovarianceM(1,0)=refTrajPtr->measurementErrors()[yIndex][xIndex];
1135  aHitCovarianceM(1,1)=refTrajPtr->measurementErrors()[yIndex][yIndex];
1136 
1137  //theHitResidualsM= new TMatrixF(2,1);
1138  aHitResidualsM(0,0)= refTrajPtr->measurements()[xIndex] - refTrajPtr->trajectoryPositions()[xIndex];
1139  aHitResidualsM(1,0)= refTrajPtr->measurements()[yIndex] - refTrajPtr->trajectoryPositions()[yIndex];
1140 
1141  // Local Derivatives into a TMatrixDSym (to use matrix operations)
1142  const AlgebraicMatrix &locDerivMatrix = refTrajPtr->derivatives();
1143  // theLocalDerivativeNumber = locDerivMatrix.num_col();
1144 
1145  //theLocalDerivativesM = new TMatrixF(2,locDerivMatrix.num_col());
1146  for (int i = 0; i < locDerivMatrix.num_col(); ++i) {
1147  aLocalDerivativesM(0,i) = locDerivMatrix[xIndex][i];
1148  aLocalDerivativesM(1,i) = locDerivMatrix[yIndex][i];
1149  }
1150 }
1151 
1152 //__________________________________________________________________________________________________
1155  unsigned int iTrajHit, const std::vector<int> &globalLabels,
1156  const std::vector<float> &globalDerivativesX,
1157  const std::vector<float> &globalDerivativesY)
1158 {
1159  const ConstRecHitPointer aRecHit(refTrajPtr->recHits()[iTrajHit]);
1160 
1161  if((aRecHit)->dimension() == 1) {
1162  return this->callMille1D(refTrajPtr, iTrajHit, globalLabels, globalDerivativesX);
1163  } else {
1164  return this->callMille2D(refTrajPtr, iTrajHit, globalLabels,
1165  globalDerivativesX, globalDerivativesY);
1166  }
1167 }
1168 
1169 
1170 //__________________________________________________________________________________________________
1173  unsigned int iTrajHit, const std::vector<int> &globalLabels,
1174  const std::vector<float> &globalDerivativesX)
1175 {
1176  const ConstRecHitPointer aRecHit(refTrajPtr->recHits()[iTrajHit]);
1177  const unsigned int xIndex = iTrajHit*2; // the even ones are local x
1178 
1179  // local derivatives
1180  const AlgebraicMatrix &locDerivMatrix = refTrajPtr->derivatives();
1181  const int nLocal = locDerivMatrix.num_col();
1182  std::vector<float> localDerivatives(nLocal);
1183  for (unsigned int i = 0; i < localDerivatives.size(); ++i) {
1184  localDerivatives[i] = locDerivMatrix[xIndex][i];
1185  }
1186 
1187  // residuum and error
1188  float residX = refTrajPtr->measurements()[xIndex] - refTrajPtr->trajectoryPositions()[xIndex];
1189  float hitErrX = TMath::Sqrt(refTrajPtr->measurementErrors()[xIndex][xIndex]);
1190 
1191  // number of global derivatives
1192  const int nGlobal = globalDerivativesX.size();
1193 
1194  // &(localDerivatives[0]) etc. are valid - as long as vector is not empty
1195  // cf. http://www.parashift.com/c++-faq-lite/containers.html#faq-34.3
1196  theMille->mille(nLocal, &(localDerivatives[0]), nGlobal, &(globalDerivativesX[0]),
1197  &(globalLabels[0]), residX, hitErrX);
1198 
1199  if (theMonitor) {
1200  theMonitor->fillDerivatives(aRecHit, &(localDerivatives[0]), nLocal,
1201  &(globalDerivativesX[0]), nGlobal, &(globalLabels[0]));
1202  theMonitor->fillResiduals(aRecHit, refTrajPtr->trajectoryStates()[iTrajHit],
1203  iTrajHit, residX, hitErrX, false);
1204  }
1205 
1206  return 1;
1207 }
1208 
1209 //__________________________________________________________________________________________________
1212  unsigned int iTrajHit, const std::vector<int> &globalLabels,
1213  const std::vector<float> &globalDerivativesx,
1214  const std::vector<float> &globalDerivativesy)
1215 {
1216  const ConstRecHitPointer aRecHit(refTrajPtr->recHits()[iTrajHit]);
1217 
1218  if((aRecHit)->dimension() != 2) {
1219  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::callMille2D"
1220  << "You try to call method for 2D hits for a "
1221  << (aRecHit)->dimension()
1222  << "D Hit. Hit gets ignored!";
1223  return -1;
1224  }
1225 
1226  TMatrixDSym aHitCovarianceM(2);
1227  TMatrixF aHitResidualsM(2,1);
1228  TMatrixF aLocalDerivativesM(2, refTrajPtr->derivatives().num_col());
1229  // below method fills above 3 matrices
1230  this->addRefTrackData2D(refTrajPtr, iTrajHit, aHitCovarianceM,aHitResidualsM,aLocalDerivativesM);
1231  TMatrixF aGlobalDerivativesM(2,globalDerivativesx.size());
1232  this->makeGlobDerivMatrix(globalDerivativesx, globalDerivativesy, aGlobalDerivativesM);
1233 
1234  // calculates correlation between Hit measurements
1235  // FIXME: Should take correlation (and resulting transformation) from original hit,
1236  // not 2x2 matrix from ReferenceTrajectory: That can come from error propagation etc.!
1237  const double corr = aHitCovarianceM(0,1) / sqrt(aHitCovarianceM(0,0) * aHitCovarianceM(1,1));
1238  if (theMonitor) theMonitor->fillCorrelations2D(corr, aRecHit);
1239  bool diag = false; // diagonalise only tracker TID, TEC
1240  switch(aRecHit->geographicalId().subdetId()) {
1241  case SiStripDetId::TID:
1242  case SiStripDetId::TEC:
1243  if (aRecHit->geographicalId().det() == DetId::Tracker && TMath::Abs(corr) > theMaximalCor2D) {
1244  this->diagonalize(aHitCovarianceM, aLocalDerivativesM, aHitResidualsM, aGlobalDerivativesM);
1245  diag = true;
1246  }
1247  break;
1248  default:;
1249  }
1250 
1251  float newResidX = aHitResidualsM(0,0);
1252  float newResidY = aHitResidualsM(1,0);
1253  float newHitErrX = TMath::Sqrt(aHitCovarianceM(0,0));
1254  float newHitErrY = TMath::Sqrt(aHitCovarianceM(1,1));
1255  float *newLocalDerivsX = aLocalDerivativesM[0].GetPtr();
1256  float *newLocalDerivsY = aLocalDerivativesM[1].GetPtr();
1257  float *newGlobDerivsX = aGlobalDerivativesM[0].GetPtr();
1258  float *newGlobDerivsY = aGlobalDerivativesM[1].GetPtr();
1259  const int nLocal = aLocalDerivativesM.GetNcols();
1260  const int nGlobal = aGlobalDerivativesM.GetNcols();
1261 
1262  if (diag && (newHitErrX > newHitErrY)) { // also for 2D hits?
1263  // measurement with smaller error is x-measurement (for !is2D do not fill y-measurement):
1264  std::swap(newResidX, newResidY);
1265  std::swap(newHitErrX, newHitErrY);
1266  std::swap(newLocalDerivsX, newLocalDerivsY);
1267  std::swap(newGlobDerivsX, newGlobDerivsY);
1268  }
1269 
1270  // &(globalLabels[0]) is valid - as long as vector is not empty
1271  // cf. http://www.parashift.com/c++-faq-lite/containers.html#faq-34.3
1272  theMille->mille(nLocal, newLocalDerivsX, nGlobal, newGlobDerivsX,
1273  &(globalLabels[0]), newResidX, newHitErrX);
1274 
1275  if (theMonitor) {
1276  theMonitor->fillDerivatives(aRecHit, newLocalDerivsX, nLocal, newGlobDerivsX, nGlobal,
1277  &(globalLabels[0]));
1278  theMonitor->fillResiduals(aRecHit, refTrajPtr->trajectoryStates()[iTrajHit],
1279  iTrajHit, newResidX, newHitErrX, false);
1280  }
1281  const bool isReal2DHit = this->is2D(aRecHit); // strip is 1D (except matched hits)
1282  if (isReal2DHit) {
1283  theMille->mille(nLocal, newLocalDerivsY, nGlobal, newGlobDerivsY,
1284  &(globalLabels[0]), newResidY, newHitErrY);
1285  if (theMonitor) {
1286  theMonitor->fillDerivatives(aRecHit, newLocalDerivsY, nLocal, newGlobDerivsY, nGlobal,
1287  &(globalLabels[0]));
1288  theMonitor->fillResiduals(aRecHit, refTrajPtr->trajectoryStates()[iTrajHit],
1289  iTrajHit, newResidY, newHitErrY, true);// true: y
1290  }
1291  }
1292 
1293  return (isReal2DHit ? 2 : 1);
1294 }
1295 
1296 //__________________________________________________________________________________________________
1298 ::addVirtualMeas(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iVirtualMeas)
1299 {
1300  TMatrixDSym aHitCovarianceM(1);
1301  TMatrixF aHitResidualsM(1,1);
1302  TMatrixF aLocalDerivativesM(1, refTrajPtr->derivatives().num_col());
1303  // below method fills above 3 'matrices'
1304  this->addRefTrackVirtualMeas1D(refTrajPtr, iVirtualMeas, aHitCovarianceM, aHitResidualsM, aLocalDerivativesM);
1305 
1306  // no global parameters (use dummy 0)
1307  TMatrixF aGlobalDerivativesM(1,1);
1308  aGlobalDerivativesM(0,0) = 0;
1309 
1310  float newResidX = aHitResidualsM(0,0);
1311  float newHitErrX = TMath::Sqrt(aHitCovarianceM(0,0));
1312  float *newLocalDerivsX = aLocalDerivativesM[0].GetPtr();
1313  float *newGlobDerivsX = aGlobalDerivativesM[0].GetPtr();
1314  const int nLocal = aLocalDerivativesM.GetNcols();
1315  const int nGlobal = 0;
1316 
1317  theMille->mille(nLocal, newLocalDerivsX, nGlobal, newGlobDerivsX,
1318  &nGlobal, newResidX, newHitErrX);
1319 }
1320 
1321 //____________________________________________________
1323  const TkFittedLasBeamCollection &lasBeams,
1324  const TsosVectorCollection &lasBeamTsoses)
1325 {
1326  TsosVectorCollection::const_iterator iTsoses = lasBeamTsoses.begin();
1327  for(TkFittedLasBeamCollection::const_iterator iBeam = lasBeams.begin(), iEnd = lasBeams.end();
1328  iBeam != iEnd; ++iBeam, ++iTsoses){ // beam/tsoses parallel!
1329 
1330  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::addLaserData"
1331  << "Beam " << iBeam->getBeamId() << " with "
1332  << iBeam->parameters().size() << " parameters and "
1333  << iBeam->getData().size() << " hits.\n There are "
1334  << iTsoses->size() << " TSOSes.";
1335 
1336  this->addLasBeam(eventInfo, *iBeam, *iTsoses);
1337  }
1338 }
1339 
1340 //____________________________________________________
1342  const TkFittedLasBeam &lasBeam,
1343  const std::vector<TrajectoryStateOnSurface> &tsoses)
1344 {
1345  AlignmentParameters *dummyPtr = 0; // for globalDerivativesHierarchy()
1346  std::vector<float> lasLocalDerivsX; // buffer for local derivatives
1347  const unsigned int beamLabel = thePedeLabels->lasBeamLabel(lasBeam.getBeamId());// for global par
1348 
1349  for (unsigned int iHit = 0; iHit < tsoses.size(); ++iHit) {
1350  if (!tsoses[iHit].isValid()) continue;
1351  // clear buffer
1352  theFloatBufferX.clear();
1353  theFloatBufferY.clear();
1354  theIntBuffer.clear();
1355  lasLocalDerivsX.clear();
1356  // get alignables and global parameters
1357  const SiStripLaserRecHit2D &hit = lasBeam.getData()[iHit];
1358  AlignableDetOrUnitPtr lasAli(theAlignableNavigator->alignableFromDetId(hit.getDetId()));
1359  this->globalDerivativesHierarchy(eventInfo,
1360  tsoses[iHit], lasAli, lasAli,
1362  // fill derivatives vector from derivatives matrix
1363  for (unsigned int nFitParams = 0;
1364  nFitParams < static_cast<unsigned int>(lasBeam.parameters().size());
1365  ++nFitParams) {
1366  const float derivative = lasBeam.derivatives()[iHit][nFitParams];
1367  if (nFitParams < lasBeam.firstFixedParameter()) { // first local beam parameters
1368  lasLocalDerivsX.push_back(derivative);
1369  } else { // now global ones
1370  const unsigned int numPar = nFitParams - lasBeam.firstFixedParameter();
1371  theIntBuffer.push_back(thePedeLabels->parameterLabel(beamLabel, numPar));
1372  theFloatBufferX.push_back(derivative);
1373  }
1374  } // end loop over parameters
1375 
1376  const float residual = hit.localPosition().x() - tsoses[iHit].localPosition().x();
1377  // error from file or assume 0.003
1378  const float error = 0.003; // hit.localPositionError().xx(); sqrt???
1379 
1380  theMille->mille(lasLocalDerivsX.size(), &(lasLocalDerivsX[0]), theFloatBufferX.size(),
1381  &(theFloatBufferX[0]), &(theIntBuffer[0]), residual, error);
1382  } // end of loop over hits
1383 
1384  theMille->end();
1385 }
1386 
1388 {
1389  // do some printing, if requested
1390  const bool doOutputOnStdout(pxbSurveyCfg.getParameter<bool>("doOutputOnStdout"));
1391  if (doOutputOnStdout) std::cout << "# Output from addPxbSurvey follows below because doOutputOnStdout is set to True" << std::endl;
1392 
1393  // instantiate a dicer object
1394  SurveyPxbDicer dicer(pxbSurveyCfg.getParameter<std::vector<edm::ParameterSet> >("toySurveyParameters"), pxbSurveyCfg.getParameter<unsigned int>("toySurveySeed"));
1395  std::ofstream outfile(pxbSurveyCfg.getUntrackedParameter<std::string>("toySurveyFile").c_str());
1396 
1397  // read data from file
1398  std::vector<SurveyPxbImageLocalFit> measurements;
1399  std::string filename(pxbSurveyCfg.getParameter<edm::FileInPath>("infile").fullPath());
1401 
1402  // loop over photographs (=measurements) and perform the fit
1403  for(std::vector<SurveyPxbImageLocalFit>::size_type i=0; i!=measurements.size(); i++)
1404  {
1405  if (doOutputOnStdout) std::cout << "Module " << i << ": ";
1406 
1407  // get the Alignables and their surfaces
1408  AlignableDetOrUnitPtr mod1(theAlignableNavigator->alignableFromDetId(measurements[i].getIdFirst()));
1409  AlignableDetOrUnitPtr mod2(theAlignableNavigator->alignableFromDetId(measurements[i].getIdSecond()));
1410  const AlignableSurface& surf1 = mod1->surface();
1411  const AlignableSurface& surf2 = mod2->surface();
1412 
1413  // the position of the fiducial points in local frame of a PXB module
1414  const LocalPoint fidpoint0(-0.91,+3.30);
1415  const LocalPoint fidpoint1(+0.91,+3.30);
1416  const LocalPoint fidpoint2(+0.91,-3.30);
1417  const LocalPoint fidpoint3(-0.91,-3.30);
1418 
1419  // We choose the local frame of the first module as reference,
1420  // so take the fidpoints of the second module and calculate their
1421  // positions in the reference frame
1422  const GlobalPoint surf2point0(surf2.toGlobal(fidpoint0));
1423  const GlobalPoint surf2point1(surf2.toGlobal(fidpoint1));
1424  const LocalPoint fidpoint0inSurf1frame(surf1.toLocal(surf2point0));
1425  const LocalPoint fidpoint1inSurf1frame(surf1.toLocal(surf2point1));
1426 
1427  // Create the vector for the fit
1429  fidpointvec.push_back(fidpoint0inSurf1frame);
1430  fidpointvec.push_back(fidpoint1inSurf1frame);
1431  fidpointvec.push_back(fidpoint2);
1432  fidpointvec.push_back(fidpoint3);
1433 
1434  // if toy survey is requested, dice the values now
1435  if (pxbSurveyCfg.getParameter<bool>("doToySurvey"))
1436  {
1437  dicer.doDice(fidpointvec,measurements[i].getIdPair(), outfile);
1438  }
1439 
1440  // do the fit
1441  measurements[i].doFit(fidpointvec, thePedeLabels->alignableLabel(mod1), thePedeLabels->alignableLabel(mod2));
1442  SurveyPxbImageLocalFit::localpars_t a; // local pars from fit
1443  a = measurements[i].getLocalParameters();
1444  const SurveyPxbImageLocalFit::value_t chi2 = measurements[i].getChi2();
1445 
1446  // do some reporting, if requested
1447  if (doOutputOnStdout)
1448  {
1449  std::cout << "a: " << a[0] << ", " << a[1] << ", " << a[2] << ", " << a[3]
1450  << " S= " << sqrt(a[2]*a[2]+a[3]*a[3])
1451  << " phi= " << atan(a[3]/a[2])
1452  << " chi2= " << chi2 << std::endl;
1453  }
1454  if (theMonitor)
1455  {
1456  theMonitor->fillPxbSurveyHistsChi2(chi2);
1457  theMonitor->fillPxbSurveyHistsLocalPars(a[0],a[1],sqrt(a[2]*a[2]+a[3]*a[3]),atan(a[3]/a[2]));
1458  }
1459 
1460  // pass the results from the local fit to mille
1462  {
1463  theMille->mille((int)measurements[i].getLocalDerivsSize(),
1464  measurements[i].getLocalDerivsPtr(j),
1465  (int)measurements[i].getGlobalDerivsSize(),
1466  measurements[i].getGlobalDerivsPtr(j),
1467  measurements[i].getGlobalDerivsLabelPtr(j),
1468  measurements[i].getResiduum(j),
1469  measurements[i].getSigma(j));
1470  }
1471  theMille->end();
1472  }
1473  outfile.close();
1474 }
1475 
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
tuple cfg
Definition: looper.py:293
void increaseHitsX(unsigned int add=1)
increase hits for x-measurement
bool read(std::vector< Alignable * > &alignables, bool setUserVars)
Definition: PedeReader.cc:52
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
void increaseHitsY(unsigned int add=1)
increase hits for y-measurement
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
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
void restoreCachedTransformations(void)
restore the previously cached position, rotation and other parameters
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
int loop
CMSSW
virtual LocalPoint localPosition() const
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)
std::pair< RunNumber, RunNumber > RunRange
#define nullptr
virtual void initialize(const edm::EventSetup &setup, AlignableTracker *tracker, AlignableMuon *muon, AlignableExtras *extras, AlignmentParameterStore *store)
Called at beginning of job.
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:57
const ConstTrajTrackPairCollection & trajTrackPairs() const
define event information passed to algorithms
tuple result
Definition: mps_fire.py:83
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
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::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
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.
CLHEP::HepMatrix AlgebraicMatrix
const std::vector< Scalar > & parameters() const
parallel to derivatives()
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
virtual bool setParametersForRunRange(const RunRange &runrange)
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
T sqrt(T t)
Definition: SSEVec.h:18
bool is2D(HitType hitType)
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
std::unique_ptr< PedeLabelerBase > thePedeLabels
const TkFittedLasBeamCollection * tkLasBeams() const
std::unique_ptr< PedeSteerer > thePedeSteer
void addGlobals(const std::vector< int > &aLabels, const TMatrixD &aDerivatives)
Add global derivatives to a point.
Definition: GblPoint.cc:298
const AlignableSurface & surface() const
Return the Surface (global position and orientation) of the object.
Definition: Alignable.h:131
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 ...
CLHEP::HepVector AlgebraicVector
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)
tuple tracks
Definition: testEve_cfg.py:39
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
T const * product() const
Definition: ESHandle.h:86
std::vector< value_t > localpars_t
virtual bool addCalibrations(const std::vector< IntegratedCalibrationBase * > &iCals)
Pass integrated calibrations to Millepede (they are not owned by Millepede!)
virtual void run(const edm::EventSetup &setup, const EventInfo &eventInfo)
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.
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
tuple filename
Definition: lut2db_cfg.py:20
TransientTrackingRecHit::ConstRecHitContainer ConstRecHitContainer
tuple cout
Definition: gather_cfg.py:145
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
void addPxbSurvey(const edm::ParameterSet &pxbSurveyCfg)
add measurement data from PXB survey
unsigned int decodeMode(const std::string &mode) const
std::string fullPath() const
Definition: FileInPath.cc:184
std::unique_ptr< AlignableNavigator > theAlignableNavigator
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:36
uint32_t dimension(pat::CandKinResolution::Parametrization parametrization)
Returns the number of free parameters in a parametrization (3 or 4)
Alignable * mother() const
Return pointer to container alignable (if any)
Definition: Alignable.h:90
tuple size
Write out results.
T get(const Candidate &c)
Definition: component.h:55
std::vector< ConstTrajTrackPair > ConstTrajTrackPairCollection
bool setAllDefault(unsigned int nParam)
set default values for all data concerning nParam (false if nParam out of range)
Point on trajectory.
Definition: GblPoint.h:46
virtual void terminate()
Called at end of job.