CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/Alignment/SurveyAnalysis/plugins/SurveyDataConverter.cc

Go to the documentation of this file.
00001 
00002 // Framework
00003 #include "FWCore/Framework/interface/MakerMacros.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 
00006 #include "Alignment/CommonAlignment/interface/AlignableModifier.h"
00007 #include "Alignment/TrackerAlignment/interface/TrackerAlignment.h"
00008 #include "Alignment/TrackerAlignment/interface/TrackerScenarioBuilder.h"
00009 #include "Alignment/CommonAlignment/interface/Alignable.h" 
00010 
00011 #include "Alignment/SurveyAnalysis/plugins/SurveyDataConverter.h"
00012 
00013 
00014 //__________________________________________________________________________________________________
00015 SurveyDataConverter::SurveyDataConverter(const edm::ParameterSet& iConfig) :
00016   theParameterSet( iConfig )
00017 {  
00018 }
00019 
00020 //__________________________________________________________________________________________________
00021 void SurveyDataConverter::analyze( const edm::Event& iEvent, const edm::EventSetup& iSetup )
00022 {
00023   
00024   edm::LogInfo("SurveyDataConverter") << "Analyzer called";
00025   applyfineinfo = theParameterSet.getParameter<bool>("applyFineInfo");
00026   applycoarseinfo = theParameterSet.getParameter<bool>("applyCoarseInfo");
00027   adderrors = theParameterSet.getParameter<bool>("applyErrors");
00028 
00029   // Read in the information from the text files
00030   edm::ParameterSet textFiles = theParameterSet.getParameter<edm::ParameterSet>( "textFileNames" );
00031 
00032   std::string textFileNames[NFILES]; 
00033   std::string fileType[NFILES];    
00034   textFileNames[0] = textFiles.getUntrackedParameter<std::string>("forTIB","NONE");  
00035   fileType[0] = "TIB";
00036   textFileNames[1] = textFiles.getUntrackedParameter<std::string>("forTID","NONE");
00037   fileType[1] = "TID";
00038 
00039   SurveyDataReader dataReader;
00040   for (int ii=0 ; ii<NFILES ;ii++) {
00041     if ( textFileNames[ii] == "NONE" )
00042       throw cms::Exception("BadConfig") << fileType[ii] << " input file not found in configuration";
00043     dataReader.readFile( textFileNames[ii], fileType[ii] );
00044   }
00045 
00046   // Get info and map
00047   const MapType& mapIdToInfo = dataReader.detIdMap();
00048   std::cout << "DATA HAS BEEN READ INTO THE MAP" << std::endl;
00049   std::cout << "DATA HAS BEEN CONVERTED IN ALIGNABLE COORDINATES" << std::endl;  
00050 
00051   TrackerAlignment tr_align( iSetup );
00052   if (applycoarseinfo) this->applyCoarseSurveyInfo(tr_align); 
00053   if (applyfineinfo) this->applyFineSurveyInfo(tr_align, mapIdToInfo);
00054   if (adderrors) this->applyAPEs(tr_align);
00055   tr_align.saveToDB();
00056 }
00057 
00058 //___________________________________
00059 //
00060 void SurveyDataConverter::applyFineSurveyInfo( TrackerAlignment& tr_align, const MapType& map ){
00061 
00062   std::cout << "Apply fine info: " << std::endl;
00063         
00064   for ( MapType::const_iterator it = map.begin(); it != map.end(); it++){
00065 
00066     const align::Scalars& align_params = (it)->second; 
00067       
00068     align::Scalars translations; 
00069     translations.push_back(align_params[0]);  
00070     translations.push_back(align_params[1]);  
00071     translations.push_back(align_params[2]); 
00072 
00073     align::RotationType bRotation(align_params[6], align_params[9], align_params[3],
00074                                   align_params[7], align_params[10], align_params[4],
00075                                   align_params[8], align_params[11], align_params[5]);
00076 
00077     align::RotationType fRotation(align_params[15], align_params[18], align_params[12],
00078                                   align_params[16], align_params[19], align_params[13],
00079                                   align_params[17], align_params[20], align_params[14]);
00080 
00081     // Use "false" for debugging only
00082     tr_align.moveAlignableTIBTIDs((it)->first, translations, bRotation, fRotation, true);
00083   }
00084 }
00085 
00086 //___________________________________
00087 //
00088 void SurveyDataConverter::applyCoarseSurveyInfo( TrackerAlignment& tr_align ){
00089         
00090   std::cout << "Apply coarse info: " << std::endl;
00091   MisalignScenario = theParameterSet.getParameter<edm::ParameterSet>( "MisalignmentScenario" );
00092 
00093   TrackerScenarioBuilder scenarioBuilder( tr_align.getAlignableTracker() );
00094   scenarioBuilder.applyScenario( MisalignScenario );
00095   
00096 }
00097 
00098 //___________________________________
00099 //
00100 void SurveyDataConverter::applyAPEs( TrackerAlignment& tr_align ) {
00101         
00102   std::cout << "Apply APEs: " << std::endl;
00103   // Neglect sensor-on-module mounting precision (10 um)
00104   // Irrelevant given other sizes ..
00105   std::vector<double> TIBerrors = theParameterSet.getParameter< std::vector<double> >("TIBerrors");
00106   std::vector<double> TOBerrors = theParameterSet.getParameter< std::vector<double> >("TOBerrors");
00107   std::vector<double> TIDerrors = theParameterSet.getParameter< std::vector<double> >("TIDerrors"); 
00108   std::vector<double> TECerrors = theParameterSet.getParameter< std::vector<double> >("TECerrors"); 
00109         
00110   if (TIBerrors.size() < 3 || TOBerrors.size() < 4 || TIDerrors.size() < 4 || TECerrors.size() < 4) {
00111     std::cout << "APE info not valid : please check test/run-converter.cfg" << std::endl;
00112     return;
00113   }
00114          
00115   AlignableModifier* theModifier = new AlignableModifier();
00116   AlignableTracker* theAlignableTracker = tr_align.getAlignableTracker() ; 
00117   align::Alignables::const_iterator iter;
00118 
00119   // TIB
00120   const align::Alignables& theTIBhb = theAlignableTracker->innerHalfBarrels();
00121   for (iter = theTIBhb.begin(); iter != theTIBhb.end(); ++iter ) 
00122     { theModifier->addAlignmentPositionErrorLocal( *iter, TIBerrors.at(0), 
00123                                                    TIBerrors.at(0), TIBerrors.at(0) ); }
00124   const align::Alignables& theTIBlayers = theAlignableTracker->innerBarrelLayers();
00125   for (iter = theTIBlayers.begin(); iter != theTIBlayers.end(); ++iter)
00126     { theModifier->addAlignmentPositionErrorLocal( *iter, TIBerrors.at(1), 
00127                                                    TIBerrors.at(1), TIBerrors.at(1) ); }
00128   const align::Alignables& theTIBgd = theAlignableTracker->innerBarrelGeomDets();
00129   for (iter = theTIBgd.begin(); iter != theTIBgd.end(); ++iter ) 
00130     { theModifier->addAlignmentPositionErrorLocal( *iter, TIBerrors.at(2), 
00131                                                    TIBerrors.at(2), TIBerrors.at(2) ); }
00132 
00133   // TOB
00134   const align::Alignables& theTOBhb = theAlignableTracker->outerHalfBarrels();
00135   for (iter = theTOBhb.begin(); iter != theTOBhb.end(); ++iter )
00136     { theModifier->addAlignmentPositionErrorLocal( *iter, TOBerrors.at(0), 
00137                                                    TOBerrors.at(0), TOBerrors.at(1) ); }
00138   const align::Alignables& theTOBrods = theAlignableTracker->outerBarrelRods();
00139   for (iter = theTOBrods.begin(); iter != theTOBrods.end(); ++iter ) 
00140     { theModifier->addAlignmentPositionErrorLocal( *iter, TOBerrors.at(2), 
00141                                                    TOBerrors.at(2), TOBerrors.at(2) ); }
00142   const align::Alignables& theTOBgd = theAlignableTracker->outerBarrelGeomDets();
00143   for (iter = theTOBgd.begin(); iter != theTOBgd.end(); ++iter )
00144     { theModifier->addAlignmentPositionErrorLocal( *iter, TOBerrors.at(3), 
00145                                                    TOBerrors.at(3), TOBerrors.at(3) ); }
00146 
00147   // TID
00148   const align::Alignables& theTIDs = theAlignableTracker->TIDs();
00149   for (iter = theTIDs.begin(); iter != theTIDs.end(); ++iter ) 
00150     { theModifier->addAlignmentPositionErrorLocal( *iter, TIDerrors.at(0), 
00151                                                    TIDerrors.at(0), TIDerrors.at(0) ); }
00152   const align::Alignables& theTIDdiscs = theAlignableTracker->TIDLayers();
00153   for (iter = theTIDdiscs.begin(); iter != theTIDdiscs.end(); ++iter )
00154     { theModifier->addAlignmentPositionErrorLocal( *iter, TIDerrors.at(1), 
00155                                                    TIDerrors.at(1), TIDerrors.at(1) ); }
00156   const align::Alignables& theTIDrings = theAlignableTracker->TIDRings();
00157   for (iter = theTIDrings.begin(); iter != theTIDrings.end(); ++iter )
00158     { theModifier->addAlignmentPositionErrorLocal( *iter, TIDerrors.at(2), 
00159                                                    TIDerrors.at(2), TIDerrors.at(2) ); } 
00160   const align::Alignables& theTIDgd = theAlignableTracker->TIDGeomDets();
00161   for (iter = theTIDgd.begin(); iter != theTIDgd.end(); ++iter )
00162     { theModifier->addAlignmentPositionErrorLocal( *iter, TIDerrors.at(3), 
00163                                                    TIDerrors.at(3), TIDerrors.at(3) ); } 
00164 
00165   // TEC
00166   const align::Alignables& theTECs = theAlignableTracker->endCaps();
00167   for (iter = theTECs.begin(); iter != theTECs.end(); ++iter ) 
00168     { theModifier->addAlignmentPositionErrorLocal( *iter, TECerrors.at(0), 
00169                                                    TECerrors.at(0), TECerrors.at(0) ); } 
00170   const align::Alignables& theTECdiscs = theAlignableTracker->endcapLayers();
00171   for (iter = theTECdiscs.begin(); iter != theTECdiscs.end(); ++iter )
00172     { theModifier->addAlignmentPositionErrorLocal( *iter, TECerrors.at(1), 
00173                                                    TECerrors.at(1), TECerrors.at(1) ); } 
00174   const align::Alignables& theTECpetals = theAlignableTracker->endcapPetals();
00175   for (iter = theTECpetals.begin(); iter != theTECpetals.end(); ++iter ) 
00176     { theModifier->addAlignmentPositionErrorLocal( *iter, TECerrors.at(2), 
00177                                                    TECerrors.at(2), TECerrors.at(2) ); }   
00178   const align::Alignables& theTECgd = theAlignableTracker->endcapGeomDets();
00179   for (iter = theTECgd.begin(); iter != theTECgd.end(); ++iter )
00180     { theModifier->addAlignmentPositionErrorLocal( *iter, TECerrors.at(3), 
00181                                                    TECerrors.at(3), TECerrors.at(3) ); }   
00182 }
00183 
00184 DEFINE_FWK_MODULE(SurveyDataConverter);
00185