CMS 3D CMS Logo

HIPAlignmentAlgorithm.cc
Go to the documentation of this file.
1 #include <fstream>
2 
3 #include "TFile.h"
4 #include "TTree.h"
5 #include "TRandom.h"
6 #include "TFormula.h"
7 #include "TMath.h"
8 
14 
17 
34 
40 
41 
43 
44 // using namespace std;
45 
46 // Constructor ----------------------------------------------------------------
47 
50  surveyResiduals_(cfg.getUntrackedParameter<std::vector<std::string> >("surveyResiduals"))
51 {
52 
53  // parse parameters
54 
55  verbose = cfg.getParameter<bool>("verbosity");
56  outpath = cfg.getParameter<std::string>("outpath");
57  outfile = cfg.getParameter<std::string>("outfile");
58  outfile2 = cfg.getParameter<std::string>("outfile2");
59  struefile = cfg.getParameter<std::string>("trueFile");
60  smisalignedfile = cfg.getParameter<std::string>("misalignedFile");
61  salignedfile = cfg.getParameter<std::string>("alignedFile");
62  siterationfile = cfg.getParameter<std::string>("iterationFile");
63  suvarfile = cfg.getParameter<std::string>("uvarFile");
64  sparameterfile = cfg.getParameter<std::string>("parameterFile");
65  ssurveyfile = cfg.getParameter<std::string>("surveyFile");
66 
67  outfile =outpath+outfile;//Eventwise tree
68  outfile2 =outpath+outfile2;//Alignablewise tree
69  struefile =outpath+struefile;
70  smisalignedfile=outpath+smisalignedfile;
71  salignedfile =outpath+salignedfile;
72  siterationfile =outpath+siterationfile;
73  suvarfile =outpath+suvarfile;
74  sparameterfile =outpath+sparameterfile;
75  ssurveyfile =outpath+ssurveyfile;
76 
77  // parameters for APE
78  theApplyAPE = cfg.getParameter<bool>("applyAPE");
79  themultiIOV = cfg.getParameter<bool>("multiIOV");
80  theAPEParameterSet = cfg.getParameter<std::vector<edm::ParameterSet> >("apeParam");
81  theIOVrangeSet = cfg.getParameter<std::vector<unsigned> >("IOVrange");
82 
83  theMaxAllowedHitPull = cfg.getParameter<double>("maxAllowedHitPull");
84  theMinimumNumberOfHits = cfg.getParameter<int>("minimumNumberOfHits");
85  theMaxRelParameterError = cfg.getParameter<double>("maxRelParameterError");
86 
87  // for collector mode (parallel processing)
88  isCollector=cfg.getParameter<bool>("collectorActive");
89  theCollectorNJobs=cfg.getParameter<int>("collectorNJobs");
90  theCollectorPath=cfg.getParameter<std::string>("collectorPath");
91  theFillTrackMonitoring=cfg.getUntrackedParameter<bool>("fillTrackMonitoring");
92 
93  if (isCollector) edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm] Collector mode";
94 
95  theEventPrescale = cfg.getParameter<int>("eventPrescale");
97 
98 //CY : added options
99  trackPs = cfg.getParameter<bool>("UsePreSelection");
100  trackWt = cfg.getParameter<bool>("UseReweighting");
101  Scale = cfg.getParameter<double>("Weight");
102  uniEta = cfg.getParameter<bool>("UniformEta");
103  IsCollision = cfg.getParameter<bool>("isCollision");
104  SetScanDet=cfg.getParameter<std::vector<double> >("setScanDet");
105  col_cut = cfg.getParameter<double>("CLAngleCut");
106  cos_cut = cfg.getParameter<double>("CSAngleCut");
107 
108  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm] constructed.";
109 
110 }
111 
112 // Call at beginning of job ---------------------------------------------------
113 
114 void
117  AlignmentParameterStore* store )
118 {
119  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm] Initializing...";
120 
121  alignableObjectId_ = std::make_unique<AlignableObjectId>
123 
124  for (const auto& level: surveyResiduals_) {
125  theLevels.push_back(alignableObjectId_->stringToId(level));
126  }
127 
128  edm::ESHandle<Alignments> globalPositionRcd;
129 
130  const edm::ValidityInterval & validity = setup.get<TrackerAlignmentRcd>().validityInterval();
131  const edm::IOVSyncValue first1 = validity.first();
132  unsigned int firstrun = first1.eventID().run();
133  if(themultiIOV){
134  if(theIOVrangeSet.size()!=1){
135  bool findMatchIOV=false;
136  for (unsigned int iovl = 0; iovl <theIOVrangeSet.size(); iovl++){
137  if(firstrun == theIOVrangeSet.at(iovl)){
138  std::string iovapp = std::to_string(firstrun);
139  iovapp.append(".root");
140  iovapp.insert(0,"_");
141  salignedfile.replace(salignedfile.end()-5, salignedfile.end(),iovapp);
142  siterationfile.replace(siterationfile.end()-5, siterationfile.end(),iovapp);
143  //sparameterfile.replace(sparameterfile.end()-5, sparameterfile.end(),iovapp);
144  if(isCollector){
145  outfile2.replace(outfile2.end()-5, outfile2.end(),iovapp);
146  ssurveyfile.replace(ssurveyfile.end()-5, ssurveyfile.end(),iovapp);
147  suvarfile.replace(suvarfile.end()-5, suvarfile.end(),iovapp);
148  }
149 
150  findMatchIOV=true;
151  break;
152  }
153  }
154  if(!findMatchIOV){
155  edm::LogWarning("Alignment") << "error! Didn't find the matched IOV file!";
156  }
157  }
158  else{
159  std::string iovapp = std::to_string(theIOVrangeSet.at(0));
160  iovapp.append(".root");
161  iovapp.insert(0,"_");
162  salignedfile.replace(salignedfile.end()-5, salignedfile.end(),iovapp);
163  siterationfile.replace(siterationfile.end()-5, siterationfile.end(),iovapp);
164  }
165  }
166 
167  // accessor Det->AlignableDet
168  if ( !muon )
170  else if ( !tracker )
172  else
173  theAlignableDetAccessor = new AlignableNavigator(tracker, muon);
174 
175  // set alignmentParameterStore
177 
178  // get alignables
180 
181  // clear theAPEParameters, if necessary
182  theAPEParameters.clear();
183 
184  // get APE parameters
185  if(theApplyAPE){
186  AlignmentParameterSelector selector(tracker, muon);
187  for (std::vector<edm::ParameterSet>::const_iterator setiter = theAPEParameterSet.begin();
188  setiter != theAPEParameterSet.end();
189  ++setiter) {
190  std::vector<Alignable*> alignables;
191 
192  selector.clear();
193  edm::ParameterSet selectorPSet = setiter->getParameter<edm::ParameterSet>("Selector");
194  std::vector<std::string> alignParams = selectorPSet.getParameter<std::vector<std::string> >("alignParams");
195  if (alignParams.size() == 1 && alignParams[0] == std::string("selected")) {
196  alignables = theAlignables;
197  }
198  else {
199  selector.addSelections(selectorPSet);
200  alignables = selector.selectedAlignables();
201  }
202 
203  std::vector<double> apeSPar = setiter->getParameter<std::vector<double> >("apeSPar");
204  std::vector<double> apeRPar = setiter->getParameter<std::vector<double> >("apeRPar");
205  std::string function = setiter->getParameter<std::string>("function");
206 
207  if (apeSPar.size() != 3 || apeRPar.size() != 3)
208  throw cms::Exception("BadConfig") << "apeSPar and apeRPar must have 3 values each" << std::endl;
209 
210  for (std::vector<double>::const_iterator i = apeRPar.begin(); i != apeRPar.end(); ++i) {
211  apeSPar.push_back(*i);
212  }
213 
214  if (function == std::string("linear")) {
215  apeSPar.push_back(0.); // c.f. note in calcAPE
216  }
217  else if (function == std::string("exponential")) {
218  apeSPar.push_back(1.); // c.f. note in calcAPE
219  }
220  else if (function == std::string("step")) {
221  apeSPar.push_back(2.); // c.f. note in calcAPE
222  }
223  else {
224  throw cms::Exception("BadConfig") << "APE function must be \"linear\" or \"exponential\"." << std::endl;
225  }
226 
227  theAPEParameters.push_back(std::pair<std::vector<Alignable*>, std::vector<double> >(alignables, apeSPar));
228  }
229  }
230 }
231 
232 // Call at new loop -------------------------------------------------------------
234 {
235 
236  // iterate over all alignables and attach user variables
237  for( std::vector<Alignable*>::const_iterator it=theAlignables.begin();
238  it!=theAlignables.end(); it++ )
239  {
240  AlignmentParameters* ap = (*it)->alignmentParameters();
241  int npar=ap->numSelected();
242  HIPUserVariables* userpar = new HIPUserVariables(npar);
243  ap->setUserVariables(userpar);
244  }
245 
246  // try to read in alignment parameters from a previous iteration
247  //
248  AlignablePositions theAlignablePositionsFromFile =
250  salignedfile.c_str(),-1,ioerr);
251 
252  int numAlignablesFromFile = theAlignablePositionsFromFile.size();
253 
254  if (numAlignablesFromFile==0) { // file not there: first iteration
255 
256  // set iteration number to 1
257  if (isCollector) theIteration=0;
258  else theIteration=1;
259  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm] File not found => iteration "<<theIteration;
260 
261  // get true (de-misaligned positions) and write to root file
262  // hardcoded iteration=1
263 // theIO.writeAlignableOriginalPositions(theAlignables,
264 // struefile.c_str(),1,false,ioerr);
265 
266  // get misaligned positions and write to root file
267  // hardcoded iteration=1
268 // theIO.writeAlignableAbsolutePositions(theAlignables,
269 // smisalignedfile.c_str(),1,false,ioerr);
270 
271  }
272 
273  else { // there have been previous iterations
274 
275  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm] Alignables Read "
276  << numAlignablesFromFile;
277 
278  // get iteration number from file
282 
283  // increase iteration
284  if(ioerr==0){
285  theIteration++;
286  edm::LogWarning("Alignment") <<"[HIPAlignmentAlgorithm] Iteration increased by one!";
287  }
288 
289  // now apply psotions of file from prev iteration
290  edm::LogInfo("Alignment") <<"[HIPAlignmentAlgorithm] Apply positions from file ...";
292  theAlignablePositionsFromFile,ioerr);
293 
294  }
295 
296  edm::LogWarning("Alignment") <<"[HIPAlignmentAlgorithm] Current Iteration number: "
297  << theIteration;
298 
299 
300  // book root trees
301  bookRoot();
302 
303  // set alignment position error
305 
306  // run collector job if we are in parallel mode
307  if (isCollector) collector();
308 
309 }
310 
311 // Call at end of job ---------------------------------------------------------
312 
314 {
315 
316  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm] Terminating";
317 
318  // calculating survey residuals
319  if (theLevels.size() > 0 )
320  {
321  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm] Using survey constraint";
322 
323  unsigned int nAlignable = theAlignables.size();
324  edm::ESHandle<TrackerTopology> tTopoHandle;
325  iSetup.get<IdealGeometryRecord>().get(tTopoHandle);
326  const TrackerTopology* const tTopo = tTopoHandle.product();
327  for (unsigned int i = 0; i < nAlignable; ++i)
328  {
329  const Alignable* ali = theAlignables[i];
330 
332 
333  HIPUserVariables* uservar =
334  dynamic_cast<HIPUserVariables*>(ap->userVariables());
335  int nhit = uservar->nhit;
336 
337  // get position
338  std::pair<int,int> tl = theAlignmentParameterStore->typeAndLayer(ali, tTopo);
339  int tmp_Type = tl.first;
340  int tmp_Layer = tl.second;
341  GlobalPoint pos = ali->surface().position();
342  float tmpz = pos.z();
343  if(nhit< 1500 || (tmp_Type==5 && tmp_Layer==4 && fabs(tmpz)>90)){
344  for (unsigned int l = 0; l < theLevels.size(); ++l)
345  {
346  SurveyResidual res(*ali, theLevels[l], true);
347 
348  if ( res.valid() )
349  {
350  AlgebraicSymMatrix invCov = res.inverseCovariance();
351 
352  // variable for tree
353  AlgebraicVector sensResid = res.sensorResidual();
354  m3_Id = ali->id();
355  m3_ObjId = theLevels[l];
356  m3_par[0] = sensResid[0]; m3_par[1] = sensResid[1]; m3_par[2] = sensResid[2];
357  m3_par[3] = sensResid[3]; m3_par[4] = sensResid[4]; m3_par[5] = sensResid[5];
358 
359  uservar->jtvj += invCov;
360  uservar->jtve += invCov * sensResid;
361 
362  theTree3->Fill();
363  }
364  }
365  }
366 
367  }
368  }
369 
370  // write user variables
372 
373  if(!isCollector) //don't store userVariable in main, to save time
375  theIteration,false,ioerr);
376 
377  // now calculate alignment corrections ...
378  int ialigned=0;
379  // iterate over alignment parameters
380  for(std::vector<Alignable*>::const_iterator
381  it=theAlignables.begin(); it!=theAlignables.end(); it++) {
382  Alignable* ali=(*it);
383  // Alignment parameters
385 
386 
387  if (SetScanDet.at(0)!=0){
388  edm::LogWarning("Alignment") << "********Starting Scan*********";
389  edm::LogWarning("Alignment") <<"det ID="<<SetScanDet.at(0)<<", starting position="<<SetScanDet.at(1)<<", step="<<SetScanDet.at(2)<<", currentDet = "<<ali->id();}
390 
391  if((SetScanDet.at(0)!=0)&&(SetScanDet.at(0)!=1)&&(ali->id()!=SetScanDet.at(0)))continue;
392 
393  bool test = calcParameters(ali,SetScanDet.at(0),SetScanDet.at(1),SetScanDet.at(2));
394  // if successful, apply parameters
395  if (test){
397  // set these parameters 'valid'
398  ali->alignmentParameters()->setValid(true);
399  // increase counter
400  ialigned++;
401  }
402  else par->setValid(false);
403 // }
404 
405  }
406 //end looping over alignables
407 
408  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::terminate] Aligned units: " << ialigned;
409 
410  // fill alignable wise root tree
411  fillRoot(iSetup);
412 
413  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm] Writing aligned parameters to file: " << theAlignables.size()<<", for Iteration "<<theIteration;
414 
415 
416  // write user variables
417  if(isCollector)
419  theIteration,false,ioerr);
420 
421  // write new absolute positions to disk
423  salignedfile.c_str(),theIteration,false,ioerr);
424 
425  // write alignment parameters to disk
426  //theIO.writeAlignmentParameters(theAlignables,
427  // sparameterfile.c_str(),theIteration,false,ioerr);
428 
429  // write iteration number to file
431 
432  // write out trees and close root file
433 
434  // eventwise tree
436  theFile->cd();
437  theTree->Write();
438  hitTree->Write();
439  delete theFile;
440  }
441 
442  if (theLevels.size() > 0){
443  theFile3->cd();
444  theTree3->Write();
445  delete theFile3;
446  }
447 
448  // alignable-wise tree is only filled once
449  if (theIteration==1) { // only for 1st iteration
450  theFile2->cd();
451  theTree2->Write();
452  delete theFile2;
453  }
454 
455 }
456 
458  const Alignable* ali,
459  const TrajectoryStateOnSurface & tsos,
461  double hitwt)
462 {
463  static const unsigned int hitDim = 1;
464 
465  // get trajectory impact point
466  LocalPoint alvec = tsos.localPosition();
467  AlgebraicVector pos(hitDim);
468  pos[0] = alvec.x();
469 
470 // CY: get trajectory impact angle
471 // LocalVector v = tsos.localDirection();
472 // double proj_z = v.dot(LocalVector(0,0,1));
473 // edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::processHit1D]detectorId="<<ali->id();
474 // edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::processHit1D]proj_Z="<<proj_z;
475 
476  // get impact point covariance
477  AlgebraicSymMatrix ipcovmat(hitDim);
478  ipcovmat[0][0] = tsos.localError().positionError().xx();
479 
480  // get hit local position and covariance
481  AlgebraicVector coor(hitDim);
482  coor[0] = hit->localPosition().x();
483 
484  AlgebraicSymMatrix covmat(hitDim);
485  covmat[0][0] = hit->localPositionError().xx();
486 
487  // add hit and impact point covariance matrices
488  covmat = covmat + ipcovmat;
489 
490  // calculate the x pull of this hit
491  double xpull = 0.;
492  if (covmat[0][0] != 0.) xpull = (pos[0] - coor[0])/sqrt(fabs(covmat[0][0]));
493 
494  // get Alignment Parameters
495  AlignmentParameters* params = ali->alignmentParameters();
496  // get derivatives
497  AlgebraicMatrix derivs2D = params->selectedDerivatives(tsos, alidet);
498  // calculate user parameters
499  int npar = derivs2D.num_row();
500  AlgebraicMatrix derivs(npar, hitDim, 0);
501 
502  for (int ipar=0;ipar<npar;ipar++) {
503  for (unsigned int idim=0;idim<hitDim;idim++) {
504  derivs[ipar][idim] = derivs2D[ipar][idim];
505  }
506  }
507  // invert covariance matrix
508  int ierr;
509  // if (covmat[1][1]>4.0) nhitDim = hitDim;
510 
511  covmat.invert(ierr);
512  if (ierr != 0) {
513  edm::LogError("Alignment") << "Matrix inversion failed!";
514  return false;
515  }
516 
517  bool useThisHit = (theMaxAllowedHitPull <= 0.);
518 
519  useThisHit |= (fabs(xpull) < theMaxAllowedHitPull);
520 
521  // bailing out
522  if (!useThisHit) return false;
523 
524  // << derivs.num_row() << " x " << derivs.num_col() << std::endl;
525 
526  // AlgebraicMatrix jtvjtmp();
527  AlgebraicMatrix covtmp(covmat);
528  AlgebraicMatrix jtvjtmp(derivs * covtmp *derivs.T());
529  AlgebraicSymMatrix thisjtvj(npar);
530  AlgebraicVector thisjtve(npar);
531  // thisjtvj = covmat.similarity(derivs);
532  thisjtvj.assign(jtvjtmp);
533  thisjtve=derivs * covmat * (pos-coor);
534 
535  AlgebraicVector hitresidual(hitDim);
536  hitresidual[0] = (pos[0] - coor[0]);
537 
538  AlgebraicMatrix hitresidualT;
539  hitresidualT = hitresidual.T();
540 
541  // access user variables (via AlignmentParameters)
542  HIPUserVariables* uservar = dynamic_cast<HIPUserVariables*>(params->userVariables());
543 
544  uservar->jtvj += hitwt*thisjtvj;
545  uservar->jtve += hitwt*thisjtve;
546  uservar->nhit++;
547  // The following variable is needed for the extended 1D/2D hit fix using
548  // matrix shrinkage and expansion
549  // uservar->hitdim = hitDim;
550 
551  //for alignable chi squared
552  float thischi2;
553  thischi2 = hitwt*(hitresidualT *covmat *hitresidual)[0];
554 
555  if ( verbose && ((thischi2/ static_cast<float>(uservar->nhit)) >10.0) ) {
556  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::processHit1D]Added to Chi2 the number " << thischi2 <<" having "
557  << uservar->nhit << " dof " << std::endl << "X-resid "
558  << hitresidual[0] << " Y-resid "
559  << hitresidual[1] << std::endl << " Cov^-1 matr (covmat): [0][0]= "
560  << covmat[0][0] << " [0][1]= "
561  << covmat[0][1] << " [1][0]= "
562  << covmat[1][0] << " [1][1]= "
563  << covmat[1][1] << std::endl;
564  }
565 
566  uservar->alichi2 += thischi2; // a bit weird, but vector.transposed * matrix * vector doesn't give a double in CMSSW's opinion
567  uservar->alindof += hitDim; // 2D hits contribute twice to the ndofs
568 
569  return true;
570 }
571 
573  const Alignable* ali,
574  const TrajectoryStateOnSurface & tsos,
576  double hitwt)
577 {
578  static const unsigned int hitDim = 2;
579 
580  // get trajectory impact point
581  LocalPoint alvec = tsos.localPosition();
582  AlgebraicVector pos(hitDim);
583  pos[0] = alvec.x();
584  pos[1] = alvec.y();
585 
586 // CY: get trajectory impact angle
587 // LocalVector v = tsos.localDirection();
588 // double proj_z = v.dot(LocalVector(0,0,1));
589 // edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::processHit2D]detectorId="<<ali->id();
590 // edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::processHit2D]proj_Z="<<proj_z;
591 
592  // get impact point covariance
593  AlgebraicSymMatrix ipcovmat(hitDim);
594  ipcovmat[0][0] = tsos.localError().positionError().xx();
595  ipcovmat[1][1] = tsos.localError().positionError().yy();
596  ipcovmat[0][1] = tsos.localError().positionError().xy();
597 
598  // get hit local position and covariance
599  AlgebraicVector coor(hitDim);
600  coor[0] = hit->localPosition().x();
601  coor[1] = hit->localPosition().y();
602 
603  AlgebraicSymMatrix covmat(hitDim);
604  covmat[0][0] = hit->localPositionError().xx();
605  covmat[1][1] = hit->localPositionError().yy();
606  covmat[0][1] = hit->localPositionError().xy();
607 
608  // add hit and impact point covariance matrices
609  covmat = covmat + ipcovmat;
610 
611  // calculate the x pull and y pull of this hit
612  double xpull = 0.;
613  double ypull = 0.;
614  if (covmat[0][0] != 0.) xpull = (pos[0] - coor[0])/sqrt(fabs(covmat[0][0]));
615  if (covmat[1][1] != 0.) ypull = (pos[1] - coor[1])/sqrt(fabs(covmat[1][1]));
616 
617  // get Alignment Parameters
618  AlignmentParameters* params = ali->alignmentParameters();
619  // get derivatives
620  AlgebraicMatrix derivs2D = params->selectedDerivatives(tsos, alidet);
621  // calculate user parameters
622  int npar = derivs2D.num_row();
623  AlgebraicMatrix derivs(npar, hitDim, 0);
624 
625  for (int ipar=0;ipar<npar;ipar++) {
626  for (unsigned int idim=0;idim<hitDim;idim++) {
627  derivs[ipar][idim] = derivs2D[ipar][idim];
628  }
629  }
630  // invert covariance matrix
631  int ierr;
632  // if (covmat[1][1]>4.0) nhitDim = hitDim;
633 
634  covmat.invert(ierr);
635  if (ierr != 0) {
636  edm::LogError("Alignment") << "Matrix inversion failed!";
637  return false;
638  }
639 
640  bool useThisHit = (theMaxAllowedHitPull <= 0.);
641 
642  useThisHit |= (fabs(xpull) < theMaxAllowedHitPull && fabs(ypull) < theMaxAllowedHitPull);
643 
644  // bailing out
645  if (!useThisHit) return false;
646 
647  // << derivs.num_row() << " x " << derivs.num_col() << std::endl;
648 
649  // AlgebraicMatrix jtvjtmp();
650  AlgebraicMatrix covtmp(covmat);
651  AlgebraicMatrix jtvjtmp(derivs * covtmp *derivs.T());
652  AlgebraicSymMatrix thisjtvj(npar);
653  AlgebraicVector thisjtve(npar);
654  // thisjtvj = covmat.similarity(derivs);
655  thisjtvj.assign(jtvjtmp);
656  thisjtve=derivs * covmat * (pos-coor);
657 
658  AlgebraicVector hitresidual(hitDim);
659  hitresidual[0] = (pos[0] - coor[0]);
660  hitresidual[1] = (pos[1] - coor[1]);
661  // if(nhitDim>1) {
662  // hitresidual[1] =0.0;
663  // nhitDim=1;
664  // }
665 
666  AlgebraicMatrix hitresidualT;
667  hitresidualT = hitresidual.T();
668  // access user variables (via AlignmentParameters)
669  HIPUserVariables* uservar = dynamic_cast<HIPUserVariables*>(params->userVariables());
670 
671  uservar->jtvj += hitwt*thisjtvj;
672  uservar->jtve += hitwt*thisjtve;
673  uservar->nhit++;
674  // The following variable is needed for the extended 1D/2D hit fix using
675  // matrix shrinkage and expansion
676  // uservar->hitdim = hitDim;
677 
678  //for alignable chi squared
679  float thischi2;
680  thischi2 = hitwt*(hitresidualT *covmat *hitresidual)[0];
681 
682  if ( verbose && ((thischi2/ static_cast<float>(uservar->nhit)) >10.0) ) {
683  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::processHit2D]Added to Chi2 the number " << thischi2 <<" having "
684  << uservar->nhit << " dof " << std::endl << "X-resid "
685  << hitresidual[0] << " Y-resid "
686  << hitresidual[1] << std::endl << " Cov^-1 matr (covmat): [0][0]= "
687  << covmat[0][0] << " [0][1]= "
688  << covmat[0][1] << " [1][0]= "
689  << covmat[1][0] << " [1][1]= "
690  << covmat[1][1] << std::endl;
691  }
692 
693  uservar->alichi2 += thischi2; // a bit weird, but vector.transposed * matrix * vector doesn't give a double in CMSSW's opinion
694  uservar->alindof += hitDim; // 2D hits contribute twice to the ndofs
695 
696  return true;
697 }
698 
699 // Run the algorithm on trajectories and tracks -------------------------------
701 {
702  if (isCollector) return;
703 
704  TrajectoryStateCombiner tsoscomb;
705 
706  // AM: not really needed
707  // AM: m_Ntracks = 0 should be sufficient
708  int itr=0;
709  m_Ntracks=0;
710 //CY : hit info
711  m_sinTheta =0;
712  m_angle = 0;
713  m_detId =0;
714  m_hitwt=1;
715 
716  // AM: what is this needed for?
717  //theFile->cd();
718 
719  // loop over tracks
721  for (ConstTrajTrackPairCollection::const_iterator it=tracks.begin();
722  it!=tracks.end();
723  ++it) {
724 
725 //CY: pre-selection
726 
727  const Trajectory* traj = (*it).first;
728  const reco::Track* track = (*it).second;
729 
730  float pt = track->pt();
731  float eta = track->eta();
732  float phi = track->phi();
733  float p = track->p();
734  float chi2n = track->normalizedChi2();
735  int nhit = track->numberOfValidHits();
736  float d0 = track->d0();
737  float dz = track->dz();
738 
739  int nhpxb = track->hitPattern().numberOfValidPixelBarrelHits();
740  int nhpxf = track->hitPattern().numberOfValidPixelEndcapHits();
741  int nhtib = track->hitPattern().numberOfValidStripTIBHits();
742  int nhtob = track->hitPattern().numberOfValidStripTOBHits();
743  int nhtid = track->hitPattern().numberOfValidStripTIDHits();
744  int nhtec = track->hitPattern().numberOfValidStripTECHits();
745 
746  if (verbose) edm::LogInfo("Alignment") << "New track pt,eta,phi,chi2n,hits: "
747  << pt << ","
748  << eta << ","
749  << phi << ","
750  << chi2n << ","
751  << nhit;
752 
753 //CY: Pre-selection
754 
755  double ihitwt = 1;
756  double trkwt = 1;
757  //eta distribution from 2015 RunD, need to change formula for other runs
758  TFormula *my_formula = new TFormula("formula","2.51469/(2.51469+4.11684*x-16.7847*pow(x,2)+46.1574*pow(x,3)-55.22*pow(x,4)+29.5591*pow(x,5)-5.39816*pow(x,6))");
759  if(uniEta){
760  trkwt = Scale*(my_formula->Eval(fabs(eta)));}
761  else trkwt=Scale;
762 
763  if (trackPs){
764  double r = gRandom->Rndm();
765  if (trkwt < r){
766  //edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::run]skip event, eta ="<<eta;
767  continue ;}
768  }
769  else if (trackWt){ihitwt=trkwt;}
770 
771 
772  //edm::LogWarning("Alignment") << "UsingReweighting="<<trackWt<<",trkwt="<<trkwt<<",hitWt="<<ihitwt;
773 
774  // fill track parameters in root tree
775  if (itr<MAXREC) {
776  m_Nhits[itr]=nhit;
777  m_Pt[itr]=pt;
778  m_P[itr]=p;
779  m_Eta[itr]=eta;
780  m_Phi[itr]=phi;
781  m_Chi2n[itr]=chi2n;
782  m_nhPXB[itr]=nhpxb;
783  m_nhPXF[itr]=nhpxf;
784  m_nhTIB[itr]=nhtib;
785  m_nhTOB[itr]=nhtob;
786  m_nhTID[itr]=nhtid;
787  m_nhTEC[itr]=nhtec;
788  m_d0[itr]=d0;
789  m_dz[itr]=dz;
790  m_wt[itr]=ihitwt;
791  itr++;
792  m_Ntracks=itr;
793  }
794  // AM: Can be simplified
795 
796  std::vector<const TransientTrackingRecHit*> hitvec;
797  std::vector<TrajectoryStateOnSurface> tsosvec;
798 
799  // loop over measurements
800  std::vector<TrajectoryMeasurement> measurements = traj->measurements();
801  for (std::vector<TrajectoryMeasurement>::iterator im=measurements.begin();
802  im!=measurements.end();
803  ++im) {
804 
805  TrajectoryMeasurement meas = *im;
806 
807  // const TransientTrackingRecHit* ttrhit = &(*meas.recHit());
808  // const TrackingRecHit *hit = ttrhit->hit();
809  const TransientTrackingRecHit* hit = &(*meas.recHit());
810 
812 
813  // this is the updated state (including the current hit)
814  // TrajectoryStateOnSurface tsos=meas.updatedState();
815  // combine fwd and bwd predicted state to get state
816  // which excludes current hit
817 
819  if (eventInfo.clusterValueMap()) {
820  // check from the PrescalingMap if the hit was taken.
821  // If not skip to the next TM
822  // bool hitTaken=false;
823  AlignmentClusterFlag myflag;
824 
825  int subDet = hit->geographicalId().subdetId();
826  //take the actual RecHit out of the Transient one
827  const TrackingRecHit *rechit=hit->hit();
828  if (subDet>2) { // AM: if possible use enum instead of hard-coded value
829  const std::type_info &type = typeid(*rechit);
830 
831  if (type == typeid(SiStripRecHit1D)) {
832 
833  const SiStripRecHit1D* stripHit1D = dynamic_cast<const SiStripRecHit1D*>(rechit);
834  if (stripHit1D) {
835  SiStripRecHit1D::ClusterRef stripclust(stripHit1D->cluster());
836  // myflag=PrescMap[stripclust];
837  myflag = (*eventInfo.clusterValueMap())[stripclust];
838  } else {
839  edm::LogError("HIPAlignmentAlgorithm")
840  << "ERROR in <HIPAlignmentAlgorithm::run>: Dynamic cast of Strip RecHit failed! "
841  << "TypeId of the RecHit: " << className(*hit) <<std::endl;
842  }
843 
844  }//end if type = SiStripRecHit1D
845  else if(type == typeid(SiStripRecHit2D)){
846 
847  const SiStripRecHit2D* stripHit2D = dynamic_cast<const SiStripRecHit2D*>(rechit);
848  if (stripHit2D) {
849  SiStripRecHit2D::ClusterRef stripclust(stripHit2D->cluster());
850  // myflag=PrescMap[stripclust];
851  myflag = (*eventInfo.clusterValueMap())[stripclust];
852  } else {
853  edm::LogError("HIPAlignmentAlgorithm")
854  << "ERROR in <HIPAlignmentAlgorithm::run>: Dynamic cast of Strip RecHit failed! "
855  // << "TypeId of the TTRH: " << className(*ttrhit) << std::endl;
856  << "TypeId of the TTRH: " << className(*hit) << std::endl;
857  }
858  } //end if type == SiStripRecHit2D
859  } //end if hit from strips
860  else {
861  const SiPixelRecHit* pixelhit= dynamic_cast<const SiPixelRecHit*>(rechit);
862  if (pixelhit) {
863  SiPixelClusterRefNew pixelclust(pixelhit->cluster());
864  // myflag=PrescMap[pixelclust];
865  myflag = (*eventInfo.clusterValueMap())[pixelclust];
866  }
867  else {
868  edm::LogError("HIPAlignmentAlgorithm")
869  << "ERROR in <HIPAlignmentAlgorithm::run>: Dynamic cast of Pixel RecHit failed! "
870  // << "TypeId of the TTRH: " << className(*ttrhit) << std::endl;
871  << "TypeId of the TTRH: " << className(*hit) << std::endl;
872  }
873  } //end 'else' it is a pixel hit
874  // bool hitTaken=myflag.isTaken();
875  if (!myflag.isTaken()) {
876  continue;
877  }
878  }//end if Prescaled Hits
880 
882  meas.backwardPredictedState());
883 
884  if(tsos.isValid()){
885  // hitvec.push_back(ttrhit);
886  hitvec.push_back(hit);
887  tsosvec.push_back(tsos);
888  }
889 
890  } //hit valid
891  }
892 
893  // transform RecHit vector to AlignableDet vector
894  std::vector <AlignableDetOrUnitPtr> alidetvec = theAlignableDetAccessor->alignablesFromHits(hitvec);
895 
896  // get concatenated alignment parameters for list of alignables
898 
899  std::vector<TrajectoryStateOnSurface>::const_iterator itsos=tsosvec.begin();
900  std::vector<const TransientTrackingRecHit*>::const_iterator ihit=hitvec.begin();
901 
902  // loop over vectors(hit,tsos)
903  while (itsos != tsosvec.end()) {
904 
905  // get AlignableDet for this hit
906  const GeomDet* det = (*ihit)->det();
907  // int subDet= (*ihit)->geographicalId().subdetId();
908  uint32_t nhitDim = (*ihit)->dimension();
909 
911 
912  // get relevant Alignable
913  Alignable* ali = aap.alignableFromAlignableDet(alidet);
914 
915  if (ali!=0) {
916 
917  const TrajectoryStateOnSurface & tsos=*itsos;
918 
919 // LocalVector v = tsos.localDirection();
920 // double proj_z = v.dot(LocalVector(0,0,1));
921 
922 //In fact, sin_theta=Abs(mom_z)
923  double mom_x = tsos.localDirection().x();
924  double mom_y = tsos.localDirection().y();
925  double mom_z = tsos.localDirection().z();
926  double sin_theta = TMath::Abs(mom_z) / sqrt(pow(mom_x,2)+pow(mom_y,2)+pow(mom_z,2) );
927  double angle = TMath::ASin(sin_theta);
928 
929 
930 //Make cut on hit impact angle, reduce collision hits perpendicular to modules
931  if(IsCollision)
932  { if (angle>col_cut)ihitwt=0;}
933  else
934  { if (angle<cos_cut)ihitwt=0;}
935  m_angle = angle;
936  m_sinTheta = sin_theta;
937  m_detId = ali->id();
938  m_hitwt = ihitwt;
939 
940  if (theFillTrackMonitoring) hitTree->Fill();
941 
942  if ((nhitDim==1)&&(ihitwt!=0)) {
943  processHit1D(alidet, ali, *itsos, *ihit, ihitwt);
944  } else if ((nhitDim==2)&&(ihitwt!=0)) {
945  processHit2D(alidet, ali, *itsos, *ihit, ihitwt);
946  }
947  }
948 
949  itsos++;
950  ihit++;
951  }
952  } // end of track loop
953 
954  // fill eventwise root tree (with prescale defined in pset)
957  //edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::run] theCurrentPrescale="<<theCurrentPrescale;
958  if (theCurrentPrescale<=0) {
959  theTree->Fill();
961  }
962  }
963 }
964 
965 // ----------------------------------------------------------------------------
966 
968 {
969  int result;
970 
971  std::ifstream inIterFile(filename.c_str(), std::ios::in);
972  if (!inIterFile) {
973  edm::LogError("Alignment") << "[HIPAlignmentAlgorithm::readIterationFile] ERROR! "
974  << "Unable to open Iteration file";
975  result = -1;
976  }
977  else {
978  inIterFile >> result;
979  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::readIterationFile] "
980  << "Read last iteration number from file: " << result;
981  }
982  inIterFile.close();
983 
984  return result;
985 }
986 
987 // ----------------------------------------------------------------------------
988 
990 {
991  std::ofstream outIterFile((filename.c_str()), std::ios::out);
992  if (!outIterFile) {
993  edm::LogError("Alignment") << "[HIPAlignmentAlgorithm::writeIterationFile] ERROR: Unable to write Iteration file";
994  }
995  else {
996  outIterFile << iter;
997  edm::LogWarning("Alignment") <<"[HIPAlignmentAlgorithm::writeIterationFile] writing iteration number to file: " << iter;
998  }
999  outIterFile.close();
1000 }
1001 
1002 
1003 // ----------------------------------------------------------------------------
1004 // set alignment position error
1005 
1007 {
1008 
1009 
1010  // Check if user wants to override APE
1011  if ( !theApplyAPE )
1012  {
1013  edm::LogWarning("Alignment") <<"[HIPAlignmentAlgorithm::setAlignmentPositionError] No APE applied";
1014  return; // NO APE APPLIED
1015  }
1016 
1017 
1018  edm::LogWarning("Alignment") <<"[HIPAlignmentAlgorithm::setAlignmentPositionError] Apply APE!";
1019 
1020  double apeSPar[3], apeRPar[3];
1021  for (std::vector<std::pair<std::vector<Alignable*>, std::vector<double> > >::const_iterator alipars = theAPEParameters.begin();
1022  alipars != theAPEParameters.end();
1023  ++alipars) {
1024  const std::vector<Alignable*> &alignables = alipars->first;
1025  const std::vector<double> &pars = alipars->second;
1026 
1027  apeSPar[0] = pars[0];
1028  apeSPar[1] = pars[1];
1029  apeSPar[2] = pars[2];
1030  apeRPar[0] = pars[3];
1031  apeRPar[1] = pars[4];
1032  apeRPar[2] = pars[5];
1033 
1034  double function = pars[6];
1035 
1036  // Printout for debug
1037  printf("APE: %u alignables\n", (unsigned int)alignables.size());
1038  for ( int i=0; i<21; ++i ) {
1039  double apelinstest=calcAPE(apeSPar,i,0.);
1040  double apeexpstest=calcAPE(apeSPar,i,1.);
1041  double apestepstest=calcAPE(apeSPar,i,2.);
1042  double apelinrtest=calcAPE(apeRPar,i,0.);
1043  double apeexprtest=calcAPE(apeRPar,i,1.);
1044  double apesteprtest=calcAPE(apeRPar,i,2.);
1045  printf("APE: iter slin sexp sstep rlin rexp rstep: %5d %12.5f %12.5f %12.5f %12.5f %12.5f %12.5f \n",
1046  i,apelinstest,apeexpstest,apestepstest, apelinrtest,apeexprtest,apesteprtest);
1047  }
1048 
1049  // set APE
1050  double apeshift=calcAPE(apeSPar,theIteration,function);
1051  double aperot =calcAPE(apeRPar,theIteration,function);
1052  theAlignmentParameterStore->setAlignmentPositionError( alignables, apeshift, aperot );
1053  }
1054 
1055 }
1056 
1057 // ----------------------------------------------------------------------------
1058 // calculate APE
1059 
1060 double
1061 HIPAlignmentAlgorithm::calcAPE(double* par, int iter, double function)
1062 {
1063  double diter=(double)iter;
1064 
1065  // The following floating-point equality check is safe because this
1066  // "0." and this "1." are generated by the compiler, in the very
1067  // same file. Whatever approximization scheme it uses to turn "1."
1068  // into 0.9999999999998 in the HIPAlignmentAlgorithm::initialize is
1069  // also used here. If I'm wrong, you'll get an assertion.
1070  if (function == 0.) {
1071  return std::max(par[1],par[0]+((par[1]-par[0])/par[2])*diter);
1072  }
1073  else if (function == 1.) {
1074  return std::max(0.,par[0]*(exp(-pow(diter,par[1])/par[2])));
1075  }
1076  else if (function == 2.) {
1077  int ipar2 = (int) par[2];
1078  int step = iter/ipar2;
1079  double dstep = (double) step;
1080  return std::max(0.0, par[0] - par[1]*dstep);
1081  }
1082  else assert(false); // should have been caught in the constructor
1083 }
1084 
1085 
1086 // ----------------------------------------------------------------------------
1087 // book root trees
1088 
1090 {
1091  TString tname="T1";
1092  char iterString[15];
1093  snprintf(iterString, sizeof(iterString), "%i",theIteration);
1094  tname.Append("_");
1095  tname.Append(iterString);
1096 
1097  // create ROOT files
1098  if (theFillTrackMonitoring) {
1099  theFile = new TFile(outfile.c_str(),"update");
1100  theFile->cd();
1101 
1102  // book event-wise ROOT Tree
1103 
1104  TString tname_hit="T1_hit";
1105  tname_hit.Append("_");
1106  tname_hit.Append(iterString);
1107 
1108  theTree = new TTree(tname,"Eventwise tree");
1109 
1110  //theTree->Branch("Run", &m_Run, "Run/I");
1111  //theTree->Branch("Event", &m_Event, "Event/I");
1112  theTree->Branch("Ntracks", &m_Ntracks, "Ntracks/I");
1113  theTree->Branch("Nhits", m_Nhits, "Nhits[Ntracks]/I");
1114  theTree->Branch("nhPXB", m_nhPXB, "nhPXB[Ntracks]/I");
1115  theTree->Branch("nhPXF", m_nhPXF, "nhPXF[Ntracks]/I");
1116  theTree->Branch("nhTIB", m_nhTIB, "nhTIB[Ntracks]/I");
1117  theTree->Branch("nhTOB", m_nhTOB, "nhTOB[Ntracks]/I");
1118  theTree->Branch("nhTID", m_nhTID, "nhTID[Ntracks]/I");
1119  theTree->Branch("nhTEC", m_nhTEC, "nhTEC[Ntracks]/I");
1120  theTree->Branch("Pt", m_Pt, "Pt[Ntracks]/F");
1121  theTree->Branch("P", m_P, "P[Ntracks]/F");
1122  theTree->Branch("Eta", m_Eta, "Eta[Ntracks]/F");
1123  theTree->Branch("Phi", m_Phi, "Phi[Ntracks]/F");
1124  theTree->Branch("Chi2n", m_Chi2n, "Chi2n[Ntracks]/F");
1125  theTree->Branch("d0", m_d0, "d0[Ntracks]/F");
1126  theTree->Branch("dz", m_dz, "dz[Ntracks]/F");
1127  theTree->Branch("wt", m_wt, "wt[Ntracks]/F");
1128 
1129  hitTree = new TTree(tname_hit,"Hitwise tree");
1130  hitTree->Branch("Id", &m_detId, "Id/i");
1131  hitTree->Branch("sinTheta", &m_sinTheta, "sinTheta/F");
1132  hitTree->Branch("hitImpactAngle", &m_angle, "hitImpactAngle/F");
1133  hitTree->Branch("wt", &m_hitwt,"wt/F");
1134  }
1135  // book Alignable-wise ROOT Tree
1136 
1137  theFile2 = new TFile(outfile2.c_str(),"update");
1138  theFile2->cd();
1139 
1140  theTree2 = new TTree("T2","Alignablewise tree");
1141 
1142  theTree2->Branch("Nhit", &m2_Nhit, "Nhit/I");
1143  theTree2->Branch("Type", &m2_Type, "Type/I");
1144  theTree2->Branch("Layer", &m2_Layer, "Layer/I");
1145  theTree2->Branch("Xpos", &m2_Xpos, "Xpos/F");
1146  theTree2->Branch("Ypos", &m2_Ypos, "Ypos/F");
1147  theTree2->Branch("Zpos", &m2_Zpos, "Zpos/F");
1148  theTree2->Branch("Eta", &m2_Eta, "Eta/F");
1149  theTree2->Branch("Phi", &m2_Phi, "Phi/F");
1150  theTree2->Branch("Id", &m2_Id, "Id/i");
1151  theTree2->Branch("ObjId", &m2_ObjId, "ObjId/I");
1152 
1153  // book survey-wise ROOT Tree only if survey is enabled
1154  if (theLevels.size() > 0){
1155 
1156  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::bookRoot] Survey trees booked.";
1157  theFile3 = new TFile(ssurveyfile.c_str(),"update");
1158  theFile3->cd();
1159  theTree3 = new TTree(tname, "Survey Tree");
1160  theTree3->Branch("Id", &m3_Id, "Id/i");
1161  theTree3->Branch("ObjId", &m3_ObjId, "ObjId/I");
1162  theTree3->Branch("Par", &m3_par, "Par[6]/F");
1163  }
1164 
1165  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::bookRoot] Root trees booked.";
1166 
1167 }
1168 
1169 // ----------------------------------------------------------------------------
1170 // fill alignable-wise root tree
1171 
1173 {
1174  using std::setw;
1175  theFile2->cd();
1176 
1177  int naligned=0;
1178 
1179  //Retrieve tracker topology from geometry
1180  edm::ESHandle<TrackerTopology> tTopoHandle;
1181 // iSetup.get<IdealGeometryRecord>().get(tTopoHandle);
1182  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
1183 
1184  const TrackerTopology* const tTopo = tTopoHandle.product();
1185 
1186  for (std::vector<Alignable*>::const_iterator it=theAlignables.begin();
1187  it!=theAlignables.end();
1188  ++it) {
1189  Alignable* ali = (*it);
1191 
1192  // consider only those parameters classified as 'valid'
1193  if (dap->isValid()) {
1194 
1195  // get number of hits from user variable
1196  HIPUserVariables* uservar = dynamic_cast<HIPUserVariables*>(dap->userVariables());
1197  m2_Nhit = uservar->nhit;
1198 
1199  // get type/layer
1200  std::pair<int,int> tl = theAlignmentParameterStore->typeAndLayer(ali, tTopo);
1201  m2_Type = tl.first;
1202  m2_Layer = tl.second;
1203 
1204  // get identifier (as for IO)
1205  m2_Id = ali->id();
1206  m2_ObjId = ali->alignableObjectId();
1207 
1208  // get position
1209  GlobalPoint pos = ali->surface().position();
1210  m2_Xpos = pos.x();
1211  m2_Ypos = pos.y();
1212  m2_Zpos = pos.z();
1213  m2_Eta = pos.eta();
1214  m2_Phi = pos.phi();
1215 
1216  AlgebraicVector pars = dap->parameters();
1217 
1218  if (verbose) {
1219  edm::LogVerbatim("Alignment")
1220  << "------------------------------------------------------------------------\n"
1221  << " ALIGNABLE: " << setw(6) << naligned
1222  << '\n'
1223  << "hits: " << setw(4) << m2_Nhit
1224  << " type: " << setw(4) << m2_Type
1225  << " layer: " << setw(4) << m2_Layer
1226  << " id: " << setw(4) << m2_Id
1227  << " objId: " << setw(4) << m2_ObjId
1228  << '\n'
1229  << std::fixed << std::setprecision(5)
1230  << "x,y,z: "
1231  << setw(12) << m2_Xpos
1232  << setw(12) << m2_Ypos
1233  << setw(12) << m2_Zpos
1234  << " eta,phi: "
1235  << setw(12) << m2_Eta
1236  << setw(12) << m2_Phi
1237  << '\n'
1238  << "params: "
1239  << setw(12) << pars[0]
1240  << setw(12) << pars[1]
1241  << setw(12) << pars[2]
1242  << setw(12) << pars[3]
1243  << setw(12) << pars[4]
1244  << setw(12) << pars[5];
1245  }
1246 
1247  naligned++;
1248  theTree2->Fill();
1249  }
1250  }
1251 }
1252 
1253 // ----------------------------------------------------------------------------
1254 
1255 bool HIPAlignmentAlgorithm::calcParameters(Alignable* ali , int setDet, double start, double step)
1256 {
1257  // Alignment parameters
1259  // access user variables
1260  HIPUserVariables* uservar = dynamic_cast<HIPUserVariables*>(par->userVariables());
1261  int nhit = uservar->nhit;
1262  // The following variable is needed for the extended 1D/2D hit fix using
1263  // matrix shrinkage and expansion
1264  // int hitdim = uservar->hitdim;
1265 
1266  edm::LogWarning("Alignment") << "Processing Detector"<<ali->id() ;
1267 
1268  if ( (setDet==0) && (nhit < theMinimumNumberOfHits)) {
1269  par->setValid(false);
1270  return false;
1271  }
1272 
1273  AlgebraicSymMatrix jtvj = uservar->jtvj;
1274  AlgebraicVector jtve = uservar->jtve;
1275 
1276 
1277  // these are the alignment corrections+covariance (for selected params)
1278 
1279 
1280  int npar = jtve.num_row();
1281  AlgebraicVector params(npar);
1282  AlgebraicVector paramerr(npar);
1283  AlgebraicSymMatrix cov(npar*npar);
1284 
1285  if (setDet!=0) {
1286  if (params.num_row()!=1){
1287  edm::LogError("Alignment") << "For scanning, please only turn on one parameter! check common_cff_py.txt";
1288  return false;
1289  }
1290  if (theIteration==1) params[0] = start;
1291  else params[0]=step;
1292  }
1293 
1294 // edm::LogWarning("Alignment") << "parameters " << params;
1295 
1296  // errors of parameters
1297 
1298  if (setDet==0){
1299  int ierr;
1300  AlgebraicSymMatrix jtvjinv = jtvj.inverse(ierr);
1301 
1302  if (ierr !=0) {
1303  edm::LogError("Alignment") << "Matrix inversion failed!";
1304  return false;
1305  }
1306  params = - (jtvjinv * jtve);
1307  cov = jtvjinv;
1308 
1309 // edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::CalcParameters]: parameters before RelErrCut= " << params;
1310 // edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::CalcParameters]: cov= " << cov;
1311  AlgebraicVector relerr(npar);
1312  for (int i=0;i<npar;i++) {
1313  if (fabs(cov[i][i])>0) paramerr[i] = sqrt(fabs(cov[i][i]));
1314  else paramerr[i] = params[i];
1315  if (params[i]!=0) relerr[i] = fabs(paramerr[i]/params[i]);
1316  else relerr[i]=0;
1317  if (relerr[i] >= theMaxRelParameterError) {
1318 // edm::LogWarning("Alignment") << "RelError="<<relerr[i]<<"too large!" ;
1319  params[i] = 0;
1320  paramerr[i]=0;
1321  }
1322  }
1323  }
1324 
1325  if (setDet!=0)
1326  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::CalcParameters]: parameters = " << params;
1327 
1328  uservar->alipar=params;
1329  uservar->alierr=paramerr;
1330 
1331  AlignmentParameters* parnew = par->cloneFromSelected(params,cov);
1332  ali->setAlignmentParameters(parnew);
1333  parnew->setValid(true);
1334 
1335  return true;
1336 }
1337 
1338 //-----------------------------------------------------------------------------
1339 
1341 {
1342  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::collector] called for iteration "
1343  << theIteration << std::endl;
1344 
1345  HIPUserVariablesIORoot HIPIO;
1346 
1347 
1348  for (int ijob=1;ijob<=theCollectorNJobs;ijob++) {
1349 
1350  edm::LogWarning("Alignment") << "reading uservar for job " << ijob;
1351 
1352  std::stringstream ss;
1353  std::string str;
1354  ss << ijob;
1355  ss >> str;
1356  std::string uvfile = theCollectorPath+"/job"+str+"/IOUserVariables.root";
1357 
1358  std::vector<AlignmentUserVariables*> uvarvec =
1359  HIPIO.readHIPUserVariables(theAlignables, uvfile.c_str(),
1360  theIteration, ioerr);
1361 
1362  if (ioerr!=0) {
1363  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::collector] could not read user variable files for job "
1364  << ijob;
1365  continue;
1366  }
1367 
1368  // add
1369  std::vector<AlignmentUserVariables*> uvarvecadd;
1370  std::vector<AlignmentUserVariables*>::const_iterator iuvarnew=uvarvec.begin();
1371  for (std::vector<Alignable*>::const_iterator it=theAlignables.begin();
1372  it!=theAlignables.end();
1373  ++it) {
1374  Alignable* ali = *it;
1376 
1377  HIPUserVariables* uvarold = dynamic_cast<HIPUserVariables*>(ap->userVariables());
1378  HIPUserVariables* uvarnew = dynamic_cast<HIPUserVariables*>(*iuvarnew);
1379 
1380  HIPUserVariables* uvar = uvarold->clone();
1381  if (uvarnew!=0) {
1382  //edm::LogWarning("Alignment") << "[collector-job"<<ijob<<"]alignables:old_nhit:"<<(uvarold->nhit)<<" new_nhit:"<<(uvarnew->nhit);
1383 
1384  uvar->nhit = (uvarold->nhit)+(uvarnew->nhit);
1385  uvar->jtvj = (uvarold->jtvj)+(uvarnew->jtvj);
1386  uvar->jtve = (uvarold->jtve)+(uvarnew->jtve);
1387  uvar->alichi2 = (uvarold->alichi2)+(uvarnew->alichi2);
1388  uvar->alindof = (uvarold->alindof)+(uvarnew->alindof);
1389  delete uvarnew;
1390  }
1391 
1392  uvarvecadd.push_back(uvar);
1393  iuvarnew++;
1394  }
1395 
1397 
1398  //fill Eventwise Tree
1399  if (theFillTrackMonitoring) {
1400  uvfile = theCollectorPath+"/job"+str+"/HIPAlignmentEvents.root";
1401  edm::LogWarning("Alignment") << "Added to the tree "
1402  << fillEventwiseTree(uvfile.c_str(), theIteration, ioerr)
1403  << "tracks";
1404  }
1405 
1406  }//end loop on jobs
1407 }
1408 
1409 //------------------------------------------------------------------------------------
1410 
1411 int HIPAlignmentAlgorithm::fillEventwiseTree(const char* filename, int iter, int ierr)
1412 {
1413  int totntrk = 0;
1414  char treeName[64];
1415  snprintf(treeName, sizeof(treeName), "T1_%d", iter);
1416 //CY: hit tree
1417  char hitTreeName[64];
1418  snprintf(hitTreeName, sizeof(hitTreeName), "T1_hit_%d", iter);
1419 
1420 
1421  //open the file "HIPAlignmentEvents.root" in the job directory
1422  TFile *jobfile = new TFile(filename, "READ");
1423  //grab the tree corresponding to this iteration
1424  TTree *jobtree = (TTree*)jobfile->Get(treeName);
1425  TTree *hittree = (TTree*)jobfile->Get(hitTreeName);
1426  //address and read the variables
1427  static const int nmaxtrackperevent = 1000;
1428  int jobNtracks, jobNhitspertrack[nmaxtrackperevent], jobnhPXB[nmaxtrackperevent], jobnhPXF[nmaxtrackperevent],jobnhTIB[nmaxtrackperevent], jobnhTOB[nmaxtrackperevent],jobnhTID[nmaxtrackperevent], jobnhTEC[nmaxtrackperevent];
1429  float jobP[nmaxtrackperevent], jobPt[nmaxtrackperevent], jobEta[nmaxtrackperevent] , jobPhi[nmaxtrackperevent];
1430  float jobd0[nmaxtrackperevent],jobwt[nmaxtrackperevent], jobdz[nmaxtrackperevent] , jobChi2n[nmaxtrackperevent];
1431  float jobsinTheta,jobHitWt,jobangle;
1432  align::ID jobDetId;
1433 
1434  jobtree->SetBranchAddress("Ntracks", &jobNtracks);
1435  jobtree->SetBranchAddress("Nhits", jobNhitspertrack);
1436  jobtree->SetBranchAddress("nhPXB", jobnhPXB);
1437  jobtree->SetBranchAddress("nhPXF", jobnhPXF);
1438  jobtree->SetBranchAddress("nhTIB", jobnhTIB);
1439  jobtree->SetBranchAddress("nhTOB", jobnhTOB);
1440  jobtree->SetBranchAddress("nhTID", jobnhTID);
1441  jobtree->SetBranchAddress("nhTEC", jobnhTEC);
1442  jobtree->SetBranchAddress("Pt", jobPt);
1443  jobtree->SetBranchAddress("P", jobP);
1444  jobtree->SetBranchAddress("d0", jobd0);
1445  jobtree->SetBranchAddress("dz", jobdz);
1446  jobtree->SetBranchAddress("Eta", jobEta);
1447  jobtree->SetBranchAddress("Phi", jobPhi);
1448  jobtree->SetBranchAddress("Chi2n", jobChi2n);
1449  jobtree->SetBranchAddress("wt", jobwt);
1450 
1451 // CY: hit info
1452  hittree->SetBranchAddress("sinTheta", &jobsinTheta);
1453  hittree->SetBranchAddress("HitImpactAngle", &jobangle);
1454  hittree->SetBranchAddress("Id", &jobDetId);
1455  hittree->SetBranchAddress("wt", &jobHitWt);
1456 
1457  int ievent = 0;
1458  for (ievent=0;ievent<jobtree->GetEntries();++ievent) {
1459  jobtree->GetEntry(ievent);
1460 
1461  //fill the collector tree with them
1462 
1463  // TO BE IMPLEMENTED: a prescale factor like in run()
1464  m_Ntracks = jobNtracks;
1465  int ntrk = 0;
1466  while (ntrk<m_Ntracks) {
1467  if (ntrk<MAXREC) {
1468  totntrk = ntrk+1;
1469  m_Nhits[ntrk] = jobNhitspertrack[ntrk];
1470  m_Pt[ntrk] = jobPt[ntrk];
1471  m_P[ntrk] = jobP[ntrk];
1472  m_nhPXB[ntrk] = jobnhPXB[ntrk];
1473  m_nhPXF[ntrk] = jobnhPXF[ntrk];
1474  m_nhTIB[ntrk] = jobnhTIB[ntrk];
1475  m_nhTOB[ntrk] = jobnhTOB[ntrk];
1476  m_nhTID[ntrk] = jobnhTID[ntrk];
1477  m_nhTEC[ntrk] = jobnhTEC[ntrk];
1478  m_Eta[ntrk] = jobEta[ntrk];
1479  m_Phi[ntrk] = jobPhi[ntrk];
1480  m_Chi2n[ntrk] = jobChi2n[ntrk];
1481  m_d0[ntrk] = jobd0[ntrk];
1482  m_dz[ntrk] = jobdz[ntrk];
1483  m_wt[ntrk] = jobwt[ntrk];
1484  }//end if j<MAXREC
1485  else{
1486  edm::LogWarning("Alignment") << "[HIPAlignmentAlgorithm::fillEventwiseTree] Number of tracks in Eventwise tree exceeds MAXREC: "
1487  << m_Ntracks << " Skipping exceeding tracks.";
1488  ntrk = m_Ntracks+1;
1489  }
1490  ++ntrk;
1491  }//end while loop
1492  theTree->Fill();
1493  }//end loop on i - entries in the job tree
1494 
1495  int ihit = 0;
1496  for (ihit=0;ihit<hittree->GetEntries();++ihit) {
1497  hittree->GetEntry(ihit);
1498  m_angle = jobangle;
1499  m_sinTheta = jobsinTheta;
1500  m_detId = jobDetId;
1501  m_hitwt=jobHitWt;
1502  hitTree->Fill();
1503  }
1504 
1505  //clean up
1506  delete jobtree;
1507  delete hittree;
1508  delete jobfile;
1509 
1510  return totntrk;
1511 }//end fillEventwiseTree
RunNumber_t run() const
Definition: EventID.h:39
ClusterRef cluster() const
Definition: start.py:1
void attachUserVariables(const align::Alignables &alivec, const std::vector< AlignmentUserVariables * > &uvarvec, int &ierr)
Attach User Variables to given alignables.
double p() const
momentum vector magnitude
Definition: TrackBase.h:610
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
align::ID id() const
Return the ID of Alignable, i.e. DetId of &#39;first&#39; component GeomDet(Unit).
Definition: Alignable.h:189
T getUntrackedParameter(std::string const &, T const &) const
float xx() const
Definition: LocalError.h:24
ConstRecHitPointer const & recHit() const
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:592
AlignableDetOrUnitPtr alignableFromGeomDet(const GeomDet *geomDet)
Returns AlignableDetOrUnitPtr corresponding to given GeomDet.
AlignmentParameterStore * theAlignmentParameterStore
std::pair< int, int > typeAndLayer(const Alignable *ali, const TrackerTopology *tTopo) const
Obtain type and layer from Alignable.
const EventID & eventID() const
Definition: IOVSyncValue.h:42
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:556
static AlignableObjectId commonObjectIdProvider(const AlignableObjectId &, const AlignableObjectId &)
TSOS combine(const TSOS &pTsos1, const TSOS &pTsos2) const
virtual AlignmentParameters * cloneFromSelected(const AlgebraicVector &par, const AlgebraicSymMatrix &cov) const =0
int fillEventwiseTree(const char *filename, int iter, int ierr)
void setAlignmentPositionError(const align::Alignables &alivec, double valshift, double valrot)
Set Alignment position error.
LocalVector localDirection() const
static char const * tname
Definition: GTSchema.h:13
uint32_t ID
Definition: Definitions.h:26
align::StructureType m2_ObjId
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
AlgebraicVector alipar
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:640
void writeIterationFile(std::string filename, int iter)
void applyParameters(void)
Obsolete: Use AlignableNavigator::alignableDetFromDetId and alignableFromAlignableDet.
AlgebraicVector jtve
Alignable * alignableFromAlignableDet(const AlignableDetOrUnitPtr &adet) const
Get relevant Alignable from AlignableDet.
HIPUserVariables * clone(void) const
void applyAlignableAbsolutePositions(const align::Alignables &alis, const AlignablePositions &newpos, int &ierr)
apply absolute positions to alignables
LocalError positionError() const
int numberOfValidStripTOBHits() const
Definition: HitPattern.h:868
AlignmentParameters * alignmentParameters() const
Get the AlignmentParameters.
Definition: Alignable.h:61
Definition: Electron.h:4
const ConstTrajTrackPairCollection & trajTrackPairs() const
bool processHit1D(const AlignableDetOrUnitPtr &alidet, const Alignable *ali, const TrajectoryStateOnSurface &tsos, const TrackingRecHit *hit, double hitwt)
define event information passed to algorithms
std::vector< unsigned > theIOVrangeSet
void writeHIPUserVariables(const Alignables &alivec, const char *filename, int iter, bool validCheck, int &ierr)
bool calcParameters(Alignable *ali, int setDet, double start, double step)
AlgebraicSymMatrix jtvj
std::vector< std::pair< std::vector< Alignable * >, std::vector< double > > > theAPEParameters
const AlgebraicVector & parameters(void) const
Get alignment parameters.
DataContainer const & measurements() const
Definition: Trajectory.h:196
AlignmentUserVariables * userVariables(void) const
Get pointer to user variables.
AlignablePositions readAlignableAbsolutePositions(const align::Alignables &alivec, const char *filename, int iter, int &ierr)
read Alignable current absolute positions
float xy() const
Definition: LocalError.h:25
void clear()
remove all selected Alignables and geometrical restrictions
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:646
double calcAPE(double *par, int iter, double function)
void startNewLoop(void)
Called at start of new loop.
int numberOfValidPixelBarrelHits() const
Definition: HitPattern.h:843
virtual StructureType alignableObjectId() const =0
Return the alignable type identifier.
std::vector< AlignmentUserVariables * > readHIPUserVariables(const Alignables &alivec, const char *filename, int iter, int &ierr)
void setAlignmentParameters(AlignmentParameters *dap)
Set the AlignmentParameters.
Definition: Alignable.cc:129
CLHEP::HepMatrix AlgebraicMatrix
HIPAlignmentAlgorithm(const edm::ParameterSet &cfg)
Constructor.
float yy() const
Definition: LocalError.h:26
void setValid(bool v)
Set validity flag.
AlgebraicVector alierr
T sqrt(T t)
Definition: SSEVec.h:18
double pt() const
track transverse momentum
Definition: TrackBase.h:616
T z() const
Definition: PV3DBase.h:64
T Abs(T a)
Definition: MathUtil.h:49
void initialize(const edm::EventSetup &setup, AlignableTracker *tracker, AlignableMuon *muon, AlignableExtras *extras, AlignmentParameterStore *store)
Call at beginning of job.
AlignableNavigator * theAlignableDetAccessor
int numberOfValidStripTIDHits() const
Definition: HitPattern.h:863
ClusterRef cluster() const
void fillRoot(const edm::EventSetup &setup)
int readIterationFile(std::string filename)
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:815
int numberOfValidStripTECHits() const
Definition: HitPattern.h:873
std::vector< double > SetScanDet
const LocalTrajectoryError & localError() const
virtual LocalPoint localPosition() const =0
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
const AlignableSurface & surface() const
Return the Surface (global position and orientation) of the object.
Definition: Alignable.h:135
TrajectoryStateOnSurface const & forwardPredictedState() const
Access to forward predicted state (from fitter or builder)
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:604
const align::Alignables & selectedAlignables() const
vector of alignables selected so far
int numSelected(void) const
Get number of selected parameters.
virtual TrackingRecHit const * hit() const
CLHEP::HepVector AlgebraicVector
void setUserVariables(AlignmentUserVariables *auv)
Set pointer to user variables.
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:445
ClusterRef cluster() const
Definition: SiPixelRecHit.h:49
const T & get() const
Definition: EventSetup.h:56
CompositeAlignmentParameters selectParameters(const std::vector< AlignableDet * > &alignabledets) const
AlgebraicSymMatrix inverseCovariance() const
Get inverse of survey covariance wrt given structure type in constructor.
bool isValid() const
std::unique_ptr< AlignableObjectId > alignableObjectId_
bool processHit2D(const AlignableDetOrUnitPtr &alidet, const Alignable *ali, const TrajectoryStateOnSurface &tsos, const TrackingRecHit *hit, double hitwt)
int numberOfValidStripTIBHits() const
Definition: HitPattern.h:858
std::vector< edm::ParameterSet > theAPEParameterSet
align::StructureType m3_ObjId
int numberOfValidPixelEndcapHits() const
Definition: HitPattern.h:848
T eta() const
Definition: PV3DBase.h:76
virtual void terminate()
Called at end of job (must be implemented in derived class)
AlgebraicVector sensorResidual() const
std::vector< AlignableDetOrUnitPtr > alignablesFromHits(const std::vector< const TransientTrackingRecHit * > &hitvec)
Returns vector AlignableDetOrUnitPtr for given vector of Hits.
const AliClusterValueMap * clusterValueMap() const
bool isValid(void) const
Get validity flag.
CLHEP::HepSymMatrix AlgebraicSymMatrix
unsigned int addSelections(const edm::ParameterSet &pSet)
const std::vector< std::string > surveyResiduals_
std::vector< AlignableAbsData > AlignablePositions
Definition: AlignableData.h:51
step
void writeAlignableAbsolutePositions(const align::Alignables &alivec, const char *filename, int iter, bool validCheck, int &ierr)
write Alignable current absolute positions
eventInfo
add run, event number and lumi section
DetId geographicalId() const
const IOVSyncValue & first() const
virtual LocalError localPositionError() const =0
bool detAndSubdetInMap(const DetId &detid) const
Given a DetId, returns true if DetIds with this detector and subdetector id are in the map (not neces...
Constructor of the full muon geometry.
Definition: AlignableMuon.h:37
T x() const
Definition: PV3DBase.h:62
const PositionType & position() const
T const * product() const
Definition: ESHandle.h:86
std::vector< align::StructureType > theLevels
std::vector< Alignable * > theAlignables
bool valid() const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
const align::Alignables & alignables(void) const
get all alignables
std::vector< ConstTrajTrackPair > ConstTrajTrackPairCollection
std::string className(const T &t)
Definition: ClassName.h:30
TrajectoryStateOnSurface const & backwardPredictedState() const
Access to backward predicted state (from smoother)
Our base class.
Definition: SiPixelRecHit.h:23
void run(const edm::EventSetup &setup, const EventInfo &eventInfo)
Run the algorithm.
virtual AlgebraicMatrix selectedDerivatives(const TrajectoryStateOnSurface &tsos, const AlignableDetOrUnitPtr &alidet) const
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11