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