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