CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AlignmentParameterStore.cc
Go to the documentation of this file.
1 
9 // This class's header should be first
11 
15 
19 
24 
25 //__________________________________________________________________________________________________
27  const edm::ParameterSet& config ) :
28  theAlignables(alis)
29 {
30  if (config.getUntrackedParameter<bool>("UseExtendedCorrelations")) {
32  (config.getParameter<edm::ParameterSet>("ExtendedCorrelationsConfig"));
33  } else {
35  }
36 
37  edm::LogInfo("Alignment") << "@SUB=AlignmentParameterStore"
38  << "Created with " << theAlignables.size() << " alignables.";
39 }
40 
41 //__________________________________________________________________________________________________
43 {
44  delete theCorrelationsStore;
45 }
46 
47 //__________________________________________________________________________________________________
49 AlignmentParameterStore::selectParameters( const std::vector<AlignableDet*>& alignabledets ) const
50 {
51  std::vector<AlignableDetOrUnitPtr> detOrUnits;
52  detOrUnits.reserve(alignabledets.size());
53 
54  std::vector<AlignableDet*>::const_iterator it, iEnd;
55  for ( it = alignabledets.begin(), iEnd = alignabledets.end(); it != iEnd; ++it)
56  detOrUnits.push_back(AlignableDetOrUnitPtr(*it));
57 
58  return this->selectParameters(detOrUnits);
59 }
60 
61 //__________________________________________________________________________________________________
63 AlignmentParameterStore::selectParameters( const std::vector<AlignableDetOrUnitPtr>& alignabledets ) const
64 {
65 
67  std::map <AlignableDetOrUnitPtr,Alignable*> alidettoalimap;
68  std::map <Alignable*,int> aliposmap;
69  std::map <Alignable*,int> alilenmap;
70  int nparam=0;
71 
72  // iterate over AlignableDet's
73  for( std::vector<AlignableDetOrUnitPtr>::const_iterator iad = alignabledets.begin();
74  iad != alignabledets.end(); ++iad )
75  {
76  Alignable* ali = alignableFromAlignableDet( *iad );
77  if ( ali )
78  {
79  alidettoalimap[ *iad ] = ali; // Add to map
80  // Check if Alignable already there, insert into vector if not
81  if ( find(alignables.begin(),alignables.end(),ali) == alignables.end() )
82  {
83  alignables.push_back(ali);
85  nparam += ap->numSelected();
86  }
87  }
88  }
89 
90  AlgebraicVector* selpar = new AlgebraicVector( nparam, 0 );
91  AlgebraicSymMatrix* selcov = new AlgebraicSymMatrix( nparam, 0 );
92 
93  // Fill in parameters and corresponding covariance matricess
94  int ipos = 1; // NOTE: .sub indices start from 1
95  align::Alignables::const_iterator it1;
96  for( it1 = alignables.begin(); it1 != alignables.end(); ++it1 )
97  {
98  AlignmentParameters* ap = (*it1)->alignmentParameters();
99  selpar->sub( ipos, ap->selectedParameters() );
100  selcov->sub( ipos, ap->selectedCovariance() );
101  int npar = ap->numSelected();
102  aliposmap[*it1]=ipos;
103  alilenmap[*it1]=npar;
104  ipos +=npar;
105  }
106 
107  // Fill in the correlations. Has to be an extra loop, because the
108  // AlignmentExtendedCorrelationsStore (if used) needs the
109  // alignables' covariance matrices already present.
110  ipos = 1;
111  for( it1 = alignables.begin(); it1 != alignables.end(); ++it1 )
112  {
113  int jpos=1;
114 
115  // Look for correlations between alignables
116  align::Alignables::const_iterator it2;
117  for( it2 = alignables.begin(); it2 != it1; ++it2 )
118  {
119  theCorrelationsStore->correlations( *it1, *it2, *selcov, ipos-1, jpos-1 );
120  jpos += (*it2)->alignmentParameters()->numSelected();
121  }
122 
123  ipos += (*it1)->alignmentParameters()->numSelected();
124  }
125 
127  CompositeAlignmentParameters aap( data, alignables, alidettoalimap, aliposmap, alilenmap );
128 
129  return aap;
130 }
131 
132 
133 //__________________________________________________________________________________________________
136 {
137 
138  align::Alignables selectedAlignables;
139  std::map <AlignableDetOrUnitPtr,Alignable*> alidettoalimap; // This map won't be filled!!!
140  std::map <Alignable*,int> aliposmap;
141  std::map <Alignable*,int> alilenmap;
142  int nparam=0;
143 
144  // iterate over Alignable's
145  align::Alignables::const_iterator ita;
146  for ( ita = alignables.begin(); ita != alignables.end(); ++ita )
147  {
148  // Check if Alignable already there, insert into vector if not
149  if ( find(selectedAlignables.begin(), selectedAlignables.end(), *ita) == selectedAlignables.end() )
150  {
151  selectedAlignables.push_back( *ita );
152  AlignmentParameters* ap = (*ita)->alignmentParameters();
153  nparam += ap->numSelected();
154  }
155  }
156 
157  AlgebraicVector* selpar = new AlgebraicVector( nparam, 0 );
158  AlgebraicSymMatrix* selcov = new AlgebraicSymMatrix( nparam, 0 );
159 
160  // Fill in parameters and corresponding covariance matrices
161  int ipos = 1; // NOTE: .sub indices start from 1
162  align::Alignables::const_iterator it1;
163  for( it1 = selectedAlignables.begin(); it1 != selectedAlignables.end(); ++it1 )
164  {
165  AlignmentParameters* ap = (*it1)->alignmentParameters();
166  selpar->sub( ipos, ap->selectedParameters() );
167  selcov->sub( ipos, ap->selectedCovariance() );
168  int npar = ap->numSelected();
169  aliposmap[*it1]=ipos;
170  alilenmap[*it1]=npar;
171  ipos +=npar;
172  }
173 
174  // Fill in the correlations. Has to be an extra loop, because the
175  // AlignmentExtendedCorrelationsStore (if used) needs the
176  // alignables' covariance matrices already present.
177  ipos = 1;
178  for( it1 = selectedAlignables.begin(); it1 != selectedAlignables.end(); ++it1 )
179  {
180  int jpos=1;
181 
182  // Look for correlations between alignables
183  align::Alignables::const_iterator it2;
184  for( it2 = selectedAlignables.begin(); it2 != it1; ++it2 )
185  {
186  theCorrelationsStore->correlations( *it1, *it2, *selcov, ipos-1, jpos-1 );
187  jpos += (*it2)->alignmentParameters()->numSelected();
188  }
189 
190  ipos += (*it1)->alignmentParameters()->numSelected();
191  }
192 
194  CompositeAlignmentParameters aap( data, selectedAlignables, alidettoalimap, aliposmap, alilenmap );
195 
196  return aap;
197 }
198 
199 
200 //__________________________________________________________________________________________________
202 {
203 
204  align::Alignables alignables = aap.components();
205  const AlgebraicVector& parameters = aap.parameters();
206  const AlgebraicSymMatrix& covariance = aap.covariance();
207 
208  int ipos = 1; // NOTE: .sub indices start from 1
209 
210  // Loop over alignables
211  for( align::Alignables::const_iterator it=alignables.begin(); it != alignables.end(); ++it )
212  {
213  // Update parameters and local covariance
214  AlignmentParameters* ap = (*it)->alignmentParameters();
215  int nsel = ap->numSelected();
216  AlgebraicVector subvec = parameters.sub( ipos, ipos+nsel-1 );
217  AlgebraicSymMatrix subcov = covariance.sub( ipos, ipos+nsel-1 );
218  AlignmentParameters* apnew = ap->cloneFromSelected( subvec, subcov );
219  (*it)->setAlignmentParameters( apnew );
220 
221  // Now update correlations between detectors
222  if ( updateCorrelations )
223  {
224  int jpos = 1;
225  for( align::Alignables::const_iterator it2 = alignables.begin(); it2 != it; ++it2 )
226  {
227  theCorrelationsStore->setCorrelations( *it, *it2, covariance, ipos-1, jpos-1 );
228  jpos += (*it2)->alignmentParameters()->numSelected();
229  }
230  }
231 
232  ipos+=nsel;
233  }
234 
235 }
236 
237 
238 //__________________________________________________________________________________________________
240 {
242  for (align::Alignables::const_iterator iali = theAlignables.begin();
243  iali != theAlignables.end(); ++iali)
244  if ( (*iali)->alignmentParameters()->isValid() ) result.push_back(*iali);
245 
246  LogDebug("Alignment") << "@SUB=AlignmentParameterStore::validAlignables"
247  << "Valid alignables: " << result.size()
248  << "out of " << theAlignables.size();
249  return result;
250 }
251 
252 //__________________________________________________________________________________________________
254 {
255  AlignableDetOrUnitPtr alignableDet = _alignableDet;
256  Alignable *mother = alignableDet;
257  while (mother) {
258  if (mother->alignmentParameters()) return mother;
259  mother = mother->mother();
260  }
261 
262  return 0;
263 }
264 
265 //__________________________________________________________________________________________________
267 {
268  align::Alignables::const_iterator iali;
269  for ( iali = theAlignables.begin(); iali != theAlignables.end(); ++iali)
270  applyParameters( *iali );
271 }
272 
273 
274 //__________________________________________________________________________________________________
276 {
277 
278  AlignmentParameters *pars = (alignable ? alignable->alignmentParameters() : 0);
279  if (!pars) {
280  throw cms::Exception("BadAlignable")
281  << "applyParameters: provided alignable does not have alignment parameters";
282  }
283  pars->apply();
284 }
285 
286 
287 //__________________________________________________________________________________________________
289 {
290  // Erase contents of correlation map
292 
293  // Iterate over alignables in the store and reset parameters
294  align::Alignables::const_iterator iali;
295  for ( iali = theAlignables.begin(); iali != theAlignables.end(); ++iali )
296  resetParameters( *iali );
297 }
298 
299 
300 //__________________________________________________________________________________________________
302 {
303  if ( ali )
304  {
305  // Get alignment parameters for this alignable
307  if ( ap )
308  {
309  int npar=ap->numSelected();
310 
311  AlgebraicVector par(npar,0);
312  AlgebraicSymMatrix cov(npar,0);
313  AlignmentParameters* apnew = ap->cloneFromSelected(par,cov);
314  ali->setAlignmentParameters(apnew);
315  apnew->setValid(false);
316  }
317  else
318  edm::LogError("BadArgument") << "@SUB=AlignmentParameterStore::resetParameters"
319  << "alignable has no alignment parameter";
320  }
321  else
322  edm::LogError("BadArgument") << "@SUB=AlignmentParameterStore::resetParameters"
323  << "argument is NULL";
324 }
325 
326 
327 //__________________________________________________________________________________________________
329 {
330  align::Alignables::const_iterator iali;
331  for ( iali = theAlignables.begin(); iali != theAlignables.end(); ++iali)
332  (*iali)->cacheTransformation();
333 }
334 
335 
336 //__________________________________________________________________________________________________
338 {
339  align::Alignables::const_iterator iali;
340  for ( iali = theAlignables.begin(); iali != theAlignables.end(); ++iali)
341  (*iali)->restoreCachedTransformation();
342 }
343 
344 //__________________________________________________________________________________________________
346 {
347 
348  unsigned int nAlignables = theAlignables.size();
349 
350  for (unsigned int i = 0; i < nAlignables; ++i)
351  {
352  Alignable* ali = theAlignables[i];
353 
355  dynamic_cast<RigidBodyAlignmentParameters*>( ali->alignmentParameters() );
356 
357  if ( !ap )
358  throw cms::Exception("BadAlignable")
359  << "acquireRelativeParameters: "
360  << "provided alignable does not have rigid body alignment parameters";
361 
362  AlgebraicVector par( ap->size(),0 );
363  AlgebraicSymMatrix cov( ap->size(), 0 );
364 
365  // Get displacement and transform global->local
366  align::LocalVector dloc = ali->surface().toLocal( ali->displacement() );
367  par[0]=dloc.x();
368  par[1]=dloc.y();
369  par[2]=dloc.z();
370 
371  // Transform to local euler angles
372  align::EulerAngles euloc = align::toAngles( ali->surface().toLocal( ali->rotation() ) );
373  par[3]=euloc(1);
374  par[4]=euloc(2);
375  par[5]=euloc(3);
376 
377  // Clone parameters
378  RigidBodyAlignmentParameters* apnew = ap->clone(par,cov);
379 
380  ali->setAlignmentParameters(apnew);
381  }
382 }
383 
384 
385 //__________________________________________________________________________________________________
386 // Get type/layer from Alignable
387 // type: -6 -5 -4 -3 -2 -1 1 2 3 4 5 6
388 // TEC- TOB- TID- TIB- PxEC- PxBR- PxBr+ PxEC+ TIB+ TID+ TOB+ TEC+
389 // Layers start from zero
390 std::pair<int,int> AlignmentParameterStore::typeAndLayer(const Alignable* ali, const TrackerTopology* tTopo) const
391 {
392  return TrackerAlignableId().typeAndLayerFromDetId( ali->id(), tTopo );
393 }
394 
395 
396 //__________________________________________________________________________________________________
399  const AlignablePositions& newpos,
400  int& ierr )
401 {
402  unsigned int nappl=0;
403  ierr=0;
404 
405  // Iterate over list of alignables
406  for (align::Alignables::const_iterator iali = alivec.begin(); iali != alivec.end(); ++iali) {
407  Alignable* ali = *iali;
408  align::ID id = ali->id();
409  align::StructureType typeId = ali->alignableObjectId();
410 
411  // Find corresponding entry in AlignablePositions
412  bool found=false;
413  for (AlignablePositions::const_iterator ipos = newpos.begin(); ipos != newpos.end(); ++ipos) {
414  if (id == ipos->id() && typeId == ipos->objId()) {
415  if (found) {
416  edm::LogError("DuplicatePosition")
417  << "New positions for alignable found more than once!";
418  } else {
419  // New position/rotation
420  const align::PositionType& pnew = ipos->pos();
421  const align::RotationType& rnew = ipos->rot();
422  // Current position / rotation
423  const align::PositionType& pold = ali->globalPosition();
424  const align::RotationType& rold = ali->globalRotation();
425 
426  // shift needed to move from current to new position
427  align::GlobalVector posDiff = pnew - pold;
428  align::RotationType rotDiff = rold.multiplyInverse(rnew);
429  align::rectify(rotDiff); // correct for rounding errors
430  ali->move( posDiff );
431  ali->rotateInGlobalFrame( rotDiff );
432  LogDebug("NewPosition") << "moving by:" << posDiff;
433  LogDebug("NewRotation") << "rotating by:\n" << rotDiff;
434 
435  // add position error
436  // AlignmentPositionError ape(shift.x(),shift.y(),shift.z());
437  // (*iali)->addAlignmentPositionError(ape);
438  // (*iali)->addAlignmentPositionErrorFromRotation(rot);
439 
440  found=true;
441  ++nappl;
442  }
443  }
444  }
445  }
446 
447  if ( nappl< newpos.size() )
448  edm::LogError("Mismatch") << "Applied only " << nappl << " new positions"
449  << " out of " << newpos.size();
450 
451  LogDebug("NewPositions") << "Applied new positions for " << nappl
452  << " out of " << alivec.size() <<" alignables.";
453 
454 }
455 
456 
457 //__________________________________________________________________________________________________
459 applyAlignableRelativePositions( const align::Alignables& alivec, const AlignableShifts& shifts, int& ierr )
460 {
461 
462  ierr=0;
463  unsigned int nappl=0;
464  unsigned int nAlignables = alivec.size();
465 
466  for (unsigned int i = 0; i < nAlignables; ++i) {
467  Alignable* ali = alivec[i];
468 
469  align::ID id = ali->id();
470  align::StructureType typeId = ali->alignableObjectId();
471 
472  // Find corresponding entry in AlignableShifts
473  bool found = false;
474  for (AlignableShifts::const_iterator ipos = shifts.begin(); ipos != shifts.end(); ++ipos) {
475  if (id == ipos->id() && typeId == ipos->objId()) {
476  if (found) {
477  edm::LogError("DuplicatePosition")
478  << "New positions for alignable found more than once!";
479  } else {
480  ali->move( ipos->pos() );
481  ali->rotateInGlobalFrame( ipos->rot() );
482 
483  // Add position error
484  //AlignmentPositionError ape(pnew.x(),pnew.y(),pnew.z());
485  //ali->addAlignmentPositionError(ape);
486  //ali->addAlignmentPositionErrorFromRotation(rnew);
487 
488  found=true;
489  ++nappl;
490  }
491  }
492  }
493  }
494 
495  if ( nappl < shifts.size() )
496  edm::LogError("Mismatch") << "Applied only " << nappl << " new positions"
497  << " out of " << shifts.size();
498 
499  LogDebug("NewPositions") << "Applied new positions for " << nappl << " alignables.";
500 }
501 
502 
503 
504 //__________________________________________________________________________________________________
506 {
508 }
509 
510 
511 
512 //__________________________________________________________________________________________________
514  const Parameters& parvec, int& ierr )
515 {
516  int ipass = 0;
517  int ifail = 0;
518  ierr = 0;
519 
520  // Iterate over alignables
521  for ( align::Alignables::const_iterator iali = alivec.begin(); iali != alivec.end(); ++iali )
522  {
523  // Iterate over Parameters
524  bool found=false;
525  for ( Parameters::const_iterator ipar = parvec.begin(); ipar != parvec.end(); ++ipar)
526  {
527  // Get new alignment parameters
528  AlignmentParameters* ap = *ipar;
529 
530  // Check if parameters belong to alignable
531  if ( ap->alignable() == (*iali) )
532  {
533  if (!found)
534  {
535  (*iali)->setAlignmentParameters(ap);
536  ++ipass;
537  found=true;
538  }
539  else edm::LogError("Alignment") << "@SUB=AlignmentParameterStore::attachAlignmentParameters"
540  << "More than one parameters for Alignable.";
541  }
542  }
543  if (!found) ++ifail;
544  }
545  if (ifail>0) ierr=-1;
546 
547  LogDebug("attachAlignmentParameters") << " Parameters, Alignables: " << parvec.size() << ","
548  << alivec.size() << "\n pass,fail: " << ipass << ","<< ifail;
549 }
550 
551 
552 //__________________________________________________________________________________________________
554  bool overwrite, int& ierr )
555 {
556  attachCorrelations( theAlignables, cormap, overwrite, ierr );
557 }
558 
559 
560 //__________________________________________________________________________________________________
562  const Correlations& cormap,
563  bool overwrite, int& ierr )
564 {
565  ierr=0;
566  int icount=0;
567 
568  // Iterate over correlations
569  for ( Correlations::const_iterator icor = cormap.begin(); icor!=cormap.end(); ++icor )
570  {
571  AlgebraicMatrix mat=(*icor).second;
572  Alignable* ali1 = (*icor).first.first;
573  Alignable* ali2 = (*icor).first.second;
574 
575  // Check if alignables exist
576  if ( find( alivec.begin(), alivec.end(), ali1 ) != alivec.end() &&
577  find( alivec.begin(), alivec.end(), ali2 ) != alivec.end() )
578  {
579  // Check if correlations already existing between these alignables
580  if ( !theCorrelationsStore->correlationsAvailable(ali1,ali2) || (overwrite) )
581  {
582  theCorrelationsStore->setCorrelations(ali1,ali2,mat);
583  ++icount;
584  }
585  else edm::LogInfo("AlreadyExists") << "Correlation existing and not overwritten";
586  }
587  else edm::LogInfo("IgnoreCorrelation") << "Ignoring correlation with no alignables!";
588  }
589 
590  LogDebug( "attachCorrelations" ) << " Alignables,Correlations: " << alivec.size() <<","<< cormap.size()
591  << "\n applied: " << icount ;
592 
593 }
594 
595 
596 //__________________________________________________________________________________________________
599  const std::vector<AlignmentUserVariables*>& uvarvec, int& ierr )
600 {
601  ierr=0;
602 
603  LogDebug("DumpArguments") << "size of alivec: " << alivec.size()
604  << "\nsize of uvarvec: " << uvarvec.size();
605 
606  std::vector<AlignmentUserVariables*>::const_iterator iuvar=uvarvec.begin();
607 
608  for ( align::Alignables::const_iterator iali=alivec.begin(); iali!=alivec.end(); ++iali, ++iuvar )
609  {
610  AlignmentParameters* ap = (*iali)->alignmentParameters();
611  AlignmentUserVariables* uvarnew = (*iuvar);
612  ap->setUserVariables(uvarnew);
613  }
614 }
615 
616 
617 //__________________________________________________________________________________________________
619  double valshift, double valrot )
620 {
621  unsigned int nAlignables = alivec.size();
622 
623  for (unsigned int i = 0; i < nAlignables; ++i)
624  {
625  Alignable* ali = alivec[i];
626 
627  // First reset APE
628  AlignmentPositionError nulApe(0,0,0);
629  ali->setAlignmentPositionError(nulApe, true);
630 
631  // Set APE from displacement
632  AlignmentPositionError ape(valshift,valshift,valshift);
633  if ( valshift > 0. ) ali->addAlignmentPositionError(ape, true);
634  else ali->setAlignmentPositionError(ape, true);
635  // GF: Resetting and setting as above does not really make sense to me,
636  // and adding to zero or setting is the same! I'd just do
637  //ali->setAlignmentPositionError(AlignmentPositionError ape(valshift,valshift,valshift),true);
638 
639  // Set APE from rotation
641  r(1)=valrot; r(2)=valrot; r(3)=valrot;
643  }
644 
645  LogDebug("StoreAPE") << "Store APE from shift: " << valshift
646  << "\nStore APE from rotation: " << valrot;
647 }
648 
649 //__________________________________________________________________________________________________
652  std::vector<std::vector<ParameterId> > &paramIdsVecOut,
653  std::vector<std::vector<double> > &factorsVecOut,
654  bool all, double epsilon) const
655 {
656  // Weak point if all = false:
657  // Ignores constraints between non-subsequent levels in case the parameter is not considered in
658  // the intermediate level, e.g. global z for dets and layers is aligned, but not for rods!
659  if (!ali || !ali->alignmentParameters()) return false;
660 
661  const std::vector<bool> &aliSel= ali->alignmentParameters()->selector();
662  paramIdsVecOut.clear();
663  factorsVecOut.clear();
664 
665  bool firstComp = true;
666  for (align::Alignables::const_iterator iComp = aliComps.begin(), iCompE = aliComps.end();
667  iComp != iCompE; ++iComp) {
668 
669  const ParametersToParametersDerivatives p2pDerivs(**iComp, *ali);
670  if (!p2pDerivs.isOK()) {
671  throw cms::Exception("BadConfig")
672  << "AlignmentParameterStore::hierarchyConstraints"
673  << " Bad match of types of AlignmentParameters classes.\n";
674  return false;
675  }
676  const std::vector<bool> &aliCompSel = (*iComp)->alignmentParameters()->selector();
677  for (unsigned int iParMast = 0, iParMastUsed = 0; iParMast < aliSel.size(); ++iParMast) {
678  if (!all && !aliSel[iParMast]) continue;// no higher level parameter & constraint deselected
679  if (firstComp) { // fill output with empty arrays
680  paramIdsVecOut.push_back(std::vector<ParameterId>());
681  factorsVecOut.push_back(std::vector<double>());
682  }
683  for (unsigned int iParComp = 0; iParComp < aliCompSel.size(); ++iParComp) {
684  if (aliCompSel[iParComp]) {
685  const double factor = p2pDerivs(iParMast, iParComp);
686  if (fabs(factor) > epsilon) {
687  paramIdsVecOut[iParMastUsed].push_back(ParameterId(*iComp, iParComp));
688  factorsVecOut[iParMastUsed].push_back(factor);
689  }
690  }
691  }
692  ++iParMastUsed;
693  }
694  firstComp = false;
695  } // end loop on components
696 
697  return true;
698 }
#define LogDebug(id)
void attachUserVariables(const align::Alignables &alivec, const std::vector< AlignmentUserVariables * > &uvarvec, int &ierr)
Attach User Variables to given alignables.
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:180
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
dictionary parameters
Definition: Parameters.py:2
void resetParameters(void)
reset parameters, correlations, user variables
std::map< std::pair< Alignable *, Alignable * >, AlgebraicMatrix > Correlations
T y() const
Cartesian y coordinate.
std::pair< int, int > typeAndLayerFromDetId(const DetId &detId, const TrackerTopology *tTopo) const
virtual RigidBodyAlignmentParameters * clone(const AlgebraicVector &parameters, const AlgebraicSymMatrix &covMatrix) const
Clone all parameters (for update of parameters)
std::pair< int, int > typeAndLayer(const Alignable *ali, const TrackerTopology *tTopo) const
Obtain type and layer from Alignable.
std::pair< Alignable *, unsigned int > ParameterId
a single alignable parameter of an Alignable
virtual void addAlignmentPositionErrorFromRotation(const RotationType &rotation, bool propagateDown)=0
AlgebraicVector selectedParameters(void) const
Get selected parameters.
void setAlignmentPositionError(const align::Alignables &alivec, double valshift, double valrot)
Set Alignment position error.
align::Alignables validAlignables(void) const
get all alignables with valid parameters
uint32_t ID
Definition: Definitions.h:26
bool isOK() const
Indicate whether able to provide the derivatives.
void attachCorrelations(const align::Alignables &alivec, const Correlations &cormap, bool overwrite, int &ierr)
Attach correlations to given alignables.
const GlobalVector & displacement() const
Return change of the global position since the creation of the object.
Definition: Alignable.h:135
void restoreCachedTransformations(void)
restore the previously cached position, rotation and other parameters
virtual ~AlignmentParameterStore()
destructor
const RotationType & globalRotation() const
Return the global orientation of the object.
Definition: Alignable.h:132
void applyParameters(void)
Obsolete: Use AlignableNavigator::alignableDetFromDetId and alignableFromAlignableDet.
virtual void move(const GlobalVector &displacement)=0
Movement with respect to the global reference frame.
const std::vector< bool > & selector(void) const
Get alignment parameter selector vector.
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
virtual void addAlignmentPositionError(const AlignmentPositionError &ape, bool propagateDown)=0
virtual bool correlationsAvailable(Alignable *ap1, Alignable *ap2) const
Check whether correlations are stored for a given pair of alignables.
void applyAlignableAbsolutePositions(const align::Alignables &alis, const AlignablePositions &newpos, int &ierr)
apply absolute positions to alignables
AlignmentParameters * alignmentParameters() const
Get the AlignmentParameters.
Definition: Alignable.h:57
void updateParameters(const CompositeAlignmentParameters &aap, bool updateCorrelations=true)
update parameters
std::vector< AlignableRelData > AlignableShifts
Definition: AlignableData.h:52
Basic3DVector< T > x() const
const AlgebraicVector & parameters() const
Get alignment parameters.
const RotationType & rotation() const
Return change of orientation since the creation of the object.
Definition: Alignable.h:138
virtual void setCorrelations(Alignable *ap1, Alignable *ap2, const AlgebraicSymMatrix &cov, int row, int col)
align::RotationType toLocal(const align::RotationType &) const
Return in local frame a rotation given in global frame.
std::vector< AlignmentParameters * > Parameters
void setAlignmentParameters(AlignmentParameters *dap)
Set the AlignmentParameters.
Definition: Alignable.cc:81
CLHEP::HepMatrix AlgebraicMatrix
void rectify(RotationType &)
Correct a rotation matrix for rounding errors.
Definition: Utilities.cc:196
void setValid(bool v)
Set validity flag.
AlignmentParameterStore(const align::Alignables &alis, const edm::ParameterSet &config)
constructor
Rot3< T > rot
bool hierarchyConstraints(const Alignable *aliMaster, const align::Alignables &aliComps, std::vector< std::vector< ParameterId > > &paramIdsVecOut, std::vector< std::vector< double > > &factorsVecOut, bool all, double epsilon) const
Components components() const
Get vector of alignable components.
void applyAlignableRelativePositions(const align::Alignables &alivec, const AlignableShifts &shifts, int &ierr)
apply relative shifts to alignables
tuple result
Definition: query.py:137
virtual StructureType alignableObjectId() const =0
Return the alignable type identifier.
const AlgebraicSymMatrix & covariance() const
Get parameter covariance matrix.
virtual void resetCorrelations(void)
Reset correlations.
(Abstract) Base class for alignment algorithm user variables
EulerAngles toAngles(const RotationType &)
Convert rotation matrix to angles about x-, y-, z-axes (frame rotation).
Definition: Utilities.cc:7
Alignable * alignable(void) const
Get pointer to corresponding alignable.
virtual void setAlignmentPositionError(const AlignmentPositionError &ape, bool propagateDown)=0
Set the alignment position error - if (!propagateDown) do not affect daughters.
const AlignableSurface & surface() const
Return the Surface (global position and orientation) of the object.
Definition: Alignable.h:126
AlgebraicSymMatrix selectedCovariance(void) const
Get covariance matrix of selected parameters.
Alignable * alignableFromAlignableDet(const AlignableDetOrUnitPtr &alignableDet) const
Obsolete: Use AlignableNavigator::alignableDetFromGeomDet and alignableFromAlignableDet.
int numSelected(void) const
Get number of selected parameters.
CLHEP::HepVector AlgebraicVector
AlgebraicVector EulerAngles
Definition: Definitions.h:36
virtual void correlations(Alignable *ap1, Alignable *ap2, AlgebraicSymMatrix &cov, int row, int col) const
void setUserVariables(AlignmentUserVariables *auv)
Set pointer to user variables.
int size(void) const
Get number of parameters.
virtual void rotateInGlobalFrame(const RotationType &rotation)=0
align::Alignables theAlignables
alignables
CompositeAlignmentParameters selectParameters(const std::vector< AlignableDet * > &alignabledets) const
std::vector< Alignable * > Alignables
Definition: Utilities.h:28
AlignmentCorrelationsStore * theCorrelationsStore
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
void cacheTransformations(void)
cache the current position, rotation and other parameters
RotationType toMatrix(const EulerAngles &)
Convert rotation angles about x-, y-, z-axes to matrix.
Definition: Utilities.cc:40
CLHEP::HepSymMatrix AlgebraicSymMatrix
std::vector< AlignableAbsData > AlignablePositions
Definition: AlignableData.h:51
const PositionType & globalPosition() const
Return the global position of the object.
Definition: Alignable.h:129
const double epsilon
Basic3DVector< T > multiplyInverse(const Basic3DVector< T > &v) const
void attachAlignmentParameters(const align::Alignables &alivec, const Parameters &parvec, int &ierr)
Attach alignment parameters to given alignables.
virtual AlignmentParameters * cloneFromSelected(const AlgebraicVector &par, const AlgebraicSymMatrix &cov) const =0
Alignable * mother() const
Return pointer to container alignable (if any)
Definition: Alignable.h:85
virtual void apply()=0
apply parameters to alignable
const align::Alignables & alignables(void) const
get all alignables