#include <Alignment/LaserAlignment/interface/LaserAlignmentAlgorithmNegTEC.h>
Public Member Functions | |
void | addLaserBeam (std::vector< double > theMeasurements, int LaserBeam, int LaserRing) |
add a LaserBeam to Millepede and do a local fit | |
void | doGlobalFit (AlignableTracker *theAlignableTracker) |
do the global fit | |
LaserAlignmentAlgorithmNegTEC (edm::ParameterSet const &theConf, int theLaserIteration) | |
constructor | |
void | resetMillepede (int UnitForIteration) |
reset Millepede | |
~LaserAlignmentAlgorithmNegTEC () | |
destructor | |
Private Member Functions | |
void | initMillepede (int UnitForIteration) |
initialize Millepede | |
Private Attributes | |
int | theFirstFixedDiskNegTEC |
float | theGlobalParametersNegTEC [27] |
float | theLocalParametersNegTEC [1] |
int | theSecondFixedDiskNegTEC |
Definition at line 19 of file LaserAlignmentAlgorithmNegTEC.h.
LaserAlignmentAlgorithmNegTEC::LaserAlignmentAlgorithmNegTEC | ( | edm::ParameterSet const & | theConf, | |
int | theLaserIteration | |||
) |
constructor
Definition at line 17 of file LaserAlignmentAlgorithmNegTEC.cc.
References initMillepede().
00018 : theFirstFixedDiskNegTEC(theConf.getUntrackedParameter<int>("FirstFixedParameterNegTEC",2)), 00019 theSecondFixedDiskNegTEC(theConf.getUntrackedParameter<int>("SecondFixedParameterNegTEC",3)), 00020 theGlobalParametersNegTEC(), theLocalParametersNegTEC() 00021 { 00022 // initialize Millepede 00023 initMillepede(theAlignmentIteration); 00024 }
LaserAlignmentAlgorithmNegTEC::~LaserAlignmentAlgorithmNegTEC | ( | ) |
void LaserAlignmentAlgorithmNegTEC::addLaserBeam | ( | std::vector< double > | theMeasurements, | |
int | LaserBeam, | |||
int | LaserRing | |||
) |
add a LaserBeam to Millepede and do a local fit
Definition at line 30 of file LaserAlignmentAlgorithmNegTEC.cc.
References funct::cos(), equloc_(), fitloc_(), i, LogDebug, funct::sin(), theGlobalParametersNegTEC, and theLocalParametersNegTEC.
Referenced by LaserAlignmentNegTEC::alignment().
00031 { 00032 LogDebug("LaserAlignmentAlgorithmNegTEC") << "<LaserAlignmentAlgorithmNegTEC::addLaserBeam()>: adding a new Laser Beam ... "; 00033 00034 float RRing = 0.0; 00035 00036 if (LaserRing == 4) 00037 { RRing = 56.4; } // radius of Ring 4 00038 else if (LaserRing == 6) 00039 { RRing = 84.0; } // radius of Ring 6 00040 00041 float LaserPhi0 = 0.392699082; // phi position of the first Laserdiode in Ring 4 00042 float LaserPhi = 0.0; // phi position of the actual beam; will be calculated 00043 int nLaserBeams = 8; // number of beams in one ring 00044 00045 float PhiMeasured[9] = {0.0}; // the measured phi positions 00046 float PhiErrors[9] = {0.0}; // errors of the measurements 00047 00048 for (int i = 0; i < 9; i++) 00049 { 00050 PhiMeasured[i] = float(theMeasurements.at(2*i)); 00051 PhiErrors[i] = float(theMeasurements.at((2*i)+1)); 00052 LogDebug("LaserAlignmentAlgorithmNegTEC") << " i = " << i << "\t PhiMeasured[" << i << "] = " << PhiMeasured[i] 00053 << "\t PhiError[" << i << "] = " << PhiErrors[i]; 00054 } 00055 00056 // calculate phi for this beam 00057 LaserPhi = LaserPhi0 + float(LaserBeam * float(float(2 * M_PI) / nLaserBeams)); 00058 00059 // loop over the discs 00060 for (int theDisk = 0; theDisk < 9; theDisk++) 00061 { 00062 theLocalParametersNegTEC[0] = 1.0; 00063 theGlobalParametersNegTEC[3*theDisk] = 1.0; // displacement by dphi 00064 theGlobalParametersNegTEC[(3*theDisk)+1] = -sin(LaserPhi)/RRing; // displacement by dx 00065 theGlobalParametersNegTEC[(3*theDisk)+2] = cos(LaserPhi)/RRing; // displacement by dy 00066 00067 // create the equation for a single measurement 00068 equloc_(theGlobalParametersNegTEC, theLocalParametersNegTEC, &PhiMeasured[theDisk], &PhiErrors[theDisk]); 00069 } 00070 00071 // local fit after one laser beam has been added 00072 fitloc_(); 00073 }
void LaserAlignmentAlgorithmNegTEC::doGlobalFit | ( | AlignableTracker * | theAlignableTracker | ) |
do the global fit
Definition at line 75 of file LaserAlignmentAlgorithmNegTEC.cc.
References Alignable::addAlignmentPositionError(), Alignable::addAlignmentPositionErrorFromRotation(), Alignable::components(), funct::cos(), muonGeometry::disk, muonGeometry::disks, GeomDetEnumerators::endcap, AlignableTracker::endCaps(), errpar_(), fitglo_(), Alignable::globalPosition(), i, j, LogDebug, Alignable::move(), p1, p2, p3, p4, p5, funct::pow(), Alignable::rotateAroundGlobalZ(), funct::sin(), funct::sqrt(), theGlobalParametersNegTEC, theLocalParametersNegTEC, and zerloc_().
Referenced by LaserAlignmentNegTEC::alignment().
00076 { 00077 edm::LogInfo("LaserAlignmentAlgorithmNegTEC") << "<LaserAlignmentAlgorithmNegTEC::doGlobalFit()>: do the global fit ... "; 00078 00079 // number of global parameters 00080 const int nGlobalParameters = 27; 00081 // array to store the results of the fit 00082 float theFittedGlobalParametersNegTEC[nGlobalParameters] = { 0.0 }; 00083 00084 // float RRing4 = 56.4; // radius of Ring 4 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 00092 thePhiCorrected.clear(); 00093 thePhiCorrectedError.clear(); 00094 00095 // do the fit 00096 fitglo_(theFittedGlobalParametersNegTEC); 00097 00098 int ep = 1; 00099 for (int i = 0; i < nGlobalParameters; i++) 00100 { 00101 LogDebug("LaserAlignmentAlgorithmNegTEC") << "Global Parameter (TEC-) " << i << " = " << theFittedGlobalParametersNegTEC[i] 00102 << " +/- " << errpar_(&ep); 00103 ep++; 00104 } 00105 00106 int p0 = 0, p1 = 1, p2 = 2, p3 = 3, p4 = 4, p5 = 5; 00107 int n1 = 1, n2 = 2, n3 = 3, n4 = 4, n5 = 5, n6 =6; 00108 for (int j = 0; j < 8; ++j) 00109 { 00110 thePhiCorrected.push_back( theFittedGlobalParametersNegTEC[p0] 00111 + (sin(LaserPhi0)/RRing6) * theFittedGlobalParametersNegTEC[p1] 00112 - (cos(LaserPhi0)/RRing6) * theFittedGlobalParametersNegTEC[p2] 00113 - ( theFittedGlobalParametersNegTEC[p3] 00114 + (sin(LaserPhi0)/RRing6) * theFittedGlobalParametersNegTEC[p4] 00115 - (cos(LaserPhi0)/RRing6) * theFittedGlobalParametersNegTEC[p5])); 00116 thePhiCorrectedError.push_back(sqrt( pow(errpar_(&n1),2) 00117 + pow(sin(LaserPhi0)/RRing6,2) * pow(errpar_(&n2),2) 00118 + pow(cos(LaserPhi0)/RRing6,2) * pow(errpar_(&n3),2) 00119 + pow(errpar_(&n4),2) 00120 + pow(sin(LaserPhi0)/RRing6,2) * pow(errpar_(&n5),2) 00121 + pow(cos(LaserPhi0)/RRing6,2) * pow(errpar_(&n6),2) ) ); 00122 00123 00124 edm::LogInfo("LaserAlignmentAlgorithmNegTEC:Results") << " Fitted Correction for TEC- in Phi[" << j << "] = " << thePhiCorrected.at(j) << " +/- " 00125 << thePhiCorrectedError.at(j); 00126 00127 p3 += 3; p4 += 3; p5 += 3; 00128 n4 += 3; n5 += 3; n6 += 3; 00129 } 00130 00131 // loop over all discs, access the AlignableTracker to move the discs 00132 // according to the calculated alignment corrections 00133 // AlignableTracker will take care to the propagation of the movements 00134 // to the lowest level of alignable objects 00135 const align::Alignables& endcaps = theAlignableTracker->endCaps(); 00136 const Alignable* endcap = endcaps[1]; 00137 if (endcap->globalPosition().z() > 0) endcap = endcaps[0]; 00138 const align::Alignables& disks = endcap->components(); 00139 00140 for (unsigned int i = 0; i < disks.size(); ++i) 00141 { 00142 int aPhi = 3*i; 00143 int aX = 3*i + 1; 00144 int aY = 3*i + 2; 00145 int ePhi = 3*i + 1; 00146 int eX = 3*i + 2; 00147 int eY = 3*i + 3; 00148 00149 align::GlobalVector translation(-1.0 * theFittedGlobalParametersNegTEC[aX], 00150 -1.0 * theFittedGlobalParametersNegTEC[aY], 00151 0.0); 00152 AlignmentPositionError positionError(errpar_(&eX), errpar_(&eY), 0.0); 00153 align::RotationType rotationError( Basic3DVector<align::Scalar>(0.0, 0.0, 1.0), errpar_(&ePhi) ); 00154 Alignable* disk = disks[i]; 00155 00156 disk->move(translation); 00157 disk->addAlignmentPositionError(positionError); 00158 disk->rotateAroundGlobalZ(-1.0 * theFittedGlobalParametersNegTEC[aPhi]); 00159 disk->addAlignmentPositionErrorFromRotation(rotationError); 00160 00161 } 00162 00163 // zero initialisation (to avoid problems with the next Millepede fit!????) 00164 zerloc_(theGlobalParametersNegTEC, theLocalParametersNegTEC); 00165 00166 edm::LogInfo("LaserAlignmentAlgorithmNegTEC") << "<LaserAlignmentAlgorithmNegTEC::doGlobalFit()>: ... done! "; 00167 }
initialize Millepede
Definition at line 169 of file LaserAlignmentAlgorithmNegTEC.cc.
References constf_(), initgl_(), initun_(), parsig_(), theFirstFixedDiskNegTEC, theGlobalParametersNegTEC, theLocalParametersNegTEC, theSecondFixedDiskNegTEC, and zerloc_().
Referenced by LaserAlignmentAlgorithmNegTEC(), and resetMillepede().
00170 { 00171 // number of global and local parameters 00172 int nGlobalParameters = 27; 00173 int nLocalParameters = 1; 00174 00175 // cut parameter for local fit (see page 11 of MILLEPEDE documentation) 00176 int theCut = 3; 00177 // verbose output 00178 int thePrintFlag = 1; 00179 00180 // SIG variable in MILLEPEDE documentation on page 14 00181 float theUseParameter = 0.0; 00182 00183 // define 0 00184 float null = 0.0; 00185 00186 // the following two parameters are needed for iterations (see page 14) 00187 int theUnitNumber = 41 + UnitForIteration; // the unit number for the data file 00188 float nIterations = 10000.0; // number of iterations 00189 00190 // the factors for the constraint of the fit (see page 15) 00191 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, 00192 1.0, 0.0, 0.0, 1.0, 0.0, 0.0, 1.0, 0.0, 0.0, 1.0, 0.0, 0.0 }; 00193 00194 // Define dimension parameters etc. (see page 11) 00195 initgl_(&nGlobalParameters,&nLocalParameters,&theCut,&thePrintFlag); 00196 00197 // Fix the parameters 00198 parsig_(&theFirstFixedDiskNegTEC,&theUseParameter); 00199 parsig_(&theSecondFixedDiskNegTEC,&theUseParameter); 00200 00201 // initialize the iterations 00202 initun_(&theUnitNumber,&nIterations); 00203 00204 // add the constraint 00205 // in this case, the sum of the rotations around phi is zero 00206 constf_(theConstraintFactors,&null); 00207 00208 // zero initialisation 00209 zerloc_(theGlobalParametersNegTEC, theLocalParametersNegTEC); 00210 00211 }
reset Millepede
Definition at line 213 of file LaserAlignmentAlgorithmNegTEC.cc.
References initMillepede().
00214 { 00215 initMillepede(UnitForIteration); 00216 }
float LaserAlignmentAlgorithmNegTEC::theGlobalParametersNegTEC[27] [private] |
Definition at line 45 of file LaserAlignmentAlgorithmNegTEC.h.
Referenced by addLaserBeam(), doGlobalFit(), and initMillepede().
float LaserAlignmentAlgorithmNegTEC::theLocalParametersNegTEC[1] [private] |
Definition at line 46 of file LaserAlignmentAlgorithmNegTEC.h.
Referenced by addLaserBeam(), doGlobalFit(), and initMillepede().