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