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 
15 
17 // in header, too
18 // end in header, too
19 
23 #include "Mille.h" // 'unpublished' interface located in src
24 #include "PedeSteerer.h" // dito
25 #include "PedeReader.h" // dito
28 
31 
34 
38 
39 // includes to make known that they inherit from Alignable:
43 
49 
52 
54 
55 #include <fstream>
56 #include <sstream>
57 #include <algorithm>
58 
59 #include <TMath.h>
60 #include <TMatrixDSymEigen.h>
64 
65 // Includes for PXB survey
66 #include <iostream>
73 
75 
76 // Constructor ----------------------------------------------------------------
77 //____________________________________________________
80  theConfig(cfg), theMode(this->decodeMode(theConfig.getUntrackedParameter<std::string>("mode"))),
81  theDir(theConfig.getUntrackedParameter<std::string>("fileDir")),
82  theAlignmentParameterStore(0), theAlignables(), theAlignableNavigator(0),
83  theMonitor(0), theMille(0), thePedeLabels(0), thePedeSteer(0),
84  theTrajectoryFactory(0),
85  theMinNumHits(cfg.getParameter<unsigned int>("minNumHits")),
86  theMaximalCor2D(cfg.getParameter<double>("max2Dcorrelation")),
87  theLastWrittenIov(0)
88 {
89  if (!theDir.empty() && theDir.find_last_of('/') != theDir.size()-1) theDir += '/';// may need '/'
90  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm" << "Start in mode '"
91  << theConfig.getUntrackedParameter<std::string>("mode")
92  << "' with output directory '" << theDir << "'.";
93  if (this->isMode(myMilleBit)) {
94  theMille = new Mille((theDir + theConfig.getParameter<std::string>("binaryFile")).c_str());// add ', false);' for text output);
95  }
96 }
97 
98 // Destructor ----------------------------------------------------------------
99 //____________________________________________________
101 {
102  delete theAlignableNavigator;
104  delete theMille;
105  theMille = 0;
106  delete theMonitor;
107  theMonitor = 0;
108  delete thePedeSteer;
109  thePedeSteer = 0;
110  delete thePedeLabels;
111  thePedeLabels = 0;
112  delete theTrajectoryFactory;
114 }
115 
116 // Call at beginning of job ---------------------------------------------------
117 //____________________________________________________
121 {
122  if (muon) {
123  edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::initialize"
124  << "Running with AlignabeMuon not yet tested.";
125  }
126 
127  theAlignableNavigator = new AlignableNavigator(extras, tracker, muon);
130 
131  edm::ParameterSet pedeLabelerCfg(theConfig.getParameter<edm::ParameterSet>("pedeLabeler"));
132  edm::VParameterSet RunRangeSelectionVPSet(theConfig.getUntrackedParameter<edm::VParameterSet>("RunRangeSelection"));
133  pedeLabelerCfg.addUntrackedParameter<edm::VParameterSet>("RunRangeSelection",
134  RunRangeSelectionVPSet);
135 
136  std::string labelerPlugin = "PedeLabeler";
137  if (RunRangeSelectionVPSet.size()>0) {
138  labelerPlugin = "RunRangeDependentPedeLabeler";
139  if (pedeLabelerCfg.exists("plugin")) {
140  std::string labelerPluginCfg = pedeLabelerCfg.getParameter<std::string>("plugin");
141  if ((labelerPluginCfg!="PedeLabeler" && labelerPluginCfg!="RunRangeDependentPedeLabeler") ||
142  pedeLabelerCfg.getUntrackedParameter<edm::VParameterSet>("parameterInstances").size()>0) {
143  throw cms::Exception("BadConfig")
144  << "MillePedeAlignmentAlgorithm::initialize"
145  << "both RunRangeSelection and generic labeler specified in config file. "
146  << "Please get rid of either one of them.\n";
147  }
148  }
149  } else {
150  if (pedeLabelerCfg.exists("plugin")) {
151  labelerPlugin = pedeLabelerCfg.getParameter<std::string>("plugin");
152  }
153  }
154 
155  if (!pedeLabelerCfg.exists("plugin")) {
156  pedeLabelerCfg.addUntrackedParameter<std::string>("plugin", labelerPlugin);
157  }
158 
159  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::initialize"
160  << "Using plugin '" << labelerPlugin << "' to generate labels.";
161 
162  thePedeLabels = PedeLabelerPluginFactory::get()->create(labelerPlugin,
163  PedeLabelerBase::TopLevelAlignables(tracker, muon, extras),
164  pedeLabelerCfg);
165 
166  // 1) Create PedeSteerer: correct alignable positions for coordinate system selection
167  edm::ParameterSet pedeSteerCfg(theConfig.getParameter<edm::ParameterSet>("pedeSteerer"));
168  thePedeSteer = new PedeSteerer(tracker, muon, extras,
170  pedeSteerCfg, theDir, !this->isMode(myPedeSteerBit));
171 
172  // 2) If requested, directly read in and apply result of previous pede run,
173  // assuming that correction from 1) was also applied to create the result:
174  const std::vector<edm::ParameterSet> mprespset
175  (theConfig.getParameter<std::vector<edm::ParameterSet> >("pedeReaderInputs"));
176  if (!mprespset.empty()) {
177  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::initialize"
178  << "Apply " << mprespset.end() - mprespset.begin()
179  << " previous MillePede constants from 'pedeReaderInputs'.";
180  }
181 
182  RunRange runrange(cond::timeTypeSpecs[cond::runnumber].beginValue,
184  for (std::vector<edm::ParameterSet>::const_iterator iSet = mprespset.begin(), iE = mprespset.end();
185  iSet != iE; ++iSet) {
186  if (!this->readFromPede((*iSet), false, runrange)) { // false: do not erase SelectionUserVariables
187  throw cms::Exception("BadConfig")
188  << "MillePedeAlignmentAlgorithm::initialize: Problems reading input constants of "
189  << "pedeReaderInputs entry " << iSet - mprespset.begin() << '.';
190  }
192  // Needed to shut up later warning from checkAliParams:
194  }
195 
196  // 3) Now create steerings with 'final' start position:
197  thePedeSteer->buildSubSteer(tracker, muon, extras);
198 
199  // After (!) 1-3 of PedeSteerer which uses the SelectionUserVariables attached to the parameters:
200  this->buildUserVariables(theAlignables); // for hit statistics and/or pede result
201 
202  if (this->isMode(myMilleBit)) {
203  if (!theConfig.getParameter<std::vector<std::string> >("mergeBinaryFiles").empty() ||
204  !theConfig.getParameter<std::vector<std::string> >("mergeTreeFiles").empty()) {
205  throw cms::Exception("BadConfig")
206  << "'vstring mergeTreeFiles' and 'vstring mergeBinaryFiles' must be empty for "
207  << "modes running mille.";
208  }
209  //theMille = new Mille((theDir + theConfig.getParameter<std::string>("binaryFile")).c_str());// add ', false);' for text output);
210  const std::string moniFile(theConfig.getUntrackedParameter<std::string>("monitorFile"));
211  if (moniFile.size()) theMonitor = new MillePedeMonitor((theDir + moniFile).c_str());
212 
213  // Get trajectory factory. In case nothing found, FrameWork will throw...
214  const edm::ParameterSet fctCfg(theConfig.getParameter<edm::ParameterSet>("TrajectoryFactory"));
215  const std::string fctName(fctCfg.getParameter<std::string>("TrajectoryFactoryName"));
216  theTrajectoryFactory = TrajectoryFactoryPlugin::get()->create(fctName, fctCfg);
217  }
218 
219  if (this->isMode(myPedeSteerBit)) {
220  // Get config for survey and set flag accordingly
221  const edm::ParameterSet pxbSurveyCfg(theConfig.getParameter<edm::ParameterSet>("surveyPixelBarrel"));
222  theDoSurveyPixelBarrel = pxbSurveyCfg.getParameter<bool>("doSurvey");
223  if (theDoSurveyPixelBarrel) this->addPxbSurvey(pxbSurveyCfg);
224  }
225 }
226 
227 //____________________________________________________
229 {
230  if (this->isMode(myPedeReadBit)) {
231  // restore initial positions, rotations and deformations
233 
234  // Needed to shut up later warning from checkAliParams:
236 
237  if (!this->readFromPede(theConfig.getParameter<edm::ParameterSet>("pedeReader"), true, runrange)) {
238  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::setParametersForRunRange"
239  << "Problems reading pede result, but applying!";
240  }
242 
243  this->doIO(++theLastWrittenIov); // pre-increment!
244  }
245 
246  return true;
247 }
248 
249 // Call at end of job ---------------------------------------------------------
250 //____________________________________________________
252 {
253  delete theMille;// delete to close binary before running pede below (flush would be enough...)
254  theMille = 0;
255 
256  std::vector<std::string> files;
257  if (this->isMode(myMilleBit) || !theConfig.getParameter<std::string>("binaryFile").empty()) {
258  files.push_back(theDir + theConfig.getParameter<std::string>("binaryFile"));
259  } else {
260  const std::vector<std::string> plainFiles
261  (theConfig.getParameter<std::vector<std::string> >("mergeBinaryFiles"));
262  for (std::vector<std::string>::const_iterator i = plainFiles.begin(), iEnd = plainFiles.end();
263  i != iEnd; ++i) {
264  files.push_back(theDir + *i);
265  }
266  }
267 
268  // cache all positions, rotations and deformations
270 
271  const std::string masterSteer(thePedeSteer->buildMasterSteer(files));// do only if myPedeSteerBit?
272  if (this->isMode(myPedeRunBit)) {
273  thePedeSteer->runPede(masterSteer);
274  }
275 
276  // parameters from pede are not yet applied,
277  // so we can still write start positions (but with hit statistics in case of mille):
278  this->doIO(0);
279  theLastWrittenIov = 0;
280 }
281 
282 // Run the algorithm on trajectories and tracks -------------------------------
283 //____________________________________________________
285 {
286  if (!this->isMode(myMilleBit)) return; // no theMille created...
288 
289  if (theMonitor) { // monitor input tracks
290  for (ConstTrajTrackPairCollection::const_iterator iTrajTrack = tracks.begin();
291  iTrajTrack != tracks.end(); ++iTrajTrack) {
292  theMonitor->fillTrack((*iTrajTrack).second);
293  }
294  }
295 
296  const RefTrajColl trajectories(theTrajectoryFactory->trajectories(setup, tracks, eventInfo.beamSpot_));
297 
298  // Now loop over ReferenceTrajectoryCollection
299  unsigned int refTrajCount = 0; // counter for track monitoring if 1 track per trajectory
300  for (RefTrajColl::const_iterator iRefTraj = trajectories.begin(), iRefTrajE = trajectories.end();
301  iRefTraj != iRefTrajE; ++iRefTraj, ++refTrajCount) {
302 
303  RefTrajColl::value_type refTrajPtr = *iRefTraj;
304  if (theMonitor) theMonitor->fillRefTrajectory(refTrajPtr);
305 
306  const std::pair<unsigned int, unsigned int> nHitXy = this->addReferenceTrajectory(eventInfo, refTrajPtr);
307 
308  if (theMonitor && (nHitXy.first || nHitXy.second)) {
309  // if track used (i.e. some hits), fill monitoring
310  // track NULL ptr if trajectories and tracks do not match
311  const reco::Track *trackPtr =
312  (trajectories.size() == tracks.size() ? tracks[refTrajCount].second : 0);
313  theMonitor->fillUsedTrack(trackPtr, nHitXy.first, nHitXy.second);
314  }
315 
316  } // end of reference trajectory and track loop
317 }
318 
319 
320 //____________________________________________________
321 std::pair<unsigned int, unsigned int>
323  const RefTrajColl::value_type &refTrajPtr)
324 {
325  std::pair<unsigned int, unsigned int> hitResultXy(0,0);
326  if (refTrajPtr->isValid()) {
327 
328  // to add hits if all fine:
329  std::vector<AlignmentParameters*> parVec(refTrajPtr->recHits().size());
330  // collect hit statistics, assuming that there are no y-only hits
331  std::vector<bool> validHitVecY(refTrajPtr->recHits().size(), false);
332  // Use recHits from ReferenceTrajectory (since they have the right order!):
333  for (unsigned int iHit = 0; iHit < refTrajPtr->recHits().size(); ++iHit) {
334  const int flagXY = this->addMeasurementData(eventInfo, refTrajPtr, iHit, parVec[iHit]);
335 
336  if (flagXY < 0) { // problem
337  hitResultXy.first = 0;
338  break;
339  } else { // hit is fine, increase x/y statistics
340  if (flagXY >= 1) ++hitResultXy.first;
341  validHitVecY[iHit] = (flagXY >= 2);
342  }
343  } // end loop on hits
344 
345  // add virtual measurements
346  for (unsigned int iVirtualMeas = 0; iVirtualMeas < refTrajPtr->numberOfVirtualMeas(); ++iVirtualMeas) {
347  this->addVirtualMeas(refTrajPtr, iVirtualMeas);
348  }
349 
350  // kill or end 'track' for mille, depends on #hits criterion
351  if (hitResultXy.first == 0 || hitResultXy.first < theMinNumHits) {
352  theMille->kill();
353  hitResultXy.first = hitResultXy.second = 0; //reset
354  } else {
355  theMille->end();
356  // add x/y hit count to MillePedeVariables of parVec,
357  // returning number of y-hits of the reference trajectory
358  hitResultXy.second = this->addHitCount(parVec, validHitVecY);
359  }
360  } // end if valid trajectory
361 
362  return hitResultXy;
363 }
364 
365 //____________________________________________________
366 unsigned int
367 MillePedeAlignmentAlgorithm::addHitCount(const std::vector<AlignmentParameters*> &parVec,
368  const std::vector<bool> &validHitVecY) const
369 {
370  // Loop on all hit information in the input arrays and count valid y-hits:
371  unsigned int nHitY = 0;
372  for (unsigned int iHit = 0; iHit < validHitVecY.size(); ++iHit) {
373  Alignable *ali = (parVec[iHit] ? parVec[iHit]->alignable() : 0);
374  // Loop upwards on hierarchy of alignables to add hits to all levels
375  // that are currently aligned. If only a non-selected alignable was hit,
376  // (i.e. flagXY == 0 in addReferenceTrajectory(..)), there is no loop at all...
377  while (ali) {
379  if (pars) { // otherwise hierarchy level not selected
380  // cast ensured by previous checks:
381  MillePedeVariables *mpVar = static_cast<MillePedeVariables*>(pars->userVariables());
382  // every hit has an x-measurement, cf. addReferenceTrajectory(..):
383  mpVar->increaseHitsX();
384  if (validHitVecY[iHit]) {
385  mpVar->increaseHitsY();
386  if (pars == parVec[iHit]) ++nHitY; // do not count hits twice
387  }
388  }
389  ali = ali->mother();
390  }
391  }
392 
393  return nHitY;
394 }
395 
396 //____________________________________________________
398  const edm::EventSetup &setup)
399 {
400  if(runInfo.tkLasBeams_ && runInfo.tkLasBeamTsoses_){
401  // LAS beam treatment
402  this->addLaserData(eventInfo, *(runInfo.tkLasBeams_), *(runInfo.tkLasBeamTsoses_));
403  }
404 }
405 
406 //____________________________________________________
409  unsigned int iHit,
410  AlignmentParameters *&params)
411 {
412  params = 0;
413  theFloatBufferX.clear();
414  theFloatBufferY.clear();
415  theIntBuffer.clear();
416 
417  const TrajectoryStateOnSurface &tsos = refTrajPtr->trajectoryStates()[iHit];
418  const ConstRecHitPointer &recHitPtr = refTrajPtr->recHits()[iHit];
419  // ignore invalid hits
420  if (!recHitPtr->isValid()) return 0;
421 
422  // get AlignableDet/Unit for this hit
423  AlignableDetOrUnitPtr alidet(theAlignableNavigator->alignableFromDetId(recHitPtr->geographicalId()));
424 
425  if (!this->globalDerivativesHierarchy(eventInfo,
426  tsos, alidet, alidet, theFloatBufferX, // 2x alidet, sic!
427  theFloatBufferY, theIntBuffer, params)) {
428  return -1; // problem
429  } else if (theFloatBufferX.empty()) {
430  return 0; // empty for X: no alignable for hit
431  } else {
432  return this->callMille(refTrajPtr, iHit, theIntBuffer, theFloatBufferX, theFloatBufferY);
433  }
434 }
435 
436 //____________________________________________________
439  const TrajectoryStateOnSurface &tsos,
440  Alignable *ali, const AlignableDetOrUnitPtr &alidet,
441  std::vector<float> &globalDerivativesX,
442  std::vector<float> &globalDerivativesY,
443  std::vector<int> &globalLabels,
444  AlignmentParameters *&lowestParams) const
445 {
446  // derivatives and labels are recursively attached
447  if (!ali) return true; // no mother might be OK
448 
449  if (false && theMonitor && alidet != ali) theMonitor->fillFrameToFrame(alidet, ali);
450 
451  AlignmentParameters *params = ali->alignmentParameters();
452 
453  if (params) {
454  if (!lowestParams) lowestParams = params; // set parameters of lowest level
455 
456  bool hasSplitParameters = thePedeLabels->hasSplitParameters(ali);
457  const unsigned int alignableLabel = thePedeLabels->alignableLabel(ali);
458 
459  if (0 == alignableLabel) { // FIXME: what about regardAllHits in Markus' code?
460  edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::globalDerivativesHierarchy"
461  << "Label not found, skip Alignable.";
462  return false;
463  }
464 
465  const std::vector<bool> &selPars = params->selector();
466  const AlgebraicMatrix derivs(params->derivatives(tsos, alidet));
467 
468  // cols: 2, i.e. x&y, rows: parameters, usually RigidBodyAlignmentParameters::N_PARAM
469  for (unsigned int iSel = 0; iSel < selPars.size(); ++iSel) {
470  if (selPars[iSel]) {
471  globalDerivativesX.push_back(derivs[iSel][kLocalX]
472  /thePedeSteer->cmsToPedeFactor(iSel));
473  if (hasSplitParameters==true) {
474  globalLabels.push_back(thePedeLabels->parameterLabel(ali, iSel, eventInfo, tsos));
475  } else {
476  globalLabels.push_back(thePedeLabels->parameterLabel(alignableLabel, iSel));
477  }
478  globalDerivativesY.push_back(derivs[iSel][kLocalY]
479  /thePedeSteer->cmsToPedeFactor(iSel));
480  }
481  }
482  // Exclude mothers if Alignable selected to be no part of a hierarchy:
483  if (thePedeSteer->isNoHiera(ali)) return true;
484  }
485  // Call recursively for mother, will stop if mother == 0:
486  return this->globalDerivativesHierarchy(eventInfo,
487  tsos, ali->mother(), alidet,
488  globalDerivativesX, globalDerivativesY,
489  globalLabels, lowestParams);
490 }
491 
492 // //____________________________________________________
493 // void MillePedeAlignmentAlgorithm
494 // ::callMille(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
495 // unsigned int iTrajHit, MeasurementDirection xOrY,
496 // const std::vector<float> &globalDerivatives, const std::vector<int> &globalLabels)
497 // {
498 // const unsigned int xyIndex = iTrajHit*2 + xOrY;
499 // // FIXME: here for residuum and sigma we could use KALMAN-Filter results
500 // const float residuum =
501 // refTrajPtr->measurements()[xyIndex] - refTrajPtr->trajectoryPositions()[xyIndex];
502 // const float covariance = refTrajPtr->measurementErrors()[xyIndex][xyIndex];
503 // const float sigma = (covariance > 0. ? TMath::Sqrt(covariance) : 0.);
504 
505 // const AlgebraicMatrix &locDerivMatrix = refTrajPtr->derivatives();
506 
507 // std::vector<float> localDerivs(locDerivMatrix.num_col());
508 // for (unsigned int i = 0; i < localDerivs.size(); ++i) {
509 // localDerivs[i] = locDerivMatrix[xyIndex][i];
510 // }
511 
512 // // &(vector[0]) is valid - as long as vector is not empty
513 // // cf. http://www.parashift.com/c++-faq-lite/containers.html#faq-34.3
514 // theMille->mille(localDerivs.size(), &(localDerivs[0]),
515 // globalDerivatives.size(), &(globalDerivatives[0]), &(globalLabels[0]),
516 // residuum, sigma);
517 // if (theMonitor) {
518 // theMonitor->fillDerivatives(refTrajPtr->recHits()[iTrajHit],localDerivs, globalDerivatives,
519 // (xOrY == kLocalY));
520 // theMonitor->fillResiduals(refTrajPtr->recHits()[iTrajHit],
521 // refTrajPtr->trajectoryStates()[iTrajHit],
522 // iTrajHit, residuum, sigma, (xOrY == kLocalY));
523 // }
524 // }
525 
526 //____________________________________________________
528 {
529  // FIXME: Check whether this is a reliable and recommended way to find out...
530 
531  if (recHit->dimension() < 2) {
532  return false; // some muon and TIB/TOB stuff really has RecHit1D
533  } else if (recHit->detUnit()) { // detunit in strip is 1D, in pixel 2D
534  return recHit->detUnit()->type().isTrackerPixel();
535  } else { // stereo strips (FIXME: endcap trouble due to non-parallel strips (wedge sensors)?)
536  if (dynamic_cast<const ProjectedSiStripRecHit2D*>(recHit->hit())) { // check persistent hit
537  // projected: 1D measurement on 'glued' module
538  return false;
539  } else {
540  return true;
541  }
542  }
543 }
544 
545 //__________________________________________________________________________________________________
546 bool MillePedeAlignmentAlgorithm::readFromPede(const edm::ParameterSet &mprespset, bool setUserVars,
547  const RunRange &runrange)
548 {
549  bool allEmpty = this->areEmptyParams(theAlignables);
550 
551  PedeReader reader(mprespset, *thePedeSteer, *thePedeLabels, runrange);
552  std::vector<Alignable*> alis;
553  bool okRead = reader.read(alis, setUserVars);
554  bool numMatch = true;
555 
556  std::stringstream out("Read ");
557  out << alis.size() << " alignables";
558  if (alis.size() != theAlignables.size()) {
559  out << " while " << theAlignables.size() << " in store";
560  numMatch = false; // FIXME: Should we check one by one? Or transfer 'alis' to the store?
561  }
562  if (!okRead) out << ", but problems in reading";
563  if (!allEmpty) out << ", possibly overwriting previous settings";
564  out << ".";
565 
566  if (okRead && allEmpty) {
567  if (numMatch) { // as many alignables with result as trying to align
568  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
569  } else if (alis.size()) { // dead module do not get hits and no pede result
570  edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
571  } else { // serious problem: no result read - and not all modules can be dead...
572  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
573  return false;
574  }
575  return true;
576  }
577  // the rest is not OK:
578  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
579  return false;
580 }
581 
582 //__________________________________________________________________________________________________
583 bool MillePedeAlignmentAlgorithm::areEmptyParams(const std::vector<Alignable*> &alignables) const
584 {
585 
586  for (std::vector<Alignable*>::const_iterator iAli = alignables.begin();
587  iAli != alignables.end(); ++iAli) {
588  const AlignmentParameters *params = (*iAli)->alignmentParameters();
589  if (params) {
590  const AlgebraicVector &parVec(params->parameters());
591  const AlgebraicMatrix &parCov(params->covariance());
592  for (int i = 0; i < parVec.num_row(); ++i) {
593  if (parVec[i] != 0.) return false;
594  for (int j = i; j < parCov.num_col(); ++j) {
595  if (parCov[i][j] != 0.) return false;
596  }
597  }
598  }
599  }
600 
601  return true;
602 }
603 
604 //__________________________________________________________________________________________________
606 {
607  unsigned int result = 0;
608 
609  const std::string outFilePlain(theConfig.getParameter<std::string>("treeFile"));
610  if (outFilePlain.empty()) {
611  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
612  << "treeFile parameter empty => skip writing for 'loop' " << loop;
613  return result;
614  }
615 
616  const std::string outFile(theDir + outFilePlain);
617 
618  AlignmentIORoot aliIO;
619  int ioerr = 0;
620  if (loop == 0) {
621  aliIO.writeAlignableOriginalPositions(theAlignables, outFile.c_str(), loop, false, ioerr);
622  if (ioerr) {
623  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
624  << "Problem " << ioerr << " in writeAlignableOriginalPositions";
625  ++result;
626  }
627  } else if (loop == 1) {
628  // only for first iov add hit counts, else 2x, 3x,... number of hits in IOV 2, 3,...
629  const std::vector<std::string> inFiles
630  (theConfig.getParameter<std::vector<std::string> >("mergeTreeFiles"));
631  const std::vector<std::string> binFiles
632  (theConfig.getParameter<std::vector<std::string> >("mergeBinaryFiles"));
633  if (inFiles.size() != binFiles.size()) {
634  edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
635  << "'vstring mergeTreeFiles' and 'vstring mergeBinaryFiles' "
636  << "differ in size";
637  }
638  this->addHitStatistics(0, outFile, inFiles); // add hit info from tree 0 in 'infiles'
639  }
640  MillePedeVariablesIORoot millePedeIO;
641  millePedeIO.writeMillePedeVariables(theAlignables, outFile.c_str(), loop, false, ioerr);
642  if (ioerr) {
643  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
644  << "Problem " << ioerr << " writing MillePedeVariables";
645  ++result;
646  }
647 
648  aliIO.writeOrigRigidBodyAlignmentParameters(theAlignables, outFile.c_str(), loop, false, ioerr);
649  if (ioerr) {
650  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO" << "Problem " << ioerr
651  << " in writeOrigRigidBodyAlignmentParameters, " << loop;
652  ++result;
653  }
654  aliIO.writeAlignableAbsolutePositions(theAlignables, outFile.c_str(), loop, false, ioerr);
655  if (ioerr) {
656  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO" << "Problem " << ioerr
657  << " in writeAlignableAbsolutePositions, " << loop;
658  ++result;
659  }
660 
661  return result;
662 }
663 
664 //__________________________________________________________________________________________________
665 void MillePedeAlignmentAlgorithm::buildUserVariables(const std::vector<Alignable*> &alis) const
666 {
667  for (std::vector<Alignable*>::const_iterator iAli = alis.begin(); iAli != alis.end(); ++iAli) {
668  AlignmentParameters *params = (*iAli)->alignmentParameters();
669  if (!params) {
670  throw cms::Exception("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::buildUserVariables"
671  << "No parameters for alignable";
672  }
673  MillePedeVariables *userVars = new MillePedeVariables(params->size(), thePedeLabels->alignableLabel(*iAli));
674  params->setUserVariables(userVars);
675  }
676 }
677 
678 //__________________________________________________________________________________________________
679 unsigned int MillePedeAlignmentAlgorithm::decodeMode(const std::string &mode) const
680 {
681  if (mode == "full") {
683  } else if (mode == "mille") {
684  return myMilleBit; // + myPedeSteerBit; // sic! Including production of steerig file. NO!
685  } else if (mode == "pede") {
687  } else if (mode == "pedeSteer") {
688  return myPedeSteerBit;
689  } else if (mode == "pedeRun") {
690  return myPedeSteerBit + myPedeRunBit + myPedeReadBit; // sic! Including steering and reading of result.
691  } else if (mode == "pedeRead") {
692  return myPedeReadBit;
693  }
694 
695  throw cms::Exception("BadConfig")
696  << "Unknown mode '" << mode
697  << "', use 'full', 'mille', 'pede', 'pedeRun', 'pedeSteer' or 'pedeRead'.";
698 
699  return 0;
700 }
701 
702 //__________________________________________________________________________________________________
703 bool MillePedeAlignmentAlgorithm::addHitStatistics(int fromIov, const std::string &outFile,
704  const std::vector<std::string> &inFiles) const
705 {
706  bool allOk = true;
707  int ierr = 0;
708  MillePedeVariablesIORoot millePedeIO;
709  for (std::vector<std::string>::const_iterator iFile = inFiles.begin();
710  iFile != inFiles.end(); ++iFile) {
711  const std::string inFile(theDir + *iFile);
712  const std::vector<AlignmentUserVariables*> mpVars =
713  millePedeIO.readMillePedeVariables(theAlignables, inFile.c_str(), fromIov, ierr);
714  if (ierr || !this->addHits(theAlignables, mpVars)) {
715  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::addHitStatistics"
716  << "Error " << ierr << " reading from " << inFile
717  << ", tree " << fromIov << ", or problems in addHits";
718  allOk = false;
719  }
720  for (std::vector<AlignmentUserVariables*>::const_iterator i = mpVars.begin();
721  i != mpVars.end(); ++i) {
722  delete *i; // clean created objects
723  }
724  }
725 
726  return allOk;
727 }
728 
729 //__________________________________________________________________________________________________
730 bool MillePedeAlignmentAlgorithm::addHits(const std::vector<Alignable*> &alis,
731  const std::vector<AlignmentUserVariables*> &mpVars) const
732 {
733  bool allOk = (mpVars.size() == alis.size());
734  std::vector<AlignmentUserVariables*>::const_iterator iUser = mpVars.begin();
735  for (std::vector<Alignable*>::const_iterator iAli = alis.begin();
736  iAli != alis.end() && iUser != mpVars.end(); ++iAli, ++iUser) {
737  MillePedeVariables *mpVarNew = dynamic_cast<MillePedeVariables*>(*iUser);
738  AlignmentParameters *ps = (*iAli)->alignmentParameters();
739  MillePedeVariables *mpVarOld = (ps ? dynamic_cast<MillePedeVariables*>(ps->userVariables()) : 0);
740  if (!mpVarNew || !mpVarOld || mpVarOld->size() != mpVarNew->size()) {
741  allOk = false;
742  continue; // FIXME error etc.?
743  }
744 
745  mpVarOld->increaseHitsX(mpVarNew->hitsX());
746  mpVarOld->increaseHitsY(mpVarNew->hitsY());
747  }
748 
749  return allOk;
750 }
751 
752 //__________________________________________________________________________________________________
753 void MillePedeAlignmentAlgorithm::makeGlobDerivMatrix(const std::vector<float> &globalDerivativesx,
754  const std::vector<float> &globalDerivativesy,
755  TMatrixF &aGlobalDerivativesM)
756 {
757 
758  for (unsigned int i = 0; i < globalDerivativesx.size(); ++i) {
759  aGlobalDerivativesM(0,i) = globalDerivativesx[i];
760  aGlobalDerivativesM(1,i) = globalDerivativesy[i];
761  }
762 }
763 
764 //__________________________________________________________________________________________________
766 (TMatrixDSym &aHitCovarianceM, TMatrixF &aLocalDerivativesM, TMatrixF &aHitResidualsM,
767  TMatrixF &aGlobalDerivativesM) const
768 {
769  TMatrixDSymEigen myDiag(aHitCovarianceM);
770  TMatrixD aTranfoToDiagonalSystem = myDiag.GetEigenVectors();
771  TMatrixD aTranfoToDiagonalSystemInv = myDiag.GetEigenVectors( );
772  TMatrixF aTranfoToDiagonalSystemInvF = myDiag.GetEigenVectors( );
773  TMatrixD aMatrix = aTranfoToDiagonalSystemInv.Invert() * aHitCovarianceM * aTranfoToDiagonalSystem;
774  // Tranformation of matrix M is done by A^T*M*A, not A^{-1}*M*A.
775  // But here A^T == A^{-1}, so we would only save CPU by Transpose()...
776  // FIXME this - I guess simply use T(), not Transpose()...
777  // TMatrixD aMatrix = aTranfoToDiagonalSystemInv.Transpose() * aHitCovarianceM
778  // * aTranfoToDiagonalSystem;
779  aHitCovarianceM = TMatrixDSym(2, aMatrix.GetMatrixArray());
780  aTranfoToDiagonalSystemInvF.Invert();
781  //edm::LogInfo("Alignment") << "NEW HIT loca in matrix"<<aLocalDerivativesM(0,0);
782  aLocalDerivativesM = aTranfoToDiagonalSystemInvF * aLocalDerivativesM;
783 
784  //edm::LogInfo("Alignment") << "NEW HIT loca in matrix after diag:"<<aLocalDerivativesM(0,0);
785  aHitResidualsM = aTranfoToDiagonalSystemInvF * aHitResidualsM;
786  aGlobalDerivativesM = aTranfoToDiagonalSystemInvF * aGlobalDerivativesM;
787 }
788 
789 //__________________________________________________________________________________________________
792  unsigned int iVirtualMeas, TMatrixDSym &aHitCovarianceM,
793  TMatrixF &aHitResidualsM, TMatrixF &aLocalDerivativesM)
794 {
795  // This Method is valid for 1D measurements only
796 
797  const unsigned int xIndex = iVirtualMeas + refTrajPtr->numberOfHitMeas();
798  // Covariance into a TMatrixDSym
799 
800  //aHitCovarianceM = new TMatrixDSym(1);
801  aHitCovarianceM(0,0)=refTrajPtr->measurementErrors()[xIndex][xIndex];
802 
803  //theHitResidualsM= new TMatrixF(1,1);
804  aHitResidualsM(0,0)= refTrajPtr->measurements()[xIndex];
805 
806  // Local Derivatives into a TMatrixDSym (to use matrix operations)
807  const AlgebraicMatrix &locDerivMatrix = refTrajPtr->derivatives();
808  // theLocalDerivativeNumber = locDerivMatrix.num_col();
809 
810  //theLocalDerivativesM = new TMatrixF(1,locDerivMatrix.num_col());
811  for (int i = 0; i < locDerivMatrix.num_col(); ++i) {
812  aLocalDerivativesM(0,i) = locDerivMatrix[xIndex][i];
813  }
814 }
815 
816 //__________________________________________________________________________________________________
819  unsigned int iTrajHit, TMatrixDSym &aHitCovarianceM,
820  TMatrixF &aHitResidualsM, TMatrixF &aLocalDerivativesM)
821 {
822  // This Method is valid for 2D measurements only
823 
824  const unsigned int xIndex = iTrajHit*2;
825  const unsigned int yIndex = iTrajHit*2+1;
826  // Covariance into a TMatrixDSym
827 
828  //aHitCovarianceM = new TMatrixDSym(2);
829  aHitCovarianceM(0,0)=refTrajPtr->measurementErrors()[xIndex][xIndex];
830  aHitCovarianceM(0,1)=refTrajPtr->measurementErrors()[xIndex][yIndex];
831  aHitCovarianceM(1,0)=refTrajPtr->measurementErrors()[yIndex][xIndex];
832  aHitCovarianceM(1,1)=refTrajPtr->measurementErrors()[yIndex][yIndex];
833 
834  //theHitResidualsM= new TMatrixF(2,1);
835  aHitResidualsM(0,0)= refTrajPtr->measurements()[xIndex] - refTrajPtr->trajectoryPositions()[xIndex];
836  aHitResidualsM(1,0)= refTrajPtr->measurements()[yIndex] - refTrajPtr->trajectoryPositions()[yIndex];
837 
838  // Local Derivatives into a TMatrixDSym (to use matrix operations)
839  const AlgebraicMatrix &locDerivMatrix = refTrajPtr->derivatives();
840  // theLocalDerivativeNumber = locDerivMatrix.num_col();
841 
842  //theLocalDerivativesM = new TMatrixF(2,locDerivMatrix.num_col());
843  for (int i = 0; i < locDerivMatrix.num_col(); ++i) {
844  aLocalDerivativesM(0,i) = locDerivMatrix[xIndex][i];
845  aLocalDerivativesM(1,i) = locDerivMatrix[yIndex][i];
846  }
847 }
848 
849 //__________________________________________________________________________________________________
852  unsigned int iTrajHit, const std::vector<int> &globalLabels,
853  const std::vector<float> &globalDerivativesX,
854  const std::vector<float> &globalDerivativesY)
855 {
856  const ConstRecHitPointer aRecHit(refTrajPtr->recHits()[iTrajHit]);
857 
858  if((aRecHit)->dimension() == 1) {
859  return this->callMille1D(refTrajPtr, iTrajHit, globalLabels, globalDerivativesX);
860  } else {
861  return this->callMille2D(refTrajPtr, iTrajHit, globalLabels,
862  globalDerivativesX, globalDerivativesY);
863  }
864 }
865 
866 
867 //__________________________________________________________________________________________________
870  unsigned int iTrajHit, const std::vector<int> &globalLabels,
871  const std::vector<float> &globalDerivativesX)
872 {
873  const ConstRecHitPointer aRecHit(refTrajPtr->recHits()[iTrajHit]);
874  const unsigned int xIndex = iTrajHit*2; // the even ones are local x
875 
876  // local derivatives
877  const AlgebraicMatrix &locDerivMatrix = refTrajPtr->derivatives();
878  const int nLocal = locDerivMatrix.num_col();
879  std::vector<float> localDerivatives(nLocal);
880  for (unsigned int i = 0; i < localDerivatives.size(); ++i) {
881  localDerivatives[i] = locDerivMatrix[xIndex][i];
882  }
883 
884  // residuum and error
885  float residX = refTrajPtr->measurements()[xIndex] - refTrajPtr->trajectoryPositions()[xIndex];
886  float hitErrX = TMath::Sqrt(refTrajPtr->measurementErrors()[xIndex][xIndex]);
887 
888  // number of global derivatives
889  const int nGlobal = globalDerivativesX.size();
890 
891  // &(localDerivatives[0]) etc. are valid - as long as vector is not empty
892  // cf. http://www.parashift.com/c++-faq-lite/containers.html#faq-34.3
893  theMille->mille(nLocal, &(localDerivatives[0]), nGlobal, &(globalDerivativesX[0]),
894  &(globalLabels[0]), residX, hitErrX);
895 
896  if (theMonitor) {
897  theMonitor->fillDerivatives(aRecHit, &(localDerivatives[0]), nLocal,
898  &(globalDerivativesX[0]), nGlobal, &(globalLabels[0]));
899  theMonitor->fillResiduals(aRecHit, refTrajPtr->trajectoryStates()[iTrajHit],
900  iTrajHit, residX, hitErrX, false);
901  }
902 
903  return 1;
904 }
905 
906 //__________________________________________________________________________________________________
909  unsigned int iTrajHit, const std::vector<int> &globalLabels,
910  const std::vector<float> &globalDerivativesx,
911  const std::vector<float> &globalDerivativesy)
912 {
913  const ConstRecHitPointer aRecHit(refTrajPtr->recHits()[iTrajHit]);
914 
915  if((aRecHit)->dimension() != 2) {
916  edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::callMille2D"
917  << "You try to call method for 2D hits for a "
918  << (aRecHit)->dimension()
919  << "D Hit. Hit gets ignored!";
920  return -1;
921  }
922 
923  TMatrixDSym aHitCovarianceM(2);
924  TMatrixF aHitResidualsM(2,1);
925  TMatrixF aLocalDerivativesM(2, refTrajPtr->derivatives().num_col());
926  // below method fills above 3 matrices
927  this->addRefTrackData2D(refTrajPtr, iTrajHit, aHitCovarianceM,aHitResidualsM,aLocalDerivativesM);
928  TMatrixF aGlobalDerivativesM(2,globalDerivativesx.size());
929  this->makeGlobDerivMatrix(globalDerivativesx, globalDerivativesy, aGlobalDerivativesM);
930 
931  // calculates correlation between Hit measurements
932  // FIXME: Should take correlation (and resulting transformation) from original hit,
933  // not 2x2 matrix from ReferenceTrajectory: That can come from error propagation etc.!
934  const double corr = aHitCovarianceM(0,1) / sqrt(aHitCovarianceM(0,0) * aHitCovarianceM(1,1));
935  if (theMonitor) theMonitor->fillCorrelations2D(corr, aRecHit);
936  bool diag = false; // diagonalise only tracker TID, TEC
937  switch(aRecHit->geographicalId().subdetId()) {
938  case SiStripDetId::TID:
939  case SiStripDetId::TEC:
940  if (aRecHit->geographicalId().det() == DetId::Tracker && TMath::Abs(corr) > theMaximalCor2D) {
941  this->diagonalize(aHitCovarianceM, aLocalDerivativesM, aHitResidualsM, aGlobalDerivativesM);
942  diag = true;
943  }
944  break;
945  default:;
946  }
947 
948  float newResidX = aHitResidualsM(0,0);
949  float newResidY = aHitResidualsM(1,0);
950  float newHitErrX = TMath::Sqrt(aHitCovarianceM(0,0));
951  float newHitErrY = TMath::Sqrt(aHitCovarianceM(1,1));
952  float *newLocalDerivsX = aLocalDerivativesM[0].GetPtr();
953  float *newLocalDerivsY = aLocalDerivativesM[1].GetPtr();
954  float *newGlobDerivsX = aGlobalDerivativesM[0].GetPtr();
955  float *newGlobDerivsY = aGlobalDerivativesM[1].GetPtr();
956  const int nLocal = aLocalDerivativesM.GetNcols();
957  const int nGlobal = aGlobalDerivativesM.GetNcols();
958 
959  if (diag && (newHitErrX > newHitErrY)) { // also for 2D hits?
960  // measurement with smaller error is x-measurement (for !is2D do not fill y-measurement):
961  std::swap(newResidX, newResidY);
962  std::swap(newHitErrX, newHitErrY);
963  std::swap(newLocalDerivsX, newLocalDerivsY);
964  std::swap(newGlobDerivsX, newGlobDerivsY);
965  }
966 
967  // &(globalLabels[0]) is valid - as long as vector is not empty
968  // cf. http://www.parashift.com/c++-faq-lite/containers.html#faq-34.3
969  theMille->mille(nLocal, newLocalDerivsX, nGlobal, newGlobDerivsX,
970  &(globalLabels[0]), newResidX, newHitErrX);
971 
972  if (theMonitor) {
973  theMonitor->fillDerivatives(aRecHit, newLocalDerivsX, nLocal, newGlobDerivsX, nGlobal,
974  &(globalLabels[0]));
975  theMonitor->fillResiduals(aRecHit, refTrajPtr->trajectoryStates()[iTrajHit],
976  iTrajHit, newResidX, newHitErrX, false);
977  }
978  const bool isReal2DHit = this->is2D(aRecHit); // strip is 1D (except matched hits)
979  if (isReal2DHit) {
980  theMille->mille(nLocal, newLocalDerivsY, nGlobal, newGlobDerivsY,
981  &(globalLabels[0]), newResidY, newHitErrY);
982  if (theMonitor) {
983  theMonitor->fillDerivatives(aRecHit, newLocalDerivsY, nLocal, newGlobDerivsY, nGlobal,
984  &(globalLabels[0]));
985  theMonitor->fillResiduals(aRecHit, refTrajPtr->trajectoryStates()[iTrajHit],
986  iTrajHit, newResidY, newHitErrY, true);// true: y
987  }
988  }
989 
990  return (isReal2DHit ? 2 : 1);
991 }
992 
993 //__________________________________________________________________________________________________
995 ::addVirtualMeas(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iVirtualMeas)
996 {
997  TMatrixDSym aHitCovarianceM(1);
998  TMatrixF aHitResidualsM(1,1);
999  TMatrixF aLocalDerivativesM(1, refTrajPtr->derivatives().num_col());
1000  // below method fills above 3 'matrices'
1001  this->addRefTrackVirtualMeas1D(refTrajPtr, iVirtualMeas, aHitCovarianceM, aHitResidualsM, aLocalDerivativesM);
1002 
1003  // no global parameters (use dummy 0)
1004  TMatrixF aGlobalDerivativesM(1,1);
1005  aGlobalDerivativesM(0,0) = 0;
1006 
1007  float newResidX = aHitResidualsM(0,0);
1008  float newHitErrX = TMath::Sqrt(aHitCovarianceM(0,0));
1009  float *newLocalDerivsX = aLocalDerivativesM[0].GetPtr();
1010  float *newGlobDerivsX = aGlobalDerivativesM[0].GetPtr();
1011  const int nLocal = aLocalDerivativesM.GetNcols();
1012  const int nGlobal = 0;
1013 
1014  theMille->mille(nLocal, newLocalDerivsX, nGlobal, newGlobDerivsX,
1015  &nGlobal, newResidX, newHitErrX);
1016 }
1017 
1018 //____________________________________________________
1020  const TkFittedLasBeamCollection &lasBeams,
1021  const TsosVectorCollection &lasBeamTsoses)
1022 {
1023  TsosVectorCollection::const_iterator iTsoses = lasBeamTsoses.begin();
1024  for(TkFittedLasBeamCollection::const_iterator iBeam = lasBeams.begin(), iEnd = lasBeams.end();
1025  iBeam != iEnd; ++iBeam, ++iTsoses){ // beam/tsoses parallel!
1026 
1027  edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::addLaserData"
1028  << "Beam " << iBeam->getBeamId() << " with "
1029  << iBeam->parameters().size() << " parameters and "
1030  << iBeam->getData().size() << " hits.\n There are "
1031  << iTsoses->size() << " TSOSes.";
1032 
1033  this->addLasBeam(eventInfo, *iBeam, *iTsoses);
1034  }
1035 }
1036 
1037 //____________________________________________________
1039  const TkFittedLasBeam &lasBeam,
1040  const std::vector<TrajectoryStateOnSurface> &tsoses)
1041 {
1042  AlignmentParameters *dummyPtr = 0; // for globalDerivativesHierarchy()
1043  std::vector<float> lasLocalDerivsX; // buffer for local derivatives
1044  const unsigned int beamLabel = thePedeLabels->lasBeamLabel(lasBeam.getBeamId());// for global par
1045 
1046  for (unsigned int iHit = 0; iHit < tsoses.size(); ++iHit) {
1047  if (!tsoses[iHit].isValid()) continue;
1048  // clear buffer
1049  theFloatBufferX.clear();
1050  theFloatBufferY.clear();
1051  theIntBuffer.clear();
1052  lasLocalDerivsX.clear();
1053  // get alignables and global parameters
1054  const SiStripLaserRecHit2D &hit = lasBeam.getData()[iHit];
1056  this->globalDerivativesHierarchy(eventInfo,
1057  tsoses[iHit], lasAli, lasAli,
1059  // fill derivatives vector from derivatives matrix
1060  for (unsigned int nFitParams = 0;
1061  nFitParams < static_cast<unsigned int>(lasBeam.parameters().size());
1062  ++nFitParams) {
1063  const float derivative = lasBeam.derivatives()[iHit][nFitParams];
1064  if (nFitParams < lasBeam.firstFixedParameter()) { // first local beam parameters
1065  lasLocalDerivsX.push_back(derivative);
1066  } else { // now global ones
1067  const unsigned int numPar = nFitParams - lasBeam.firstFixedParameter();
1068  theIntBuffer.push_back(thePedeLabels->parameterLabel(beamLabel, numPar));
1069  theFloatBufferX.push_back(derivative);
1070  }
1071  } // end loop over parameters
1072 
1073  const float residual = hit.localPosition().x() - tsoses[iHit].localPosition().x();
1074  // error from file or assume 0.003
1075  const float error = 0.003; // hit.localPositionError().xx(); sqrt???
1076 
1077  theMille->mille(lasLocalDerivsX.size(), &(lasLocalDerivsX[0]), theFloatBufferX.size(),
1078  &(theFloatBufferX[0]), &(theIntBuffer[0]), residual, error);
1079  } // end of loop over hits
1080 
1081  theMille->end();
1082 }
1083 
1085 {
1086  // do some printing, if requested
1087  const bool doOutputOnStdout(pxbSurveyCfg.getParameter<bool>("doOutputOnStdout"));
1088  if (doOutputOnStdout) std::cout << "# Output from addPxbSurvey follows below because doOutputOnStdout is set to True" << std::endl;
1089 
1090  // instantiate a dicer object
1091  SurveyPxbDicer dicer(pxbSurveyCfg.getParameter<std::vector<edm::ParameterSet> >("toySurveyParameters"), pxbSurveyCfg.getParameter<unsigned int>("toySurveySeed"));
1092  std::ofstream outfile(pxbSurveyCfg.getUntrackedParameter<std::string>("toySurveyFile").c_str());
1093 
1094  // read data from file
1095  std::vector<SurveyPxbImageLocalFit> measurements;
1096  std::string filename(pxbSurveyCfg.getParameter<edm::FileInPath>("infile").fullPath());
1098 
1099  // loop over photographs (=measurements) and perform the fit
1100  for(std::vector<SurveyPxbImageLocalFit>::size_type i=0; i!=measurements.size(); i++)
1101  {
1102  if (doOutputOnStdout) std::cout << "Module " << i << ": ";
1103 
1104  // get the Alignables and their surfaces
1105  AlignableDetOrUnitPtr mod1(theAlignableNavigator->alignableFromDetId(measurements[i].getIdFirst()));
1106  AlignableDetOrUnitPtr mod2(theAlignableNavigator->alignableFromDetId(measurements[i].getIdSecond()));
1107  const AlignableSurface& surf1 = mod1->surface();
1108  const AlignableSurface& surf2 = mod2->surface();
1109 
1110  // the position of the fiducial points in local frame of a PXB module
1111  const LocalPoint fidpoint0(-0.91,+3.30);
1112  const LocalPoint fidpoint1(+0.91,+3.30);
1113  const LocalPoint fidpoint2(+0.91,-3.30);
1114  const LocalPoint fidpoint3(-0.91,-3.30);
1115 
1116  // We choose the local frame of the first module as reference,
1117  // so take the fidpoints of the second module and calculate their
1118  // positions in the reference frame
1119  const GlobalPoint surf2point0(surf2.toGlobal(fidpoint0));
1120  const GlobalPoint surf2point1(surf2.toGlobal(fidpoint1));
1121  const LocalPoint fidpoint0inSurf1frame(surf1.toLocal(surf2point0));
1122  const LocalPoint fidpoint1inSurf1frame(surf1.toLocal(surf2point1));
1123 
1124  // Create the vector for the fit
1126  fidpointvec.push_back(fidpoint0inSurf1frame);
1127  fidpointvec.push_back(fidpoint1inSurf1frame);
1128  fidpointvec.push_back(fidpoint2);
1129  fidpointvec.push_back(fidpoint3);
1130 
1131  // if toy survey is requested, dice the values now
1132  if (pxbSurveyCfg.getParameter<bool>("doToySurvey"))
1133  {
1134  dicer.doDice(fidpointvec,measurements[i].getIdPair(), outfile);
1135  }
1136 
1137  // do the fit
1138  measurements[i].doFit(fidpointvec, thePedeLabels->alignableLabel(mod1), thePedeLabels->alignableLabel(mod2));
1139  SurveyPxbImageLocalFit::localpars_t a; // local pars from fit
1140  a = measurements[i].getLocalParameters();
1141  const SurveyPxbImageLocalFit::value_t chi2 = measurements[i].getChi2();
1142 
1143  // do some reporting, if requested
1144  if (doOutputOnStdout)
1145  {
1146  std::cout << "a: " << a[0] << ", " << a[1] << ", " << a[2] << ", " << a[3]
1147  << " S= " << sqrt(a[2]*a[2]+a[3]*a[3])
1148  << " phi= " << atan(a[3]/a[2])
1149  << " chi2= " << chi2 << std::endl;
1150  }
1151  if (theMonitor)
1152  {
1154  theMonitor->fillPxbSurveyHistsLocalPars(a[0],a[1],sqrt(a[2]*a[2]+a[3]*a[3]),atan(a[3]/a[2]));
1155  }
1156 
1157  // pass the results from the local fit to mille
1159  {
1160  theMille->mille((int)measurements[i].getLocalDerivsSize(),
1161  measurements[i].getLocalDerivsPtr(j),
1162  (int)measurements[i].getGlobalDerivsSize(),
1163  measurements[i].getGlobalDerivsPtr(j),
1164  measurements[i].getGlobalDerivsLabelPtr(j),
1165  measurements[i].getResiduum(j),
1166  measurements[i].getSigma(j));
1167  }
1168  theMille->end();
1169  }
1170  outfile.close();
1171 }
1172 
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
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
void increaseHitsX(unsigned int add=1)
increase hits for x-measurement
bool read(std::vector< Alignable * > &alignables, bool setUserVars)
Definition: PedeReader.cc:50
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
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
const ConstTrajTrackPairCollection & trajTrackPairs_
const TsosVectorCollection * tkLasBeamTsoses_
might be null!
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
void fillPxbSurveyHistsChi2(const float &chi2)
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
tuple files
Definition: linker.py:146
void addLaserData(const EventInfo &eventInfo, const TkFittedLasBeamCollection &tkLasBeams, const TsosVectorCollection &tkLasBeamTsoses)
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)
virtual bool setParametersForRunRange(const RunRange &runrange)
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
T sqrt(T t)
Definition: SSEVec.h:46
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:759
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
bool is2D(const TransientTrackingRecHit::ConstRecHitPointer &recHit) const
true if hit belongs to 2D detector (currently tracker specific)
void diagonalize(TMatrixDSym &aHitCovarianceM, TMatrixF &aLocalDerivativesM, TMatrixF &aHitResidualsM, TMatrixF &theGlobalDerivativesM) const
int j
Definition: DBlmapReader.cc:9
std::string buildMasterSteer(const std::vector< std::string > &binaryFiles)
construct (and return name of) master steering file from config, binaryFiles etc. ...
Definition: PedeSteerer.cc:719
Container::value_type value_type
int addMeasurementData(const EventInfo &eventInfo, const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iHit, AlignmentParameters *&params)
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:9
bool addHits(const std::vector< Alignable * > &alis, const std::vector< AlignmentUserVariables * > &mpVars) const
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
std::vector< ConstRecHitPointer > ConstRecHitContainer
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::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 setUserVariables(AlignmentUserVariables *auv)
Set pointer to user variables.
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:671
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:209
Class to hold one picture of the BPix survey.
const TkFittedLasBeamCollection * tkLasBeams_
std::vector< value_t > localpars_t
void fillTrack(const reco::Track *track)
virtual void run(const edm::EventSetup &setup, const EventInfo &eventInfo)
Run the algorithm on trajectories and tracks.
define run information passed to algorithms (in endRun)
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.
std::vector< TkFittedLasBeam > TkFittedLasBeamCollection
double a
Definition: hdecay.h:121
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
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:171
std::pair< unsigned int, unsigned int > addReferenceTrajectory(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) ...
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
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)
T x() const
Definition: PV3DBase.h:61
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:56
const align::Alignables & alignables(void) const
get all alignables
std::vector< ConstTrajTrackPair > ConstTrajTrackPairCollection
A derivative(const A &_)
Definition: Derivative.h:18
define event information passed to algorithms
virtual void terminate()
Call at end of job.