CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/Alignment/LaserAlignment/plugins/LaserAlignment.cc

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