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