CMS 3D CMS Logo

LaserAlignmentAlgorithmPosTEC.cc

Go to the documentation of this file.
00001 
00009 #include "Alignment/LaserAlignment/interface/LaserAlignmentAlgorithmPosTEC.h"
00010 #include "Alignment/LaserAlignment/interface/Millepede.h" 
00011 #include "FWCore/ParameterSet/interface/ParameterSet.h" 
00012 
00013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00014 
00015 #include "DataFormats/TrackingRecHit/interface/AlignmentPositionError.h"
00016 
00017 LaserAlignmentAlgorithmPosTEC::LaserAlignmentAlgorithmPosTEC(edm::ParameterSet const& theConf, int theAlignmentIteration) 
00018   : theFirstFixedDiskPosTEC(theConf.getUntrackedParameter<int>("FirstFixedParameterPosTEC",2)), 
00019     theSecondFixedDiskPosTEC(theConf.getUntrackedParameter<int>("SecondFixedParameterPosTEC",3)),
00020     theGlobalParametersPosTEC(), theLocalParametersPosTEC()
00021 {
00022 
00023   // initialize Millepede
00024   initMillepede(theAlignmentIteration);
00025 }
00026 
00027 LaserAlignmentAlgorithmPosTEC::~LaserAlignmentAlgorithmPosTEC()
00028 {
00029 }
00030 
00031 void LaserAlignmentAlgorithmPosTEC::addLaserBeam(std::vector<double> theMeasurements, int LaserBeam, int LaserRing)
00032 {
00033   LogDebug("LaserAlignmentAlgorithmPosTEC") << "<LaserAlignmentAlgorithmPosTEC::addLaserBeam()>: adding a new Laser Beam ... ";
00034   
00035   float RRing = 0.0;
00036 
00037   if (LaserRing == 4) 
00038     { RRing = 56.4; }             // radius of Ring 4 
00039   else if (LaserRing == 6)
00040     { RRing = 84.0; }             // radius of Ring 6
00041 
00042   float LaserPhi0 = 0.392699082;  // phi position of the first Laserdiode in Ring 4
00043   float LaserPhi = 0.0;           // phi position of the actual beam; will be calculated
00044   int nLaserBeams = 8;            // number of beams in one ring
00045 
00046   float PhiMeasured[9] = {0.0};   // the measured phi positions
00047   float PhiErrors[9] = {0.0};     // errors of the measurements
00048 
00049   for (int i = 0; i < 9; i++)
00050     {
00051       PhiMeasured[i] = float(theMeasurements.at(2*i));
00052       PhiErrors[i] = float(theMeasurements.at((2*i)+1));
00053       LogDebug("AlignmentAlgorithmPosTEC") << " i = " << i << "\t PhiMeasured[" << i << "] = " << PhiMeasured[i]
00054                                            << "\t PhiError[" << i << "] = " << PhiErrors[i];
00055     }
00056 
00057   // calculate phi for this beam
00058   LaserPhi = LaserPhi0 + float(LaserBeam * float(float(2 * M_PI) / nLaserBeams));
00059 
00060   // loop over the discs
00061   for (int theDisk = 0; theDisk < 9; theDisk++)
00062     {
00063       theLocalParametersPosTEC[0] = 1.0;
00064       theGlobalParametersPosTEC[3*theDisk] = 1.0;                        // displacement by dphi
00065       theGlobalParametersPosTEC[(3*theDisk)+1] = -sin(LaserPhi)/RRing;  // displacement by dx
00066       theGlobalParametersPosTEC[(3*theDisk)+2] = cos(LaserPhi)/RRing;   // displacement by dy
00067 
00068       // create the equation for a single measurement
00069       equloc_(theGlobalParametersPosTEC, theLocalParametersPosTEC, &PhiMeasured[theDisk], &PhiErrors[theDisk]);
00070     }
00071   
00072   // local fit after one laser beam has been added
00073   fitloc_();
00074 }
00075 
00076 void LaserAlignmentAlgorithmPosTEC::doGlobalFit(AlignableTracker * theAlignableTracker)
00077 {
00078   edm::LogInfo("AlignmentAlgorithmPosTEC") << "<AlignmentAlgorithmPosTEC::doGlobalFit()>: do the global fit ... ";
00079   
00080   // number of global parameters
00081   const int nGlobalParameters = 27;
00082   // array to store the results of the fit
00083   float theFittedGlobalParametersPosTEC[nGlobalParameters] = { 0.0 };
00084 
00085   float RRing6 = 84.0;            // radius of Ring 6
00086   float LaserPhi0 = 0.392699082;  // phi position of the first Laserdiode in Ring 4
00087 
00088   // corrected global Phi + error
00089   std::vector<float> thePhiCorrected;
00090   std::vector<float> thePhiCorrectedError;
00091   std::vector<float> theAbsPhiCorrected;
00092   std::vector<float> theAbsPhiCorrectedError;
00093 
00094   thePhiCorrected.clear();
00095   thePhiCorrectedError.clear();
00096   theAbsPhiCorrected.clear();
00097   theAbsPhiCorrectedError.clear();
00098 
00099   // do the fit
00100   fitglo_(theFittedGlobalParametersPosTEC);
00101 
00102   int ep = 1;
00103   for (int i = 0; i < nGlobalParameters; i++)
00104     {
00105       LogDebug("AlignmentAlgorithmPosTEC") << "Global Parameter (TEC+) " << i << " = " << theFittedGlobalParametersPosTEC[i] 
00106                                            << " +/- " << errpar_(&ep);
00107       ep++;
00108     }
00109 
00110   int p0 = 0, p1 = 1, p2 = 2, p3 = 3, p4 = 4, p5 = 5;
00111   int n1 = 1, n2 = 2, n3 = 3, n4 = 4, n5 = 5, n6 =6;
00112 
00113   for (int j = 0; j < 8; ++j)
00114     {
00115       // calculate the corrections and the errors
00116       thePhiCorrected.push_back( theFittedGlobalParametersPosTEC[p0] 
00117                                 + (sin(LaserPhi0)/RRing6) * theFittedGlobalParametersPosTEC[p1]
00118                                 - (cos(LaserPhi0)/RRing6) * theFittedGlobalParametersPosTEC[p2]
00119                                 - ( theFittedGlobalParametersPosTEC[p3]
00120                                    + (sin(LaserPhi0)/RRing6) * theFittedGlobalParametersPosTEC[p4]
00121                                    - (cos(LaserPhi0)/RRing6) * theFittedGlobalParametersPosTEC[p5] ));
00122       thePhiCorrectedError.push_back(sqrt(pow(errpar_(&n1),2) 
00123                                           + pow(sin(LaserPhi0)/RRing6,2) * pow(errpar_(&n2),2) 
00124                                           + pow(cos(LaserPhi0)/RRing6,2) * pow(errpar_(&n3),2)
00125                                           + pow(errpar_(&n4),2)
00126                                           + pow(sin(LaserPhi0)/RRing6,2) * pow(errpar_(&n5),2)
00127                                           + pow(cos(LaserPhi0)/RRing6,2) * pow(errpar_(&n6),2) ) );
00128 
00129       // for debugging
00130       edm::LogInfo("LaserAlignmentAlgorithmPosTEC:Results") << " Fitted relative Correction for TEC+ in Phi[" << j << "] = " << thePhiCorrected.at(j) << " +/- "
00131                                                        << thePhiCorrectedError.at(j);
00132 
00133 
00134       p3 += 3; p4 += 3; p5 += 3;
00135       n4 += 3; n5 += 3; n6 += 3;
00136 
00137     }
00138 
00139   for (int j = 0; j < 9; ++j)
00140     {
00141       // calculate the correction for each disk (not relative to disk one)
00142       int e1 = 3*j+1;
00143       int e2 = 3*j+2;
00144       int e3 = 3*j+3;
00145       theAbsPhiCorrected.push_back(-1.0 * theFittedGlobalParametersPosTEC[j*3] 
00146                                    + (sin(LaserPhi0)/RRing6) * theFittedGlobalParametersPosTEC[j*3 + 1]
00147                                    - (cos(LaserPhi0)/RRing6) * theFittedGlobalParametersPosTEC[j*3 + 2]);
00148       theAbsPhiCorrectedError.push_back(sqrt(pow(errpar_(&e1),2) 
00149                                              + pow(sin(LaserPhi0)/RRing6,2) * pow(errpar_(&e2),2) 
00150                                              + pow(cos(LaserPhi0)/RRing6,2) * pow(errpar_(&e3),2)));
00151       
00152 
00153       // for debugging
00154       edm::LogInfo("LaserAlignmentAlgorithmPosTEC:Results") << " Fitted Correction for TEC+ in Phi[" << j << "] = " << theAbsPhiCorrected.at(j) << " +/- "
00155                                                        << theAbsPhiCorrectedError.at(j);
00156 
00157     }
00158 
00159   // loop over all discs, access the AlignableTracker to move the discs 
00160   // according to the calculated alignment corrections
00161   // AlignableTracker will take care to the propagation of the movements
00162   // to the lowest level of alignable objects
00163   const align::Alignables& endcaps = theAlignableTracker->endCaps();
00164   const Alignable* endcap = endcaps[0];
00165   if (endcap->globalPosition().z() < 0) endcap = endcaps[1];
00166   const align::Alignables& disks = endcap->components();
00167 
00168   for (unsigned int i = 0; i < disks.size(); ++i)
00169     {
00170       int aPhi = 3*i;
00171       int aX   = 3*i + 1;
00172       int aY   = 3*i + 2;
00173       int ePhi = 3*i + 1;
00174       int eX = 3*i + 2;
00175       int eY = 3*i + 3;
00176 
00177       align::GlobalVector translation(-1.0 * theFittedGlobalParametersPosTEC[aX], 
00178                                       -1.0 * theFittedGlobalParametersPosTEC[aY],
00179                                       0.0);
00180       AlignmentPositionError positionError(errpar_(&eX),errpar_(&eY), 0.0);
00181       align::RotationType rotationError( Basic3DVector<float>(0.0, 0.0, 1.0), errpar_(&ePhi) );
00182       Alignable* disk = disks[i];
00183 
00184       disk->move(translation);
00185       disk->addAlignmentPositionError(positionError);
00186       disk->rotateAroundGlobalZ(-1.0 * theFittedGlobalParametersPosTEC[aPhi]);
00187       disk->addAlignmentPositionErrorFromRotation(rotationError);
00188     }
00189 
00190 
00191   // zero initialisation (to avoid problems with the fit of NegTEC!????)
00192   zerloc_(theGlobalParametersPosTEC, theLocalParametersPosTEC);
00193   
00194   edm::LogInfo("LaserAlignmentAlgorithPosTEC") << "<LaserAlignmentAlgorithmPosTEC::doGlobalFit()>: ... done! ";
00195 }
00196 
00197 void LaserAlignmentAlgorithmPosTEC::initMillepede(int UnitForIteration)
00198 {
00199   // number of global and local parameters
00200   int nGlobalParameters = 27;
00201   int nLocalParameters = 1;
00202 
00203   // cut parameter for local fit (see page 11 of MILLEPEDE documentation)
00204   int theCut = 3;
00205   // verbose output
00206   int thePrintFlag = 1;
00207   
00208   // SIG variable in MILLEPEDE documentation on page 14
00209   float theUseParameter = 0.0;
00210 
00211   // define 0
00212   float null = 0.0;
00213 
00214   // the following two parameters are needed for iterations (see page 14)
00215   int theUnitNumber = 11 + UnitForIteration;        // the unit number for the data file
00216   float nIterations = 10000.0;   // number of iterations
00217 
00218   // the factors for the constraint of the fit (see page 15)
00219   float theConstraintFactors[27] = { 1.0, 0.0, 0.0, 1.0, 0.0, 0.0, 1.0, 0.0, 0.0, 1.0, 0.0, 0.0, 1.0, 0.0, 0.0, 
00220                                      1.0, 0.0, 0.0, 1.0, 0.0, 0.0, 1.0, 0.0, 0.0, 1.0, 0.0, 0.0 };
00221 
00222   // Define dimension parameters etc. (see page 11)
00223   initgl_(&nGlobalParameters,&nLocalParameters,&theCut,&thePrintFlag);
00224 
00225   // Fix the parameters
00226   parsig_(&theFirstFixedDiskPosTEC,&theUseParameter);
00227   parsig_(&theSecondFixedDiskPosTEC,&theUseParameter);
00228 
00229   // initialize the iterations
00230   initun_(&theUnitNumber,&nIterations);
00231 
00232   // add the constraint
00233   // in this case, the sum of the rotations around phi is zero
00234   constf_(theConstraintFactors,&null);
00235 
00236   // zero initialisation
00237   zerloc_(theGlobalParametersPosTEC, theLocalParametersPosTEC);
00238 }
00239 
00240 void LaserAlignmentAlgorithmPosTEC::resetMillepede(int UnitForIteration)
00241 {
00242   initMillepede(UnitForIteration);
00243 }

Generated on Tue Jun 9 17:24:08 2009 for CMSSW by  doxygen 1.5.4