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