CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/src/Alignment/CommonAlignmentParametrization/src/TwoBowedSurfacesAlignmentParameters.cc

Go to the documentation of this file.
00001 
00008 #include "FWCore/Utilities/interface/Exception.h"
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 
00011 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00012 
00013 #include "Alignment/CommonAlignment/interface/Alignable.h"
00014 #include "Alignment/CommonAlignment/interface/AlignableDetOrUnitPtr.h"
00015 #include "Alignment/CommonAlignmentParametrization/interface/KarimakiAlignmentDerivatives.h"
00016 #include "Alignment/CommonAlignmentParametrization/interface/AlignmentParametersFactory.h"
00017 #include "CondFormats/Alignment/interface/Definitions.h"
00018 //#include "DataFormats/Alignment/interface/SurfaceDeformation.h"
00019 #include "Geometry/CommonTopologies/interface/TwoBowedSurfacesDeformation.h"
00020 
00021 // This class's header 
00022 #include "Alignment/CommonAlignmentParametrization/interface/TwoBowedSurfacesAlignmentParameters.h"
00023 
00024 #include <math.h>
00025 #include <iostream>
00026 //_________________________________________________________________________________________________
00027 TwoBowedSurfacesAlignmentParameters::TwoBowedSurfacesAlignmentParameters(Alignable *ali) :
00028   AlignmentParameters(ali, AlgebraicVector(N_PARAM), AlgebraicSymMatrix(N_PARAM, 0)),
00029   ySplit_(this->ySplitFromAlignable(ali))
00030 {
00031 }
00032 
00033 //_________________________________________________________________________________________________
00034 TwoBowedSurfacesAlignmentParameters
00035 ::TwoBowedSurfacesAlignmentParameters(Alignable *alignable, 
00036                                       const AlgebraicVector &parameters, 
00037                                       const AlgebraicSymMatrix &covMatrix) :
00038   AlignmentParameters(alignable, parameters, covMatrix),
00039   ySplit_(this->ySplitFromAlignable(alignable))
00040 {
00041   if (parameters.num_row() != N_PARAM) {
00042     throw cms::Exception("BadParameters") << "in TwoBowedSurfacesAlignmentParameters(): "
00043                                           << parameters.num_row() << " instead of " << N_PARAM 
00044                                           << " parameters.";
00045   }
00046 }
00047 
00048 //_________________________________________________________________________________________________
00049 TwoBowedSurfacesAlignmentParameters
00050 ::TwoBowedSurfacesAlignmentParameters(Alignable *alignable, 
00051                                       const AlgebraicVector &parameters, 
00052                                       const AlgebraicSymMatrix &covMatrix,
00053                                       const std::vector<bool> &selection) :
00054   AlignmentParameters(alignable, parameters, covMatrix, selection),
00055   ySplit_(this->ySplitFromAlignable(alignable))
00056 {
00057   if (parameters.num_row() != N_PARAM) {
00058     throw cms::Exception("BadParameters") << "in TwoBowedSurfacesAlignmentParameters(): "
00059                                           << parameters.num_row() << " instead of " << N_PARAM 
00060                                           << " parameters.";
00061   }
00062 }
00063 
00064 //_________________________________________________________________________________________________
00065 TwoBowedSurfacesAlignmentParameters* 
00066 TwoBowedSurfacesAlignmentParameters::clone(const AlgebraicVector &parameters, 
00067                                            const AlgebraicSymMatrix &covMatrix) const 
00068 {
00069   TwoBowedSurfacesAlignmentParameters* rbap = 
00070     new TwoBowedSurfacesAlignmentParameters(this->alignable(), parameters, covMatrix, selector());
00071 
00072   if (this->userVariables()) rbap->setUserVariables(this->userVariables()->clone());
00073   rbap->setValid(this->isValid());
00074 
00075   return rbap;
00076 }
00077 
00078 //_________________________________________________________________________________________________
00079 TwoBowedSurfacesAlignmentParameters* 
00080 TwoBowedSurfacesAlignmentParameters::cloneFromSelected(const AlgebraicVector &parameters,
00081                                                        const AlgebraicSymMatrix &covMatrix) const
00082 {
00083   return this->clone(this->expandVector(parameters, this->selector()),
00084                      this->expandSymMatrix(covMatrix, this->selector()));
00085 }
00086 
00087 //_________________________________________________________________________________________________
00088 AlgebraicMatrix 
00089 TwoBowedSurfacesAlignmentParameters::derivatives(const TrajectoryStateOnSurface &tsos,
00090                                                  const AlignableDetOrUnitPtr &alidet) const
00091 {
00092   const Alignable *ali = this->alignable(); // Alignable of these parameters
00093   AlgebraicMatrix result(N_PARAM, 2); // initialised with zeros
00094 
00095   if (ali == alidet) {
00096     const AlignableSurface &surf = ali->surface();
00097  
00098     // matrix of dimension BowedDerivs::N_PARAM x 2 
00099     const AlgebraicMatrix derivs(BowedDerivs()(tsos, surf.width(), surf.length(),
00100                                                true, ySplit_)); // split at ySplit_!
00101 
00102     // Parameters belong to surface part with y < ySplit_ or y >= ySplit_?
00103     const double localY = tsos.localParameters().mixedFormatVector()[4];
00104     const unsigned int indexOffset = (localY < ySplit_ ? 0 : dx2);
00105     // Copy derivatives to relevant part of result
00106     for (unsigned int i = BowedDerivs::dx; i < BowedDerivs::N_PARAM; ++i) {
00107       result[indexOffset + i][0] = derivs[i][0];
00108       result[indexOffset + i][1] = derivs[i][1];
00109     } 
00110   } else {
00111     // The following is even more difficult for TwoBowedSurfacesAlignmentParameters
00112     // than for BowedSurfaceAlignmentParameters where this text comes from:
00113     //
00114     // We could give this a meaning by applying frame-to-frame derivatives 
00115     // to the rigid body part of the parameters (be careful that alpha ~= dslopeY
00116     // and beta ~= -dslopeX, but with changed scale!)
00117     // and keep the surface structure parameters untouched in local meaning.
00118     // In this way we could do higher level alignment and determine 'average'
00119     // surface structures for the components.
00120     throw cms::Exception("MisMatch")
00121       << "TwoBowedSurfacesAlignmentParameters::derivatives: The hit alignable must match the "
00122       << "aligned one (i.e. bowed surface parameters cannot be used for composed alignables)\n";
00123   }
00124 
00125   return result;
00126 }
00127 
00128 //_________________________________________________________________________________________________
00129 void TwoBowedSurfacesAlignmentParameters::apply()
00130 {
00131   Alignable *alignable = this->alignable();
00132   if (!alignable) {
00133     throw cms::Exception("BadParameters") 
00134       << "TwoBowedSurfacesAlignmentParameters::apply: parameters without alignable";
00135   }
00136   
00137   // Some repeatedly needed variables
00138   const AlignableSurface &surface = alignable->surface();
00139   const double halfLength  = surface.length() * 0.5; // full module
00140   const double halfLength1 = (halfLength + ySplit_) * 0.5;        // low-y surface
00141   const double halfLength2 = (halfLength - ySplit_) * 0.5;        // high-y surface
00142 
00143   // first copy the parameters into separate parts for the two surfaces
00144   const AlgebraicVector &params = theData->parameters();
00145   std::vector<double> rigidBowPar1(BowedDerivs::N_PARAM); // 1st surface (y <  ySplit_)
00146   std::vector<double> rigidBowPar2(BowedDerivs::N_PARAM); // 2nd surface (y >= ySplit_)
00147   for (unsigned int i = 0; i < BowedDerivs::N_PARAM; ++i) {
00148     rigidBowPar1[i] = params[i];
00149     rigidBowPar2[i] = params[i + BowedDerivs::N_PARAM];
00150   }
00151   // Now adjust slopes to angles, note that dslopeX <-> -beta & dslopeY <-> alpha,
00152   // see BowedSurfaceAlignmentParameters::rotation(): FIXME: use atan?
00153   rigidBowPar1[3] =  params[dslopeY1] / halfLength1; // alpha1
00154   rigidBowPar2[3] =  params[dslopeY2] / halfLength2; // alpha2 
00155   rigidBowPar1[4] = -params[dslopeX1] / (surface.width() * 0.5); // beta1
00156   rigidBowPar2[4] = -params[dslopeX2] / (surface.width() * 0.5); // beta2
00157   // gamma is simply scaled
00158   const double gammaScale1 = BowedDerivs::gammaScale(surface.width(), 2.0*halfLength1);
00159   rigidBowPar1[5] = params[drotZ1] / gammaScale1;
00160 //   const double gammaScale2 = std::sqrt(halfLength2 * halfLength2
00161 //                                     + surface.width() * surface.width()/4.);
00162   const double gammaScale2 = BowedDerivs::gammaScale(surface.width(), 2.0*halfLength2);
00163   rigidBowPar2[5] = params[drotZ2] / gammaScale2;
00164 
00165   // Get rigid body rotations of full module as mean of the two surfaces:
00166   align::EulerAngles angles(3); // to become 'common' rotation in local frame
00167   for (unsigned int i = 0; i < 3; ++i) {
00168     angles[i] = (rigidBowPar1[i+3] + rigidBowPar2[i+3]) * 0.5;
00169   }
00170   // Module rotations are around other axes than the one we determined,
00171   // so we have to correct that the surfaces are shifted by the rotation around
00172   // the module axis - in linear approximation just an additional shift:
00173   const double yMean1 = -halfLength + halfLength1;// y of alpha1 rotation axis in module frame
00174   const double yMean2 =  halfLength - halfLength2;// y of alpha2 rotation axis in module frame
00175   rigidBowPar1[dz1] -= angles[0] * yMean1; // correct w1 for alpha
00176   rigidBowPar2[dz1] -= angles[0] * yMean2; // correct w2 for alpha
00177   // Nothing for beta1/2 since anyway both around the y-axis of the module.
00178   rigidBowPar1[dx1] += angles[2] * yMean1; // correct x1 for gamma
00179   rigidBowPar2[dx1] += angles[2] * yMean2; // correct x1 for gamma
00180 
00181   // Get rigid body shifts of full module as mean of the two surfaces:
00182   const align::LocalVector shift((rigidBowPar1[dx1] + rigidBowPar2[dx1]) * 0.5, // dx1!
00183                                  (rigidBowPar1[dy1] + rigidBowPar2[dy1]) * 0.5, // dy1!
00184                                  (rigidBowPar1[dz1] + rigidBowPar2[dz1]) * 0.5);// dz1!
00185   // Apply module shift and rotation:
00186   alignable->move(surface.toGlobal(shift));
00187   // original code:
00188   //  alignable->rotateInLocalFrame( align::toMatrix(angles) );
00189   // correct for rounding errors:
00190   align::RotationType rot(surface.toGlobal(align::toMatrix(angles)));
00191   align::rectify(rot);
00192   alignable->rotateInGlobalFrame(rot);
00193 
00194   // Fill surface structures with mean bows and half differences for all parameters:
00195   std::vector<align::Scalar> deformations; deformations.reserve(13);
00196   // first part: average bows
00197   deformations.push_back((params[dsagittaX1 ] + params[dsagittaX2 ]) * 0.5);
00198   deformations.push_back((params[dsagittaXY1] + params[dsagittaXY2]) * 0.5);
00199   deformations.push_back((params[dsagittaY1 ] + params[dsagittaY2 ]) * 0.5);
00200   // second part: half difference of all corrections
00201   for (unsigned int i = 0; i < BowedDerivs::N_PARAM; ++i) { 
00202     // sign means that we have to apply e.g.
00203     // - sagittaX for sensor 1: deformations[0] + deformations[9]
00204     // - sagittaX for sensor 2: deformations[0] - deformations[9]
00205     // - additional dx for sensor 1:  deformations[3]
00206     // - additional dx for sensor 2: -deformations[3]
00207     deformations.push_back((rigidBowPar1[i] - rigidBowPar2[i]) * 0.5);
00208   }
00209   // finally: keep track of where we have split the module
00210   deformations.push_back(ySplit_); // index is 12
00211 
00212   //  const SurfaceDeformation deform(SurfaceDeformation::kTwoBowedSurfaces, deformations);
00213   const TwoBowedSurfacesDeformation deform(deformations);
00214 
00215   // FIXME: true to propagate down?
00216   //        Needed for hierarchy with common deformation parameter,
00217   //        but that is not possible now anyway.
00218   alignable->addSurfaceDeformation(&deform, false);
00219 }
00220 
00221 //_________________________________________________________________________________________________
00222 int TwoBowedSurfacesAlignmentParameters::type() const
00223 {
00224   return AlignmentParametersFactory::kTwoBowedSurfaces;
00225 }
00226 
00227 //_________________________________________________________________________________________________
00228 void TwoBowedSurfacesAlignmentParameters::print() const
00229 {
00230 
00231   std::cout << "Contents of TwoBowedSurfacesAlignmentParameters:"
00232             << "\nParameters: " << theData->parameters()
00233             << "\nCovariance: " << theData->covariance() << std::endl;
00234 }
00235 
00236 //_________________________________________________________________________________________________
00237 double TwoBowedSurfacesAlignmentParameters::ySplitFromAlignable(const Alignable *ali) const
00238 {
00239   if (!ali) return 0.;
00240 
00241   const align::PositionType pos(ali->globalPosition());  
00242   const double r = pos.perp();
00243 
00244   // The returned numbers for TEC are calculated as stated below from
00245   // what is found in CMS-NOTE 2003/20.
00246   // Note that at that time it was planned to use ST sensors for the outer TEC
00247   // while in the end there are only a few of them in the tracker - the others are HPK.
00248   // No idea whether there are subtle changes in geometry.
00249   // The numbers have been cross checked with y-residuals in data, see 
00250   // https://hypernews.cern.ch/HyperNews/CMS/get/recoTracking/1018/1/1/2/1/1/1/2/1.html.
00251 
00252   if (r < 58.) { // Pixel, TIB, TID or TEC ring 1-4
00253     edm::LogError("Alignment") << "@SUB=TwoBowedSurfacesAlignmentParameters::ySplitFromAlignable"
00254                                << "There are no split modules for radii < 58, but r = " << r;
00255     return 0.;
00256   } else if (fabs(pos.z()) < 118.) { // TOB
00257     return 0.;
00258   } else if (r > 90.) { // TEC ring 7
00259     // W7a Height active= 106.900mm (Tab. 2) (but 106.926 mm p. 40!?)
00260     // W7a R min active = 888.400mm (Paragraph 5.5.7)
00261     // W7a R max active = W7a R min active + W7a Height active = 
00262     //                  = 888.400mm + 106.900mm = 995.300mm
00263     // W7b Height active=  94.900mm (Tab. 2)  (but 94.876 mm p. 41!?)
00264     // W7b R min active = 998.252mm (Paragraph 5.5.8)
00265     // W7b R max active = 998.252mm + 94.900mm = 1093.152mm
00266     // mean radius module = 0.5*(1093.152mm + 888.400mm) = 990.776mm
00267     // mean radius gap = 0.5*(998.252mm + 995.300mm) = 996.776mm
00268     // ySplit = (mean radius gap - mean radius module) // local y and r have same directions!
00269     //        = 996.776mm - 990.776mm = 6mm
00270     return 0.6;
00271   } else if (r > 75.) { // TEC ring 6
00272     // W6a Height active= 96.100mm (Tab. 2) (but 96.136 mm p. 38!?)
00273     // W6a R min active = 727.000mm (Paragraph 5.5.5)
00274     // W6a R max active = W6a R min active + W6a Height active = 
00275     //                  = 727.000mm + 96.100mm = 823.100mm
00276     // W6b Height active= 84.900mm (Tab. 2) (but 84.936 mm p. 39!?)
00277     // W6b R min active = 826.060mm (Paragraph 5.5.6)
00278     // W6b R max active = 826.060mm + 84.900mm = 910.960mm
00279     // mean radius module = 0.5*(910.960mm + 727.000mm) = 818.980mm
00280     // mean radius gap = 0.5*(826.060mm + 823.100mm) = 824.580mm
00281     //          -1: local y and r have opposite directions!
00282     // ySplit = -1*(mean radius gap - mean radius module)
00283     //        = -1*(824.580mm - 818.980mm) = -5.6mm
00284     return -0.56;
00285   } else {              // TEC ring 5 - smaller radii alreay excluded before
00286     // W5a Height active= 81.200mm (Tab. 2) (but 81.169 mm p. 36!?)
00287     // W5a R min active = 603.200mm (Paragraph 5.5.3)
00288     // W5a R max active = W5a R min active + W5a Height active = 
00289     //                  = 603.200mm + 81.200mm = 684.400mm
00290     // W5b Height active= 63.200mm (Tab. 2) (63.198 mm on p. 37)
00291     // W5b R min active = 687.293mm (Abschnitt 5.5.4 der note)
00292     // W5b R max active = 687.293mm + 63.200mm = 750.493mm
00293     // mean radius module = 0.5*(750.493mm + 603.200mm) = 676.8465mm
00294     // mean radius gap = 0.5*(687.293mm + 684.400mm) = 685.8465mm
00295     //          -1: local y and r have opposite directions!
00296     // ySplit = -1*(mean radius gap - mean radius module)
00297     //        = -1*(685.8465mm - 676.8465mm) = -9mm
00298     return -0.9;
00299   }
00300   //  return 0.;
00301 }