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