00001
00010 #include "Alignment/LaserAlignment/plugins/LaserAlignment.h"
00011 #include "FWCore/Framework/interface/Run.h"
00012
00013
00014
00015
00019 LaserAlignment::LaserAlignment( edm::ParameterSet const& theConf ) :
00020 theEvents(0),
00021 theDoPedestalSubtraction( theConf.getUntrackedParameter<bool>( "SubtractPedestals", true ) ),
00022 theUseMinuitAlgorithm( theConf.getUntrackedParameter<bool>( "RunMinuitAlignmentTubeAlgorithm", false ) ),
00023 theApplyBeamKinkCorrections( theConf.getUntrackedParameter<bool>( "ApplyBeamKinkCorrections", true ) ),
00024 peakFinderThreshold( theConf.getUntrackedParameter<double>( "PeakFinderThreshold", 10. ) ),
00025 enableJudgeZeroFilter( theConf.getUntrackedParameter<bool>( "EnableJudgeZeroFilter", true ) ),
00026 judgeOverdriveThreshold( theConf.getUntrackedParameter<unsigned int>( "JudgeOverdriveThreshold", 220 ) ),
00027 updateFromInputGeometry( theConf.getUntrackedParameter<bool>( "UpdateFromInputGeometry", false ) ),
00028 misalignedByRefGeometry( theConf.getUntrackedParameter<bool>( "MisalignedByRefGeometry", false ) ),
00029 theStoreToDB ( theConf.getUntrackedParameter<bool>( "SaveToDbase", false ) ),
00030 theDigiProducersList( theConf.getParameter<std::vector<edm::ParameterSet> >( "DigiProducersList" ) ),
00031 theSaveHistograms( theConf.getUntrackedParameter<bool>( "SaveHistograms", false ) ),
00032 theCompression( theConf.getUntrackedParameter<int>( "ROOTFileCompression", 1 ) ),
00033 theFileName( theConf.getUntrackedParameter<std::string>( "ROOTFileName", "test.root" ) ),
00034 theMaskTecModules( theConf.getUntrackedParameter<std::vector<unsigned int> >( "MaskTECModules" ) ),
00035 theMaskAtModules( theConf.getUntrackedParameter<std::vector<unsigned int> >( "MaskATModules" ) ),
00036 theSetNominalStrips( theConf.getUntrackedParameter<bool>( "ForceFitterToNominalStrips", false ) ),
00037 theLasConstants( theConf.getUntrackedParameter<std::vector<edm::ParameterSet> >( "LaserAlignmentConstants" ) ),
00038 theFile(),
00039 theAlignableTracker(),
00040 theAlignRecordName( "TrackerAlignmentRcd" ),
00041 theErrorRecordName( "TrackerAlignmentErrorRcd" ),
00042 firstEvent_(true) {
00043
00044
00045 std::cout << std::endl;
00046 std::cout << "=============================================================="
00047 << "\n=== LaserAlignment module configuration ==="
00048 << "\n"
00049 << "\n Write histograms to file = " << (theSaveHistograms?"true":"false")
00050 << "\n Histogram file name = " << theFileName
00051 << "\n Histogram file compression = " << theCompression
00052 << "\n Subtract pedestals = " << (theDoPedestalSubtraction?"true":"false")
00053 << "\n Run Minuit AT algorithm = " << (theUseMinuitAlgorithm?"true":"false")
00054 << "\n Apply beam kink corrections = " << (theApplyBeamKinkCorrections?"true":"false")
00055 << "\n Peak Finder Threshold = " << peakFinderThreshold
00056 << "\n EnableJudgeZeroFilter = " << (enableJudgeZeroFilter?"true":"false")
00057 << "\n JudgeOverdriveThreshold = " << judgeOverdriveThreshold
00058 << "\n Update from input geometry = " << (updateFromInputGeometry?"true":"false")
00059 << "\n Misalignment from ref geometry = " << (misalignedByRefGeometry?"true":"false")
00060 << "\n Number of TEC modules masked = " << theMaskTecModules.size() << " (s. below list if > 0)"
00061 << "\n Number of AT modules masked = " << theMaskAtModules.size() << " (s. below list if > 0)"
00062 << "\n Store to database = " << (theStoreToDB?"true":"false")
00063 << "\n ----------------------------------------------- ----------"
00064 << (theSetNominalStrips?"\n Set strips to nominal = true":"\n")
00065 << "\n=============================================================" << std::endl;
00066
00067
00068
00069 if( theMaskTecModules.size() ) {
00070 std::cout << " ===============================================================================================\n" << std::flush;
00071 std::cout << " The following " << theMaskTecModules.size() << " TEC modules have been masked out and will not be considered by the TEC algorithm:\n " << std::flush;
00072 for( std::vector<unsigned int>::iterator moduleIt = theMaskTecModules.begin(); moduleIt != theMaskTecModules.end(); ++moduleIt ) {
00073 std::cout << *moduleIt << (moduleIt!=--theMaskTecModules.end()?", ":"") << std::flush;
00074 }
00075 std::cout << std::endl << std::flush;
00076 std::cout << " ===============================================================================================\n\n" << std::flush;
00077 }
00078 if( theMaskAtModules.size() ) {
00079 std::cout << " ===============================================================================================\n" << std::flush;
00080 std::cout << " The following " << theMaskAtModules.size() << " AT modules have been masked out and will not be considered by the AT algorithm:\n " << std::flush;
00081 for( std::vector<unsigned int>::iterator moduleIt = theMaskAtModules.begin(); moduleIt != theMaskAtModules.end(); ++moduleIt ) {
00082 std::cout << *moduleIt << (moduleIt!=--theMaskAtModules.end()?", ":"") << std::flush;
00083 }
00084 std::cout << std::endl << std::flush;
00085 std::cout << " ===============================================================================================\n\n" << std::flush;
00086 }
00087
00088
00089
00090
00091 std::string alias ( theConf.getParameter<std::string>("@module_label") );
00092
00093
00094 produces<TkLasBeamCollection, edm::InRun>( "tkLaserBeams" ).setBranchAlias( alias + "TkLasBeamCollection" );
00095
00096
00097 judge.EnableZeroFilter( enableJudgeZeroFilter );
00098
00099
00100 judge.SetOverdriveThreshold( judgeOverdriveThreshold );
00101
00102 }
00103
00104
00105
00106
00107
00111 LaserAlignment::~LaserAlignment() {
00112
00113 if ( theSaveHistograms ) theFile->Write();
00114 if ( theFile ) { delete theFile; }
00115 if ( theAlignableTracker ) { delete theAlignableTracker; }
00116
00117 }
00118
00119
00120
00121
00122
00126 void LaserAlignment::beginJob() {
00127
00128
00129
00130 if( theSaveHistograms ) {
00131
00132
00133 theFile = new TFile( theFileName.c_str(), "RECREATE", "CMS ROOT file" );
00134
00135
00136 if ( theFile ) {
00137 theFile->SetCompressionLevel(theCompression);
00138 singleModulesDir = theFile->mkdir( "single modules" );
00139 } else
00140 throw cms::Exception( " [LaserAlignment::beginJob]") << " ** ERROR: could not open file:"
00141 << theFileName.c_str() << " for writing." << std::endl;
00142
00143 }
00144
00145
00146 fillDetectorId();
00147
00148
00149
00150
00151
00152
00153 std::stringstream nameBuilder;
00154
00155
00156 int det, ring, beam, disk, pos;
00157
00158
00159 det = 0; ring = 0; beam = 0; disk = 0;
00160 do {
00161
00162 pedestalProfiles.GetTECEntry( det, ring, beam, disk ).SetAllValuesTo( 0. );
00163 currentDataProfiles.GetTECEntry( det, ring, beam, disk ).SetAllValuesTo( 0. );
00164 collectedDataProfiles.GetTECEntry( det, ring, beam, disk ).SetAllValuesTo( 0. );
00165
00166
00167 isAcceptedProfile.SetTECEntry( det, ring, beam, disk, 0 );
00168 numberOfAcceptedProfiles.SetTECEntry( det, ring, beam, disk, 0 );
00169
00170
00171 nameBuilder.clear();
00172 nameBuilder.str( "" );
00173 nameBuilder << "TEC";
00174 if( det == 0 ) nameBuilder << "+"; else nameBuilder << "-";
00175 nameBuilder << "_Ring";
00176 if( ring == 0 ) nameBuilder << "4"; else nameBuilder << "6";
00177 nameBuilder << "_Beam" << beam;
00178 nameBuilder << "_Disk" << disk;
00179 theProfileNames.SetTECEntry( det, ring, beam, disk, nameBuilder.str() );
00180
00181
00182 if( theSaveHistograms ) {
00183 nameBuilder << "_Histo";
00184 summedHistograms.SetTECEntry( det, ring, beam, disk, new TH1D( nameBuilder.str().c_str(), nameBuilder.str().c_str(), 512, 0, 512 ) );
00185 summedHistograms.GetTECEntry( det, ring, beam, disk )->SetDirectory( singleModulesDir );
00186 }
00187
00188 } while ( moduleLoop.TECLoop( det, ring, beam, disk ) );
00189
00190
00191
00192 det = 2; beam = 0; pos = 0;
00193 do {
00194
00195 pedestalProfiles.GetTIBTOBEntry( det, beam, pos ).SetAllValuesTo( 0. );
00196 currentDataProfiles.GetTIBTOBEntry( det, beam, pos ).SetAllValuesTo( 0. );
00197 collectedDataProfiles.GetTIBTOBEntry( det, beam, pos ).SetAllValuesTo( 0. );
00198
00199
00200 isAcceptedProfile.SetTIBTOBEntry( det, beam, pos, 0 );
00201 numberOfAcceptedProfiles.SetTIBTOBEntry( det, beam, pos, 0 );
00202
00203
00204 nameBuilder.clear();
00205 nameBuilder.str( "" );
00206 if( det == 2 ) nameBuilder << "TIB"; else nameBuilder << "TOB";
00207 nameBuilder << "_Beam" << beam;
00208 nameBuilder << "_Zpos" << pos;
00209
00210 theProfileNames.SetTIBTOBEntry( det, beam, pos, nameBuilder.str() );
00211
00212
00213 if( theSaveHistograms ) {
00214 nameBuilder << "_Histo";
00215 summedHistograms.SetTIBTOBEntry( det, beam, pos, new TH1D( nameBuilder.str().c_str(), nameBuilder.str().c_str(), 512, 0, 512 ) );
00216 summedHistograms.GetTIBTOBEntry( det, beam, pos )->SetDirectory( singleModulesDir );
00217 }
00218
00219 } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
00220
00221
00222
00223 det = 0; beam = 0; disk = 0;
00224 do {
00225
00226 pedestalProfiles.GetTEC2TECEntry( det, beam, disk ).SetAllValuesTo( 0. );
00227 currentDataProfiles.GetTEC2TECEntry( det, beam, disk ).SetAllValuesTo( 0. );
00228 collectedDataProfiles.GetTEC2TECEntry( det, beam, disk ).SetAllValuesTo( 0. );
00229
00230
00231 isAcceptedProfile.SetTEC2TECEntry( det, beam, disk, 0 );
00232 numberOfAcceptedProfiles.SetTEC2TECEntry( det, beam, disk, 0 );
00233
00234
00235 nameBuilder.clear();
00236 nameBuilder.str( "" );
00237 nameBuilder << "TEC(AT)";
00238 if( det == 0 ) nameBuilder << "+"; else nameBuilder << "-";
00239 nameBuilder << "_Beam" << beam;
00240 nameBuilder << "_Disk" << disk;
00241 theProfileNames.SetTEC2TECEntry( det, beam, disk, nameBuilder.str() );
00242
00243
00244 if( theSaveHistograms ) {
00245 nameBuilder << "_Histo";
00246 summedHistograms.SetTEC2TECEntry( det, beam, disk, new TH1D( nameBuilder.str().c_str(), nameBuilder.str().c_str(), 512, 0, 512 ) );
00247 summedHistograms.GetTEC2TECEntry( det, beam, disk )->SetDirectory( singleModulesDir );
00248 }
00249
00250 } while( moduleLoop.TEC2TECLoop( det, beam, disk ) );
00251
00252 firstEvent_ = true;
00253 }
00254
00255
00256
00257
00258
00259
00263 void LaserAlignment::produce(edm::Event& theEvent, edm::EventSetup const& theSetup) {
00264
00265 if (firstEvent_) {
00266
00267
00268 theSetup.get<TrackerDigiGeometryRecord>().get( theTrackerGeometry );
00269 theSetup.get<IdealGeometryRecord>().get( gD );
00270
00271
00272 edm::ESHandle<SiStripPedestals> pedestalsHandle;
00273 if( theDoPedestalSubtraction ) {
00274 theSetup.get<SiStripPedestalsRcd>().get( pedestalsHandle );
00275 fillPedestalProfiles( pedestalsHandle );
00276 }
00277
00278
00279
00280 theSetup.get<TrackerDigiGeometryRecord>().getRecord<GlobalPositionRcd>().get( theGlobalPositionRcd );
00281
00282
00283 if( !updateFromInputGeometry ) {
00284
00285 edm::ESHandle<GeometricDet> theGeometricDet;
00286 theSetup.get<IdealGeometryRecord>().get(theGeometricDet);
00287 TrackerGeomBuilderFromGeometricDet trackerBuilder;
00288 TrackerGeometry* theRefTracker = trackerBuilder.build(&*theGeometricDet);
00289 theAlignableTracker = new AlignableTracker(&(*theRefTracker));
00290 }
00291 else {
00292
00293 theAlignableTracker = new AlignableTracker( &(*theTrackerGeometry) );
00294 }
00295
00296 firstEvent_ = false;
00297 }
00298
00299 LogDebug("LaserAlignment") << "==========================================================="
00300 << "\n Private analysis of event #"<< theEvent.id().event()
00301 << " in run #" << theEvent.id().run();
00302
00303
00304
00305 fillDataProfiles( theEvent, theSetup );
00306
00307
00308 int det, ring, beam, disk, pos;
00309
00310
00311
00312
00313
00314
00315
00316
00317
00318 ring = 1;
00319 for( det = 0; det < 2; ++det ) {
00320 for( beam = 0; beam < 8; ++ beam ) {
00321 for( disk = 0; disk < 9; ++disk ) {
00322 if( judge.IsSignalIn( currentDataProfiles.GetTECEntry( det, ring, beam, disk ) - pedestalProfiles.GetTECEntry( det, ring, beam, disk ), 0 ) ) {
00323 isAcceptedProfile.SetTECEntry( det, ring, beam, disk, 1 );
00324 }
00325 else {
00326 isAcceptedProfile.SetTECEntry( det, ring, beam, disk, 0 );
00327 }
00328 }
00329 }
00330 }
00331
00332
00333 det = 2; beam = 0; pos = 0;
00334 do {
00335
00336 if( judge.IsSignalIn( currentDataProfiles.GetTIBTOBEntry( det, beam, pos ) - pedestalProfiles.GetTIBTOBEntry( det, beam, pos ), getTIBTOBNominalBeamOffset( det, beam, pos ) ) ) {
00337 isAcceptedProfile.SetTIBTOBEntry( det, beam, pos, 1 );
00338 }
00339 else {
00340 isAcceptedProfile.SetTIBTOBEntry( det, beam, pos, 0 );
00341 }
00342
00343 } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
00344
00345
00346
00347
00348
00349 bool isTECMode = isTECBeam();
00350
00351 std::cout << " [LaserAlignment::produce] -- LaserAlignment::isTECBeam declares this event " << ( isTECMode ? "" : "NOT " ) << "a TEC event." << std::endl;
00352
00353 bool isATMode = isATBeam();
00354
00355 std::cout << " [LaserAlignment::produce] -- LaserAlignment::isATBeam declares this event " << ( isATMode ? "" : "NOT " ) << "an AT event." << std::endl;
00356
00357
00358
00359
00360
00361
00362
00363
00364
00365
00366
00367
00368 det = 0; ring = 0; beam = 0; disk = 0;
00369 do {
00370
00371 LogDebug( "[LaserAlignment::produce]" ) << "Profile is: " << theProfileNames.GetTECEntry( det, ring, beam, disk ) << "." << std::endl;
00372
00373
00374
00375 if( ring == 0 && find( tecDoubleHitDetId.begin(), tecDoubleHitDetId.end(), detectorId.GetTECEntry( det, ring, beam, disk ) ) != tecDoubleHitDetId.end() ) {
00376
00377 if( isTECMode ) {
00378
00379 if( judge.JudgeProfile( currentDataProfiles.GetTECEntry( det, ring, beam, disk ) - pedestalProfiles.GetTECEntry( det, ring, beam, disk ), 0 ) ) {
00380 collectedDataProfiles.GetTECEntry( det, ring, beam, disk ) += currentDataProfiles.GetTECEntry( det, ring, beam, disk ) - pedestalProfiles.GetTECEntry( det, ring, beam, disk );
00381 numberOfAcceptedProfiles.GetTECEntry( det, ring, beam, disk )++;
00382 }
00383 }
00384 }
00385
00386 else {
00387
00388 if( judge.JudgeProfile( currentDataProfiles.GetTECEntry( det, ring, beam, disk ) - pedestalProfiles.GetTECEntry( det, ring, beam, disk ), 0 ) ) {
00389 collectedDataProfiles.GetTECEntry( det, ring, beam, disk ) += currentDataProfiles.GetTECEntry( det, ring, beam, disk ) - pedestalProfiles.GetTECEntry( det, ring, beam, disk );
00390 numberOfAcceptedProfiles.GetTECEntry( det, ring, beam, disk )++;
00391 }
00392 }
00393
00394 } while( moduleLoop.TECLoop( det, ring, beam, disk ) );
00395
00396
00397
00398
00399
00400
00401 det = 2; beam = 0; pos = 0;
00402 do {
00403
00404 LogDebug( "[LaserAlignment::produce]" ) << "Profile is: " << theProfileNames.GetTIBTOBEntry( det, beam, pos ) << "." << std::endl;
00405
00406
00407 if( judge.JudgeProfile( currentDataProfiles.GetTIBTOBEntry( det, beam, pos ) - pedestalProfiles.GetTIBTOBEntry( det, beam, pos ), getTIBTOBNominalBeamOffset( det, beam, pos ) ) ) {
00408 collectedDataProfiles.GetTIBTOBEntry( det, beam, pos ) += currentDataProfiles.GetTIBTOBEntry( det, beam, pos ) - pedestalProfiles.GetTIBTOBEntry( det, beam, pos );
00409 numberOfAcceptedProfiles.GetTIBTOBEntry( det, beam, pos )++;
00410 }
00411
00412 } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
00413
00414
00415
00416
00417 det = 0; beam = 0; disk = 0;
00418 do {
00419
00420 LogDebug( "[LaserAlignment::produce]" ) << "Profile is: " << theProfileNames.GetTEC2TECEntry( det, beam, disk ) << "." << std::endl;
00421
00422
00423
00424 if( ring == 0 && find( tecDoubleHitDetId.begin(), tecDoubleHitDetId.end(), detectorId.GetTECEntry( det, ring, beam, disk ) ) != tecDoubleHitDetId.end() ) {
00425
00426 if( isATMode ) {
00427
00428 if( judge.JudgeProfile( currentDataProfiles.GetTEC2TECEntry( det, beam, disk ) - pedestalProfiles.GetTEC2TECEntry( det, beam, disk ), 0 ) ) {
00429 collectedDataProfiles.GetTEC2TECEntry( det, beam, disk ) += currentDataProfiles.GetTEC2TECEntry( det, beam, disk ) - pedestalProfiles.GetTEC2TECEntry( det, beam, disk );
00430 numberOfAcceptedProfiles.GetTEC2TECEntry( det, beam, disk )++;
00431 }
00432 }
00433
00434 }
00435
00436 else {
00437
00438 if( judge.JudgeProfile( currentDataProfiles.GetTEC2TECEntry( det, beam, disk ) - pedestalProfiles.GetTEC2TECEntry( det, beam, disk ), 0 ) ) {
00439 collectedDataProfiles.GetTEC2TECEntry( det, beam, disk ) += currentDataProfiles.GetTEC2TECEntry( det, beam, disk ) - pedestalProfiles.GetTEC2TECEntry( det, beam, disk );
00440 numberOfAcceptedProfiles.GetTEC2TECEntry( det, beam, disk )++;
00441 }
00442 }
00443
00444
00445 } while( moduleLoop.TEC2TECLoop( det, beam, disk ) );
00446
00447
00448
00449
00450 theEvents++;
00451
00452 }
00453
00454
00455
00456
00457
00461 void LaserAlignment::endRun( edm::Run& theRun, const edm::EventSetup& theSetup ) {
00462
00463
00464 std::cout << " [LaserAlignment::endRun] -- Total number of events processed: " << theEvents << std::endl;
00465
00466
00467 DumpHitmaps( numberOfAcceptedProfiles );
00468
00469
00470 int det, ring, beam, disk, pos;
00471
00472
00473 LASGlobalData<LASCoordinateSet> measuredCoordinates;
00474
00475
00476 LASGlobalData<std::pair<float,float> > measuredStripPositions;
00477
00478
00479 LASPeakFinder peakFinder;
00480 peakFinder.SetAmplitudeThreshold( peakFinderThreshold );
00481 std::pair<double,double> peakFinderResults;
00482 bool isGoodFit;
00483
00484
00485 const TrackerGeometry& theTracker( *theTrackerGeometry );
00486
00487
00488 CalculateNominalCoordinates();
00489
00490
00491
00492
00493
00494
00495
00496
00497
00498 det = 0; ring = 0; beam = 0; disk = 0;
00499 do {
00500
00501
00502 isGoodFit = peakFinder.FindPeakIn( collectedDataProfiles.GetTECEntry( det, ring, beam, disk ), peakFinderResults,
00503 summedHistograms.GetTECEntry( det, ring, beam, disk ), 0 );
00504
00505
00506 if( !isGoodFit ) std::cout << " [LaserAlignment::endRun] ** WARNING: Fit failed for TEC det: "
00507 << det << ", ring: " << ring << ", beam: " << beam << ", disk: " << disk
00508 << " (id: " << detectorId.GetTECEntry( det, ring, beam, disk ) << ")." << std::endl;
00509
00510
00511
00512
00513
00514
00515 const DetId theDetId( detectorId.GetTECEntry( det, ring, beam, disk ) );
00516 const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
00517
00518 if (theStripDet) {
00519
00520 measuredCoordinates.SetTECEntry( det, ring, beam, disk, nominalCoordinates.GetTECEntry( det, ring, beam, disk ) );
00521
00522 if( isGoodFit ) {
00523
00524 measuredStripPositions.GetTECEntry( det, ring, beam, disk ) = peakFinderResults;
00525 const float positionInStrips = theSetNominalStrips ? 256. : peakFinderResults.first;
00526 const GlobalPoint& globalPoint = theStripDet->surface().toGlobal( theStripDet->specificTopology().localPosition( positionInStrips ) );
00527 measuredCoordinates.GetTECEntry( det, ring, beam, disk ).SetPhi( ConvertAngle( globalPoint.barePhi() ) );
00528
00529
00530
00531 measuredCoordinates.GetTECEntry( det, ring, beam, disk ).SetPhiError( 0.00046 );
00532
00533 }
00534 else {
00535 measuredStripPositions.GetTECEntry( det, ring, beam, disk ) = std::pair<float,float>( 256., 1000. );
00536 const GlobalPoint& globalPoint = theStripDet->surface().toGlobal( theStripDet->specificTopology().localPosition( 256. ) );
00537 measuredCoordinates.GetTECEntry( det, ring, beam, disk ).SetPhi( ConvertAngle( globalPoint.barePhi() ) );
00538 measuredCoordinates.GetTECEntry( det, ring, beam, disk ).SetPhiError( 1000. );
00539 }
00540 }
00541
00542 } while( moduleLoop.TECLoop( det, ring, beam, disk ) );
00543
00544
00545
00546
00547
00548 det = 2; beam = 0; pos = 0;
00549 do {
00550
00551
00552 isGoodFit = peakFinder.FindPeakIn( collectedDataProfiles.GetTIBTOBEntry( det, beam, pos ), peakFinderResults,
00553 summedHistograms.GetTIBTOBEntry( det, beam, pos ), getTIBTOBNominalBeamOffset( det, beam, pos ) );
00554
00555
00556 if( !isGoodFit ) std::cout << " [LaserAlignment::endJob] ** WARNING: Fit failed for TIB/TOB det: "
00557 << det << ", beam: " << beam << ", pos: " << pos
00558 << " (id: " << detectorId.GetTIBTOBEntry( det, beam, pos ) << ")." << std::endl;
00559
00560
00561
00562
00563
00564 const DetId theDetId( detectorId.GetTIBTOBEntry( det, beam, pos ) );
00565 const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
00566
00567 if (theStripDet) {
00568
00569 measuredCoordinates.SetTIBTOBEntry( det, beam, pos, nominalCoordinates.GetTIBTOBEntry( det, beam, pos ) );
00570
00571 if( isGoodFit ) {
00572 measuredStripPositions.GetTIBTOBEntry( det, beam, pos ) = peakFinderResults;
00573 const float positionInStrips = theSetNominalStrips ? 256. + getTIBTOBNominalBeamOffset( det, beam, pos ) : peakFinderResults.first;
00574 const GlobalPoint& globalPoint = theStripDet->surface().toGlobal( theStripDet->specificTopology().localPosition( positionInStrips ) );
00575 measuredCoordinates.GetTIBTOBEntry( det, beam, pos ).SetPhi( ConvertAngle( globalPoint.barePhi() ) );
00576 measuredCoordinates.GetTIBTOBEntry( det, beam, pos ).SetPhiError( 0.00028 );
00577 }
00578 else {
00579 measuredStripPositions.GetTIBTOBEntry( det, beam, pos ) = std::pair<float,float>( 256. + getTIBTOBNominalBeamOffset( det, beam, pos ), 1000. );
00580 const GlobalPoint& globalPoint = theStripDet->surface().toGlobal( theStripDet->specificTopology().localPosition( 256. + getTIBTOBNominalBeamOffset( det, beam, pos ) ) );
00581 measuredCoordinates.GetTIBTOBEntry( det, beam, pos ).SetPhi( ConvertAngle( globalPoint.barePhi() ) );
00582 measuredCoordinates.GetTIBTOBEntry( det, beam, pos ).SetPhiError( 1000. );
00583 }
00584 }
00585
00586 } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
00587
00588
00589
00590
00591
00592 det = 0; beam = 0; disk = 0;
00593 do {
00594
00595
00596 isGoodFit = peakFinder.FindPeakIn( collectedDataProfiles.GetTEC2TECEntry( det, beam, disk ), peakFinderResults,
00597 summedHistograms.GetTEC2TECEntry( det, beam, disk ), getTEC2TECNominalBeamOffset( det, beam, disk ) );
00598
00599 if( !isGoodFit ) std::cout << " [LaserAlignment::endRun] ** WARNING: Fit failed for TEC2TEC det: "
00600 << det << ", beam: " << beam << ", disk: " << disk
00601 << " (id: " << detectorId.GetTEC2TECEntry( det, beam, disk ) << ")." << std::endl;
00602
00603
00604
00605
00606
00607 const DetId theDetId( detectorId.GetTEC2TECEntry( det, beam, disk ) );
00608 const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
00609
00610 if (theStripDet) {
00611
00612 measuredCoordinates.SetTEC2TECEntry( det, beam, disk, nominalCoordinates.GetTEC2TECEntry( det, beam, disk ) );
00613
00614 if( isGoodFit ) {
00615 measuredStripPositions.GetTEC2TECEntry( det, beam, disk ) = peakFinderResults;
00616 const float positionInStrips = theSetNominalStrips ? 256. + getTEC2TECNominalBeamOffset( det, beam, disk ) : peakFinderResults.first;
00617 const GlobalPoint& globalPoint = theStripDet->surface().toGlobal( theStripDet->specificTopology().localPosition( positionInStrips ) );
00618 measuredCoordinates.GetTEC2TECEntry( det, beam, disk ).SetPhi( ConvertAngle( globalPoint.barePhi() ) );
00619 measuredCoordinates.GetTEC2TECEntry( det, beam, disk ).SetPhiError( 0.00047 );
00620 }
00621 else {
00622 measuredStripPositions.GetTEC2TECEntry( det, beam, disk ) = std::pair<float,float>( 256. + getTEC2TECNominalBeamOffset( det, beam, disk ), 1000. );
00623 const GlobalPoint& globalPoint = theStripDet->surface().toGlobal( theStripDet->specificTopology().localPosition( 256. + getTEC2TECNominalBeamOffset( det, beam, disk ) ) );
00624 measuredCoordinates.GetTEC2TECEntry( det, beam, disk ).SetPhi( ConvertAngle( globalPoint.barePhi() ) );
00625 measuredCoordinates.GetTEC2TECEntry( det, beam, disk ).SetPhiError( 1000. );
00626 }
00627 }
00628
00629 } while( moduleLoop.TEC2TECLoop( det, beam, disk ) );
00630
00631
00632
00633
00634
00635
00636
00637
00638
00639
00640
00641
00642
00643
00644
00645
00646
00647
00648 LASGeometryUpdater geometryUpdater( nominalCoordinates, theLasConstants );
00649
00650
00651 if( theApplyBeamKinkCorrections ) geometryUpdater.ApplyBeamKinkCorrections( measuredCoordinates );
00652
00653
00654
00655 if( updateFromInputGeometry ) geometryUpdater.SetReverseDirection( true );
00656
00657
00658
00659 if( misalignedByRefGeometry ) geometryUpdater.SetMisalignmentFromRefGeometry( true );
00660
00661
00662 LASEndcapAlgorithm endcapAlgorithm;
00663 LASEndcapAlignmentParameterSet endcapParameters;
00664
00665
00666
00667
00668 if( theMaskTecModules.size() ) {
00669 ApplyEndcapMaskingCorrections( measuredCoordinates, nominalCoordinates, endcapParameters );
00670 }
00671
00672
00673 endcapParameters = endcapAlgorithm.CalculateParameters( measuredCoordinates, nominalCoordinates );
00674
00675
00676
00677
00678
00679
00680
00681 if( theMaskTecModules.size() ) {
00682
00683 const unsigned int nIterations = 30;
00684 for( unsigned int iteration = 0; iteration < nIterations; ++iteration ) {
00685
00686
00687
00688 ApplyEndcapMaskingCorrections( measuredCoordinates, nominalCoordinates, endcapParameters );
00689
00690
00691 endcapParameters = endcapAlgorithm.CalculateParameters( measuredCoordinates, nominalCoordinates );
00692
00693 }
00694
00695 }
00696
00697
00698 endcapParameters.Print();
00699
00700
00701
00702
00703
00704
00705 geometryUpdater.EndcapUpdate( endcapParameters, measuredCoordinates );
00706
00707
00708
00709 LASBarrelAlignmentParameterSet alignmentTubeParameters;
00710
00711 LASBarrelAlgorithm barrelAlgorithm;
00712
00713 LASAlignmentTubeAlgorithm alignmentTubeAlgorithm;
00714
00715
00716
00717
00718 if( theMaskAtModules.size() ) {
00719 ApplyATMaskingCorrections( measuredCoordinates, nominalCoordinates, alignmentTubeParameters );
00720 }
00721
00722 if( theUseMinuitAlgorithm ) {
00723
00724 alignmentTubeParameters = barrelAlgorithm.CalculateParameters( measuredCoordinates, nominalCoordinates );
00725 }
00726 else {
00727
00728 alignmentTubeParameters = alignmentTubeAlgorithm.CalculateParameters( measuredCoordinates, nominalCoordinates );
00729 }
00730
00731
00732
00733
00734
00735
00736
00737
00738
00739 if( theMaskAtModules.size() ) {
00740
00741 const unsigned int nIterations = 30;
00742 for( unsigned int iteration = 0; iteration < nIterations; ++iteration ) {
00743
00744
00745
00746 ApplyATMaskingCorrections( measuredCoordinates, nominalCoordinates, alignmentTubeParameters );
00747
00748
00749 if( theUseMinuitAlgorithm ) {
00750 alignmentTubeParameters = barrelAlgorithm.CalculateParameters( measuredCoordinates, nominalCoordinates );
00751 }
00752 else {
00753 alignmentTubeParameters = alignmentTubeAlgorithm.CalculateParameters( measuredCoordinates, nominalCoordinates );
00754 }
00755
00756 }
00757
00758 }
00759
00760
00761
00762 alignmentTubeParameters.Print();
00763
00764
00765
00766 geometryUpdater.TrackerUpdate( endcapParameters, alignmentTubeParameters, *theAlignableTracker );
00767
00768
00769
00770
00771
00772
00773
00778
00779
00780
00781 std::auto_ptr<TkLasBeamCollection> laserBeams( new TkLasBeamCollection );
00782
00783
00784
00785 for( det = 0; det < 2; ++det ) {
00786 for( ring = 0; ring < 2; ++ring ) {
00787 for( beam = 0; beam < 8; ++beam ) {
00788
00789
00790 TkLasBeam currentBeam( 100 * det + 10 * beam + ring );
00791
00792
00793 const int firstDisk = det==0 ? 0 : 8;
00794 const int lastDisk = det==0 ? 8 : 0;
00795
00796
00797 for( disk = firstDisk; det==0 ? disk <= lastDisk : disk >= lastDisk; det==0 ? ++disk : --disk ) {
00798
00799
00800 const SiStripDetId theDetId( detectorId.GetTECEntry( det, ring, beam, disk ) );
00801
00802
00803 const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
00804
00805
00806 const SiStripLaserRecHit2D currentHit(
00807 theStripDet->specificTopology().localPosition( measuredStripPositions.GetTECEntry( det, ring, beam, disk ).first ),
00808 theStripDet->specificTopology().localError( measuredStripPositions.GetTECEntry( det, ring, beam, disk ).first, measuredStripPositions.GetTECEntry( det, ring, beam, disk ).second ),
00809 theDetId
00810 );
00811
00812 currentBeam.push_back( currentHit );
00813
00814 }
00815
00816 laserBeams->push_back( currentBeam );
00817
00818 }
00819 }
00820 }
00821
00822
00823
00824
00825
00826
00827 for( beam = 0; beam < 8; ++beam ) {
00828
00829
00830 TkLasBeam currentBeam( 100 * 2 + 10 * beam + 0 );
00831
00832
00833
00834 det = 1;
00835 for( disk = 4; disk >= 0; --disk ) {
00836
00837
00838 const SiStripDetId theDetId( detectorId.GetTEC2TECEntry( det, beam, disk ) );
00839
00840
00841 const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
00842
00843
00844 const SiStripLaserRecHit2D currentHit(
00845 theStripDet->specificTopology().localPosition( measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).first ),
00846 theStripDet->specificTopology().localError( measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).first, measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).second ),
00847 theDetId
00848 );
00849
00850 currentBeam.push_back( currentHit );
00851
00852 }
00853
00854
00855
00856 for( det = 2; det < 4; ++det ) {
00857 for( pos = 5; pos >= 0; --pos ) {
00858
00859
00860 const SiStripDetId theDetId( detectorId.GetTIBTOBEntry( det, beam, pos ) );
00861
00862
00863 const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
00864
00865
00866 const SiStripLaserRecHit2D currentHit(
00867 theStripDet->specificTopology().localPosition( measuredStripPositions.GetTIBTOBEntry( det, beam, pos ).first ),
00868 theStripDet->specificTopology().localError( measuredStripPositions.GetTIBTOBEntry( det, beam, pos ).first, measuredStripPositions.GetTIBTOBEntry( det, beam, pos ).second ),
00869 theDetId
00870 );
00871
00872 currentBeam.push_back( currentHit );
00873
00874 }
00875 }
00876
00877
00878
00879 det = 0;
00880 for( disk = 0; disk < 5; ++disk ) {
00881
00882
00883 const SiStripDetId theDetId( detectorId.GetTEC2TECEntry( det, beam, disk ) );
00884
00885
00886 const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
00887
00888
00889 const SiStripLaserRecHit2D currentHit(
00890 theStripDet->specificTopology().localPosition( measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).first ),
00891 theStripDet->specificTopology().localError( measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).first, measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).second ),
00892 theDetId
00893 );
00894
00895 currentBeam.push_back( currentHit );
00896
00897 }
00898
00899
00900
00901
00902 laserBeams->push_back( currentBeam );
00903
00904 }
00905
00906
00907
00908 theRun.put( laserBeams, "tkLaserBeams" );
00909
00910
00911
00912
00913
00914
00915
00916 Alignments* alignments = theAlignableTracker->alignments();
00917 AlignmentErrors* alignmentErrors = theAlignableTracker->alignmentErrors();
00918
00919 if ( theStoreToDB ) {
00920
00921 std::cout << " [LaserAlignment::endRun] -- Storing the calculated alignment parameters to the DataBase:" << std::endl;
00922
00923
00924 edm::Service<cond::service::PoolDBOutputService> poolDbService;
00925 if( !poolDbService.isAvailable() )
00926 throw cms::Exception( "NotAvailable" ) << "PoolDBOutputService not available";
00927
00928
00929
00930
00931
00932
00933
00934
00935
00936 poolDbService->writeOne<Alignments>( alignments, poolDbService->beginOfTime(), theAlignRecordName );
00937
00938
00939
00940
00941
00942
00943
00944 poolDbService->writeOne<AlignmentErrors>( alignmentErrors, poolDbService->beginOfTime(), theErrorRecordName );
00945
00946 std::cout << " [LaserAlignment::endRun] -- Storing done." << std::endl;
00947
00948 }
00949
00950 }
00951
00952
00953
00954
00955
00959 void LaserAlignment::endJob() {
00960 }
00961
00962
00963
00964
00965
00970 void LaserAlignment::fillDataProfiles( edm::Event const& theEvent, edm::EventSetup const& theSetup ) {
00971
00972
00973 edm::Handle< edm::DetSetVector<SiStripRawDigi> > theStripRawDigis;
00974 edm::Handle< edm::DetSetVector<SiStripDigi> > theStripDigis;
00975
00976 bool isRawDigi = false;
00977
00978
00979 int det = 0, ring = 0, beam = 0, disk = 0, pos = 0;
00980
00981
00982 for ( std::vector<edm::ParameterSet>::iterator itDigiProducersList = theDigiProducersList.begin(); itDigiProducersList != theDigiProducersList.end(); ++itDigiProducersList ) {
00983
00984 std::string digiProducer = itDigiProducersList->getParameter<std::string>( "DigiProducer" );
00985 std::string digiLabel = itDigiProducersList->getParameter<std::string>( "DigiLabel" );
00986 std::string digiType = itDigiProducersList->getParameter<std::string>( "DigiType" );
00987
00988
00989
00990 if( digiType == "Raw" ) {
00991 theEvent.getByLabel( digiProducer, digiLabel, theStripRawDigis );
00992 isRawDigi = true;
00993 }
00994 else if( digiType == "Processed" ) {
00995 theEvent.getByLabel( digiProducer, digiLabel, theStripDigis );
00996 isRawDigi = false;
00997 }
00998 else {
00999 throw cms::Exception( " [LaserAlignment::fillDataProfiles]") << " ** ERROR: Invalid digi type: \"" << digiType << "\" specified in configuration." << std::endl;
01000 }
01001
01002
01003
01004
01005 det = 0; ring = 0; beam = 0; disk = 0;
01006 do {
01007
01008
01009 currentDataProfiles.GetTECEntry( det, ring, beam, disk ).SetAllValuesTo( 0. );
01010
01011
01012 const int detRawId = detectorId.GetTECEntry( det, ring, beam, disk );
01013
01014 if( isRawDigi ) {
01015
01016
01017 edm::DetSetVector<SiStripRawDigi>::const_iterator detSetIter = theStripRawDigis->find( detRawId );
01018 if( detSetIter == theStripRawDigis->end() ) {
01019 throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: No raw DetSet found for det: " << detRawId << "." << std::endl;
01020 }
01021
01022
01023 edm::DetSet<SiStripRawDigi>::const_iterator digiRangeIterator = detSetIter->data.begin();
01024 edm::DetSet<SiStripRawDigi>::const_iterator digiRangeStart = digiRangeIterator;
01025
01026
01027 for (; digiRangeIterator != detSetIter->data.end(); ++digiRangeIterator ) {
01028 const SiStripRawDigi& digi = *digiRangeIterator;
01029 const int channel = distance( digiRangeStart, digiRangeIterator );
01030 if ( channel >= 0 && channel < 512 ) currentDataProfiles.GetTECEntry( det, ring, beam, disk ).SetValue( channel, digi.adc() );
01031 else throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: raw digi channel: " << channel << " out of range for det: " << detRawId << "." << std::endl;
01032 }
01033
01034 }
01035
01036 else {
01037
01038
01039 edm::DetSetVector<SiStripDigi>::const_iterator detSetIter = theStripDigis->find( detRawId );
01040
01041
01042 if( detSetIter == theStripDigis->end() ) continue;
01043
01044
01045 edm::DetSet<SiStripDigi>::const_iterator digiRangeIterator = detSetIter->data.begin();
01046
01047 for(; digiRangeIterator != detSetIter->data.end(); ++digiRangeIterator ) {
01048 const SiStripDigi& digi = *digiRangeIterator;
01049 if ( digi.strip() < 512 ) currentDataProfiles.GetTECEntry( det, ring, beam, disk ).SetValue( digi.strip(), digi.adc() );
01050 else throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: digi strip: " << digi.strip() << " out of range for det: " << detRawId << "." << std::endl;
01051 }
01052
01053 }
01054
01055
01056 } while( moduleLoop.TECLoop( det, ring, beam, disk ) );
01057
01058
01059
01060
01061
01062
01063 det = 2; beam = 0; pos = 0;
01064 do {
01065
01066
01067 currentDataProfiles.GetTIBTOBEntry( det, beam, pos ).SetAllValuesTo( 0. );
01068
01069
01070 const int detRawId = detectorId.GetTIBTOBEntry( det, beam, pos );
01071
01072 if( isRawDigi ) {
01073
01074
01075 edm::DetSetVector<SiStripRawDigi>::const_iterator detSetIter = theStripRawDigis->find( detRawId );
01076 if( detSetIter == theStripRawDigis->end() ) {
01077 throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: No raw DetSet found for det: " << detRawId << "." << std::endl;
01078 }
01079
01080
01081 edm::DetSet<SiStripRawDigi>::const_iterator digiRangeIterator = detSetIter->data.begin();
01082 edm::DetSet<SiStripRawDigi>::const_iterator digiRangeStart = digiRangeIterator;
01083
01084
01085 for (; digiRangeIterator != detSetIter->data.end(); ++digiRangeIterator ) {
01086 const SiStripRawDigi& digi = *digiRangeIterator;
01087 const int channel = distance( digiRangeStart, digiRangeIterator );
01088 if ( channel >= 0 && channel < 512 ) currentDataProfiles.GetTIBTOBEntry( det, beam, pos ).SetValue( channel, digi.adc() );
01089 else throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: raw digi channel: " << channel << " out of range for det: " << detRawId << "." << std::endl;
01090 }
01091
01092 }
01093
01094 else {
01095
01096
01097 edm::DetSetVector<SiStripDigi>::const_iterator detSetIter = theStripDigis->find( detRawId );
01098
01099
01100 if( detSetIter == theStripDigis->end() ) continue;
01101
01102
01103 edm::DetSet<SiStripDigi>::const_iterator digiRangeIterator = detSetIter->data.begin();
01104
01105 for(; digiRangeIterator != detSetIter->data.end(); ++digiRangeIterator ) {
01106 const SiStripDigi& digi = *digiRangeIterator;
01107 if ( digi.strip() < 512 ) currentDataProfiles.GetTIBTOBEntry( det, beam, pos ).SetValue( digi.strip(), digi.adc() );
01108 else throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: digi strip: " << digi.strip() << " out of range for det: " << detRawId << "." << std::endl;
01109 }
01110
01111 }
01112
01113 } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
01114
01115
01116
01117
01118 det = 0; beam = 0; disk = 0;
01119 do {
01120
01121
01122 currentDataProfiles.GetTEC2TECEntry( det, beam, disk ).SetAllValuesTo( 0. );
01123
01124
01125 const int detRawId = detectorId.GetTEC2TECEntry( det, beam, disk );
01126
01127 if( isRawDigi ) {
01128
01129
01130 edm::DetSetVector<SiStripRawDigi>::const_iterator detSetIter = theStripRawDigis->find( detRawId );
01131 if( detSetIter == theStripRawDigis->end() ) {
01132 throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: No raw DetSet found for det: " << detRawId << "." << std::endl;
01133 }
01134
01135
01136 edm::DetSet<SiStripRawDigi>::const_iterator digiRangeIterator = detSetIter->data.begin();
01137 edm::DetSet<SiStripRawDigi>::const_iterator digiRangeStart = digiRangeIterator;
01138
01139
01140 for (; digiRangeIterator != detSetIter->data.end(); ++digiRangeIterator ) {
01141 const SiStripRawDigi& digi = *digiRangeIterator;
01142 const int channel = distance( digiRangeStart, digiRangeIterator );
01143 if ( channel >= 0 && channel < 512 ) currentDataProfiles.GetTEC2TECEntry( det, beam, disk ).SetValue( channel, digi.adc() );
01144 else throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: raw digi channel: " << channel << " out of range for det: " << detRawId << "." << std::endl;
01145 }
01146
01147 }
01148
01149 else {
01150
01151
01152 edm::DetSetVector<SiStripDigi>::const_iterator detSetIter = theStripDigis->find( detRawId );
01153
01154
01155 if( detSetIter == theStripDigis->end() ) continue;
01156
01157
01158 edm::DetSet<SiStripDigi>::const_iterator digiRangeIterator = detSetIter->data.begin();
01159
01160 for(; digiRangeIterator != detSetIter->data.end(); ++digiRangeIterator ) {
01161 const SiStripDigi& digi = *digiRangeIterator;
01162 if ( digi.strip() < 512 ) currentDataProfiles.GetTEC2TECEntry( det, beam, disk ).SetValue( digi.strip(), digi.adc() );
01163 else throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: digi strip: " << digi.strip() << " out of range for det: " << detRawId << "." << std::endl;
01164 }
01165
01166 }
01167
01168 } while( moduleLoop.TEC2TECLoop( det, beam, disk ) );
01169
01170 }
01171
01172 }
01173
01174
01175
01176
01177
01186 void LaserAlignment::fillPedestalProfiles( edm::ESHandle<SiStripPedestals>& pedestalsHandle ) {
01187
01188 int det, ring, beam, disk, pos;
01189
01190
01191 det = 0; ring = 0; beam = 0; disk = 0;
01192 do {
01193 SiStripPedestals::Range pedRange = pedestalsHandle->getRange( detectorId.GetTECEntry( det, ring, beam, disk ) );
01194 for( int strip = 0; strip < 512; ++strip ) {
01195 int thePedestal = int( pedestalsHandle->getPed( strip, pedRange ) );
01196 if( thePedestal > 895 ) thePedestal -= 1024;
01197 pedestalProfiles.GetTECEntry( det, ring, beam, disk ).SetValue( strip, thePedestal );
01198 }
01199 } while ( moduleLoop.TECLoop( det, ring, beam, disk ) );
01200
01201
01202
01203 det = 2; beam = 0; pos = 0;
01204 do {
01205 SiStripPedestals::Range pedRange = pedestalsHandle->getRange( detectorId.GetTIBTOBEntry( det, beam, pos ) );
01206 for( int strip = 0; strip < 512; ++strip ) {
01207 int thePedestal = int( pedestalsHandle->getPed( strip, pedRange ) );
01208 if( thePedestal > 895 ) thePedestal -= 1024;
01209 pedestalProfiles.GetTIBTOBEntry( det, beam, pos ).SetValue( strip, thePedestal );
01210 }
01211 } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
01212
01213
01214
01215 det = 0; beam = 0; disk = 0;
01216 do {
01217 SiStripPedestals::Range pedRange = pedestalsHandle->getRange( detectorId.GetTEC2TECEntry( det, beam, disk ) );
01218 for( int strip = 0; strip < 512; ++strip ) {
01219 int thePedestal = int( pedestalsHandle->getPed( strip, pedRange ) );
01220 if( thePedestal > 895 ) thePedestal -= 1024;
01221 pedestalProfiles.GetTEC2TECEntry( det, beam, disk ).SetValue( strip, thePedestal );
01222 }
01223 } while( moduleLoop.TEC2TECLoop( det, beam, disk ) );
01224
01225 }
01226
01227
01228
01229
01230
01236 bool LaserAlignment::isTECBeam( void ) {
01237
01238 int numberOfProfiles = 0;
01239
01240 int ring = 1;
01241 for( int det = 0; det < 2; ++det ) {
01242 for( int beam = 0; beam < 8; ++ beam ) {
01243 for( int disk = 0; disk < 9; ++disk ) {
01244 if( isAcceptedProfile.GetTECEntry( det, ring, beam, disk ) == 1 ) numberOfProfiles++;
01245 }
01246 }
01247 }
01248
01249 LogDebug( "[LaserAlignment::isTECBeam]" ) << " Found: " << numberOfProfiles << "hits." << std::endl;
01250 std::cout << " [LaserAlignment::isTECBeam] -- Found: " << numberOfProfiles << " hits." << std::endl;
01251
01252 if( numberOfProfiles > 10 ) return( true );
01253 return( false );
01254
01255 }
01256
01257
01258
01259
01260
01266
01267 bool LaserAlignment::isATBeam( void ) {
01268
01269 int numberOfProfiles = 0;
01270
01271 int det = 2; int beam = 0; int pos = 0;
01272 do {
01273 if( isAcceptedProfile.GetTIBTOBEntry( det, beam, pos ) == 1 ) numberOfProfiles++;
01274 } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
01275
01276 LogDebug( "[LaserAlignment::isATBeam]" ) << " Found: " << numberOfProfiles << "hits." << std::endl;
01277 std::cout << " [LaserAlignment::isATBeam] -- Found: " << numberOfProfiles << " hits." << std::endl;
01278
01279 if( numberOfProfiles > 10 ) return( true );
01280 return( false );
01281
01282 }
01283
01284
01285
01286
01287
01296 double LaserAlignment::getTIBTOBNominalBeamOffset( unsigned int det, unsigned int beam, unsigned int pos ) {
01297
01298 if( det < 2 || det > 3 || beam > 7 || pos > 5 ) {
01299 throw cms::Exception( "[LaserAlignment::getTIBTOBNominalBeamOffset]" ) << " ERROR ** Called with nonexisting parameter set: det " << det << " beam " << beam << " pos " << pos << "." << std::endl;
01300 }
01301
01302 const double nominalOffsetsTIB[8] = { 0.00035, 2.10687, -2.10827, -0.00173446, 2.10072, -0.00135114, 2.10105, -2.10401 };
01303
01304
01305
01306
01307 const int orientationPattern[6] = { -1, 1, 1, -1, -1, 1 };
01308 const double nominalOffsetsTOB[8] = { 0.00217408, 1.58678, 117.733, 119.321, 120.906, 119.328, 117.743, 1.58947 };
01309
01310
01311 if( det == 2 ) return( -1. * nominalOffsetsTIB[beam] );
01312
01313 else {
01314 if( beam == 0 or beam > 4 ) return( nominalOffsetsTOB[beam] * orientationPattern[pos] );
01315 else return( -1. * nominalOffsetsTOB[beam] * orientationPattern[pos] );
01316 }
01317
01318 }
01319
01320
01321
01322
01331 double LaserAlignment::getTEC2TECNominalBeamOffset( unsigned int det, unsigned int beam, unsigned int disk ) {
01332
01333 if( det > 1 || beam > 7 || disk > 5 ) {
01334 throw cms::Exception( "[LaserAlignment::getTEC2TECNominalBeamOffset]" ) << " ERROR ** Called with nonexisting parameter set: det " << det << " beam " << beam << " disk " << disk << "." << std::endl;
01335 }
01336
01337 const double nominalOffsets[8] = { 0., 2.220, -2.221, 0., 2.214, 0., 2.214, -2.217 };
01338
01339 if( det == 0 ) return -1. * nominalOffsets[beam];
01340 else return nominalOffsets[beam];
01341
01342 }
01343
01344
01345
01346
01347
01351 void LaserAlignment::CalculateNominalCoordinates( void ) {
01352
01353
01354
01355
01356
01357
01358 const double tecPhiPositions[8] = { 0.392699, 1.178097, 1.963495, 2.748894, 3.534292, 4.319690, 5.105088, 5.890486 };
01359 const double atPhiPositions[8] = { 0.392699, 1.289799, 1.851794, 2.748894, 3.645995, 4.319690, 5.216791, 5.778784 };
01360
01361
01362 const double tobRPosition = 600.;
01363 const double tibRPosition = 514.;
01364 const double tecRPosition[2] = { 564., 840. };
01365
01366
01367 const double tobZPosition[6] = { 1040., 580., 220., -140., -500., -860. };
01368 const double tibZPosition[6] = { 620., 380., 180., -100., -340., -540. };
01369
01370
01371 const double tecZPosition[9] = { 1322.5, 1462.5, 1602.5, 1742.5, 1882.5, 2057.5, 2247.5, 2452.5, 2667.5 };
01372
01373
01374
01375
01376
01377
01378
01379
01380 LASGlobalLoop moduleLoop;
01381 int det, ring, beam, disk, pos;
01382
01383
01384
01385 det = 0; ring = 0, beam = 0; disk = 0;
01386 do {
01387
01388 if( det == 0 ) {
01389 nominalCoordinates.SetTECEntry( det, ring, beam, disk, LASCoordinateSet( tecPhiPositions[beam], 0., tecRPosition[ring], 0., tecZPosition[disk], 0. ) );
01390 }
01391 else {
01392 nominalCoordinates.SetTECEntry( det, ring, beam, disk, LASCoordinateSet( tecPhiPositions[beam], 0., tecRPosition[ring], 0., -1. * tecZPosition[disk], 0. ) );
01393 }
01394
01395 } while( moduleLoop.TECLoop( det, ring, beam, disk ) );
01396
01397
01398
01399
01400 det = 2; beam = 0; pos = 0;
01401 do {
01402 if( det == 2 ) {
01403 nominalCoordinates.SetTIBTOBEntry( det, beam, pos, LASCoordinateSet( atPhiPositions[beam], 0., tibRPosition, 0., tibZPosition[pos], 0. ) );
01404 }
01405 else {
01406 nominalCoordinates.SetTIBTOBEntry( det, beam, pos, LASCoordinateSet( atPhiPositions[beam], 0., tobRPosition, 0., tobZPosition[pos], 0. ) );
01407 }
01408
01409 } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
01410
01411
01412
01413
01414
01415 det = 0; beam = 0; disk = 0;
01416 do {
01417
01418 if( det == 0 ) {
01419 nominalCoordinates.SetTEC2TECEntry( det, beam, disk, LASCoordinateSet( atPhiPositions[beam], 0., tecRPosition[0], 0., tecZPosition[disk], 0. ) );
01420 }
01421 else {
01422 nominalCoordinates.SetTEC2TECEntry( det, beam, disk, LASCoordinateSet( atPhiPositions[beam], 0., tecRPosition[0], 0., -1. * tecZPosition[disk], 0. ) );
01423 }
01424
01425 } while( moduleLoop.TEC2TECLoop( det, beam, disk ) );
01426
01427
01428 }
01429
01430
01431
01432
01433
01438 double LaserAlignment::ConvertAngle( double angle ) {
01439
01440 if( angle < -1. * M_PI || angle > M_PI ) {
01441 throw cms::Exception(" [LaserAlignment::ConvertAngle] ") << "** ERROR: Called with illegal input angle: " << angle << "." << std::endl;
01442 }
01443
01444 if( angle >= 0. ) return angle;
01445 else return( angle + 2. * M_PI );
01446
01447 }
01448
01449
01450
01451
01452
01456 void LaserAlignment::DumpPosFileSet( LASGlobalData<LASCoordinateSet>& coordinates ) {
01457
01458 LASGlobalLoop loop;
01459 int det, ring, beam, disk, pos;
01460
01461 std:: cout << std::endl << " [LaserAlignment::DumpPosFileSet] -- Dump: " << std::endl;
01462
01463
01464 det = 0; ring = 0; beam = 0; disk = 0;
01465 do {
01466 std::cout << "POS " << det << "\t" << beam << "\t" << disk << "\t" << ring << "\t" << coordinates.GetTECEntry( det, ring, beam, disk ).GetPhi() << "\t" << coordinates.GetTECEntry( det, ring, beam, disk ).GetPhiError() << std::endl;
01467 } while ( loop.TECLoop( det, ring, beam, disk ) );
01468
01469
01470 det = 2; beam = 0; pos = 0;
01471 do {
01472 std::cout << "POS " << det << "\t" << beam << "\t" << pos << "\t" << "-1" << "\t" << coordinates.GetTIBTOBEntry( det, beam, pos ).GetPhi() << "\t" << coordinates.GetTIBTOBEntry( det, beam, pos ).GetPhiError() << std::endl;
01473 } while( loop.TIBTOBLoop( det, beam, pos ) );
01474
01475
01476 det = 0; beam = 0; disk = 0;
01477 do {
01478 std::cout << "POS " << det << "\t" << beam << "\t" << disk << "\t" << "-1" << "\t" << coordinates.GetTEC2TECEntry( det, beam, disk ).GetPhi() << "\t" << coordinates.GetTEC2TECEntry( det, beam, disk ).GetPhiError() << std::endl;
01479 } while( loop.TEC2TECLoop( det, beam, disk ) );
01480
01481 std:: cout << std::endl << " [LaserAlignment::DumpPosFileSet] -- End dump: " << std::endl;
01482
01483 }
01484
01485
01486
01487
01488
01492 void LaserAlignment::DumpStripFileSet( LASGlobalData<std::pair<float,float> >& measuredStripPositions ) {
01493
01494 LASGlobalLoop loop;
01495 int det, ring, beam, disk, pos;
01496
01497 std:: cout << std::endl << " [LaserAlignment::DumpStripFileSet] -- Dump: " << std::endl;
01498
01499
01500 det = 0; ring = 0; beam = 0; disk = 0;
01501 do {
01502 std::cout << "STRIP " << det << "\t" << beam << "\t" << disk << "\t" << ring << "\t" << measuredStripPositions.GetTECEntry( det, ring, beam, disk ).first
01503 << "\t" << measuredStripPositions.GetTECEntry( det, ring, beam, disk ).second << std::endl;
01504 } while ( loop.TECLoop( det, ring, beam, disk ) );
01505
01506
01507 det = 2; beam = 0; pos = 0;
01508 do {
01509 std::cout << "STRIP " << det << "\t" << beam << "\t" << pos << "\t" << "-1" << "\t" << measuredStripPositions.GetTIBTOBEntry( det, beam, pos ).first
01510 << "\t" << measuredStripPositions.GetTIBTOBEntry( det, beam, pos ).second << std::endl;
01511 } while( loop.TIBTOBLoop( det, beam, pos ) );
01512
01513
01514 det = 0; beam = 0; disk = 0;
01515 do {
01516 std::cout << "STRIP " << det << "\t" << beam << "\t" << disk << "\t" << "-1" << "\t" << measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).first
01517 << "\t" << measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).second << std::endl;
01518 } while( loop.TEC2TECLoop( det, beam, disk ) );
01519
01520 std:: cout << std::endl << " [LaserAlignment::DumpStripFileSet] -- End dump: " << std::endl;
01521
01522
01523 }
01524
01525
01526
01527
01528
01532 void LaserAlignment::DumpHitmaps( LASGlobalData<int> &numberOfAcceptedProfiles ) {
01533
01534 std::cout << " [LaserAlignment::DumpHitmaps] -- Dumping hitmap for TEC+:" << std::endl;
01535 std::cout << " [LaserAlignment::DumpHitmaps] -- Ring4:" << std::endl;
01536 std::cout << " disk0 disk1 disk2 disk3 disk4 disk5 disk6 disk7 disk8" << std::endl;
01537
01538 for( int beam = 0; beam < 8; ++beam ) {
01539 std::cout << " beam" << beam << ":";
01540 for( int disk = 0; disk < 9; ++disk ) {
01541 std::cout << "\t" << numberOfAcceptedProfiles.GetTECEntry( 0, 0, beam, disk );
01542 }
01543 std::cout << std::endl;
01544 }
01545
01546 std::cout << " [LaserAlignment::DumpHitmaps] -- Ring6:" << std::endl;
01547 std::cout << " disk0 disk1 disk2 disk3 disk4 disk5 disk6 disk7 disk8" << std::endl;
01548
01549 for( int beam = 0; beam < 8; ++beam ) {
01550 std::cout << " beam" << beam << ":";
01551 for( int disk = 0; disk < 9; ++disk ) {
01552 std::cout << "\t" << numberOfAcceptedProfiles.GetTECEntry( 0, 1, beam, disk );
01553 }
01554 std::cout << std::endl;
01555 }
01556
01557 std::cout << " [LaserAlignment::DumpHitmaps] -- Dumping hitmap for TEC-:" << std::endl;
01558 std::cout << " [LaserAlignment::DumpHitmaps] -- Ring4:" << std::endl;
01559 std::cout << " disk0 disk1 disk2 disk3 disk4 disk5 disk6 disk7 disk8" << std::endl;
01560
01561 for( int beam = 0; beam < 8; ++beam ) {
01562 std::cout << " beam" << beam << ":";
01563 for( int disk = 0; disk < 9; ++disk ) {
01564 std::cout << "\t" << numberOfAcceptedProfiles.GetTECEntry( 1, 0, beam, disk );
01565 }
01566 std::cout << std::endl;
01567 }
01568
01569 std::cout << " [LaserAlignment::DumpHitmaps] -- Ring6:" << std::endl;
01570 std::cout << " disk0 disk1 disk2 disk3 disk4 disk5 disk6 disk7 disk8" << std::endl;
01571
01572 for( int beam = 0; beam < 8; ++beam ) {
01573 std::cout << " beam" << beam << ":";
01574 for( int disk = 0; disk < 9; ++disk ) {
01575 std::cout << "\t" << numberOfAcceptedProfiles.GetTECEntry( 1, 1, beam, disk );
01576 }
01577 std::cout << std::endl;
01578 }
01579
01580 std::cout << " [LaserAlignment::DumpHitmaps] -- End of dump." << std::endl << std::endl;
01581
01582 }
01583
01584
01585
01586
01587
01592 void LaserAlignment::ApplyEndcapMaskingCorrections( LASGlobalData<LASCoordinateSet>& measuredCoordinates, LASGlobalData<LASCoordinateSet>& nominalCoordinates, LASEndcapAlignmentParameterSet& endcapParameters ) {
01593
01594
01595 for( std::vector<unsigned int>::iterator moduleIt = theMaskTecModules.begin(); moduleIt != theMaskTecModules.end(); ++moduleIt ) {
01596
01597
01598 LASGlobalLoop moduleLoop;
01599 int det, ring, beam, disk;
01600
01601
01602 LASEndcapAlgorithm endcapAlgorithm;
01603
01604
01605 det = 0; ring = 0; beam = 0; disk = 0;
01606 do {
01607
01608
01609 if( detectorId.GetTECEntry( det, ring, beam, disk ) == *moduleIt ) {
01610
01611
01612 const double nominalPhi = nominalCoordinates.GetTECEntry( det, ring, beam, disk ).GetPhi();
01613
01614
01615 const double phiCorrection = endcapAlgorithm.GetAlignmentParameterCorrection( det, ring, beam, disk, nominalCoordinates, endcapParameters );
01616
01617
01618 measuredCoordinates.GetTECEntry( det, ring, beam, disk ).SetPhi( nominalPhi - phiCorrection );
01619
01620 }
01621
01622 } while ( moduleLoop.TECLoop( det, ring, beam, disk ) );
01623
01624 }
01625
01626 }
01627
01628
01629
01630
01631
01636 void LaserAlignment::ApplyATMaskingCorrections( LASGlobalData<LASCoordinateSet>& measuredCoordinates, LASGlobalData<LASCoordinateSet>& nominalCoordinates, LASBarrelAlignmentParameterSet& atParameters ) {
01637
01638
01639 for( std::vector<unsigned int>::iterator moduleIt = theMaskAtModules.begin(); moduleIt != theMaskAtModules.end(); ++moduleIt ) {
01640
01641
01642 LASGlobalLoop moduleLoop;
01643 int det, beam, disk, pos;
01644
01645
01646 LASAlignmentTubeAlgorithm atAlgorithm;
01647
01648
01649
01650
01651
01652 det = 2; beam = 0; pos = 0;
01653 do {
01654
01655
01656 if( detectorId.GetTIBTOBEntry( det, beam, pos ) == *moduleIt ) {
01657
01658
01659 const double nominalPhi = nominalCoordinates.GetTIBTOBEntry( det, beam, pos ).GetPhi();
01660
01661
01662 const double phiCorrection = atAlgorithm.GetTIBTOBAlignmentParameterCorrection( det, beam, pos, nominalCoordinates, atParameters );
01663
01664
01665 measuredCoordinates.GetTIBTOBEntry( det, beam, pos ).SetPhi( nominalPhi - phiCorrection );
01666
01667 }
01668
01669 } while ( moduleLoop.TIBTOBLoop( det, beam, pos ) );
01670
01671
01672
01673
01674 det = 0; beam = 0; disk = 0;
01675 do {
01676
01677
01678 if( detectorId.GetTEC2TECEntry( det, beam, disk ) == *moduleIt ) {
01679
01680
01681 const double nominalPhi = nominalCoordinates.GetTEC2TECEntry( det, beam, disk ).GetPhi();
01682
01683
01684 const double phiCorrection = atAlgorithm.GetTEC2TECAlignmentParameterCorrection( det, beam, disk, nominalCoordinates, atParameters );
01685
01686
01687 measuredCoordinates.GetTEC2TECEntry( det, beam, disk ).SetPhi( nominalPhi - phiCorrection );
01688
01689 }
01690
01691 } while ( moduleLoop.TEC2TECLoop( det, beam, disk ) );
01692
01693 }
01694
01695 }
01696
01697
01698
01699
01700
01705 void LaserAlignment::testRoutine( void ) {
01706
01707
01708
01709 const TrackerGeometry& theTracker( *theTrackerGeometry );
01710
01711 const double atPhiPositions[8] = { 0.392699, 1.289799, 1.851794, 2.748894, 3.645995, 4.319690, 5.216791, 5.778784 };
01712 const double tecPhiPositions[8] = { 0.392699, 1.178097, 1.963495, 2.748894, 3.534292, 4.319690, 5.105088, 5.890486 };
01713 const double zPositions[9] = { 125.0, 139.0, 153.0, 167.0, 181.0, 198.5, 217.5, 238.0, 259.5 };
01714 const double zPositionsTIB[6] = { 62.0, 38.0, 18.0, -10.0, -34.0, -54.0 };
01715 const double zPositionsTOB[6] = { 104.0, 58.0, 22.0, -14.0, -50.0, -86.0 };
01716
01717 int det, beam, disk, pos, ring;
01718
01719
01720 det = 0; ring = 0; beam = 0; disk = 0;
01721 do {
01722
01723 const double radius = ring?84.0:56.4;
01724
01725
01726 const DetId theDetId( detectorId.GetTECEntry( det, ring, beam, disk ) );
01727 const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
01728
01729 if (theStripDet) {
01730 const GlobalPoint gp( GlobalPoint::Cylindrical( radius, tecPhiPositions[beam], zPositions[disk] ) );
01731
01732 const LocalPoint lp( theStripDet->surface().toLocal( gp ) );
01733 std::cout << "__TEC: " << 256. - theStripDet->specificTopology().strip( lp ) << std::endl;
01734 }
01735
01736 } while( moduleLoop.TECLoop( det, ring, beam, disk ) );
01737
01738
01739
01740 det = 2; beam = 0; pos = 0;
01741 do {
01742
01743 const double radius = (det==2?51.4:58.4);
01744 const double theZ = (det==2?zPositionsTIB[pos]:zPositionsTOB[pos]);
01745
01746
01747 const DetId theDetId( detectorId.GetTIBTOBEntry( det, beam, pos ) );
01748 const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
01749
01750 if (theStripDet) {
01751 const GlobalPoint gp( GlobalPoint::Cylindrical( radius, atPhiPositions[beam], theZ ) );
01752
01753 const LocalPoint lp( theStripDet->surface().toLocal( gp ) );
01754 std::cout << "__TIBTOB det " << det << " beam " << beam << " pos " << pos << " " << 256. - theStripDet->specificTopology().strip( lp );
01755 std::cout << " " << theStripDet->position().perp()<< std::endl;
01756 }
01757
01758 } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
01759
01760
01761
01762 det = 0; beam = 0; disk = 0;
01763 do {
01764
01765 const double radius = 56.4;
01766
01767
01768 const DetId theDetId( detectorId.GetTEC2TECEntry( det, beam, disk ) );
01769 const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
01770
01771 if (theStripDet) {
01772 const GlobalPoint gp( GlobalPoint::Cylindrical( radius, atPhiPositions[beam], zPositions[disk] ) );
01773
01774 const LocalPoint lp( theStripDet->surface().toLocal( gp ) );
01775 std::cout << "__TEC2TEC det " << det << " beam " << beam << " disk " << disk << " " << 256. - theStripDet->specificTopology().strip( lp ) << std::endl;
01776 }
01777
01778 } while( moduleLoop.TEC2TECLoop( det, beam, disk ) );
01779
01780
01781 }
01782
01783
01784
01785
01786
01787
01788
01789
01790
01791
01792 #include "FWCore/Framework/interface/MakerMacros.h"
01793
01794 DEFINE_FWK_MODULE(LaserAlignment);
01795
01796
01797
01798
01799
01800