#include <Alignment/SurveyDataConverter/interface/SurveyDataConverter.h>
Public Member Functions | |
virtual void | analyze (const edm::Event &iEvent, const edm::EventSetup &iSetup) |
virtual void | endJob () |
SurveyDataConverter (const edm::ParameterSet &iConfig) | |
Private Types | |
typedef SurveyDataReader::MapType | MapType |
typedef SurveyDataReader::MapTypeOr | MapTypeOr |
typedef SurveyDataReader::PairType | PairType |
typedef SurveyDataReader::PairTypeOr | PairTypeOr |
Private Member Functions | |
void | applyAPEs (TrackerAlignment &tr_align) |
void | applyCoarseSurveyInfo (TrackerAlignment &tr_align) |
void | applyFineSurveyInfo (TrackerAlignment &tr_align, const MapType &map) |
Private Attributes | |
bool | adderrors |
bool | applycoarseinfo |
bool | applyfineinfo |
edm::ParameterSet | MisalignScenario |
edm::ParameterSet | theParameterSet |
Static Private Attributes | |
static const int | NFILES = 2 |
Description: Reads survey corrections and applies them to the geometry.
Implementation: <Notes on="" implementation>="">
Definition at line 25 of file SurveyDataConverter.h.
typedef SurveyDataReader::MapType SurveyDataConverter::MapType [private] |
Definition at line 28 of file SurveyDataConverter.h.
typedef SurveyDataReader::MapTypeOr SurveyDataConverter::MapTypeOr [private] |
Definition at line 30 of file SurveyDataConverter.h.
typedef SurveyDataReader::PairType SurveyDataConverter::PairType [private] |
Definition at line 29 of file SurveyDataConverter.h.
typedef SurveyDataReader::PairTypeOr SurveyDataConverter::PairTypeOr [private] |
Definition at line 31 of file SurveyDataConverter.h.
SurveyDataConverter::SurveyDataConverter | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 15 of file SurveyDataConverter.cc.
: theParameterSet( iConfig ) { }
void SurveyDataConverter::analyze | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements edm::EDAnalyzer.
Definition at line 21 of file SurveyDataConverter.cc.
References adderrors, applyAPEs(), applycoarseinfo, applyCoarseSurveyInfo(), applyfineinfo, applyFineSurveyInfo(), gather_cfg::cout, SurveyDataReader::detIdMap(), Exception, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), NFILES, SurveyDataReader::readFile(), TrackerAlignment::saveToDB(), and theParameterSet.
{ edm::LogInfo("SurveyDataConverter") << "Analyzer called"; applyfineinfo = theParameterSet.getParameter<bool>("applyFineInfo"); applycoarseinfo = theParameterSet.getParameter<bool>("applyCoarseInfo"); adderrors = theParameterSet.getParameter<bool>("applyErrors"); // Read in the information from the text files edm::ParameterSet textFiles = theParameterSet.getParameter<edm::ParameterSet>( "textFileNames" ); std::string textFileNames[NFILES]; std::string fileType[NFILES]; textFileNames[0] = textFiles.getUntrackedParameter<std::string>("forTIB","NONE"); fileType[0] = "TIB"; textFileNames[1] = textFiles.getUntrackedParameter<std::string>("forTID","NONE"); fileType[1] = "TID"; SurveyDataReader dataReader; for (int ii=0 ; ii<NFILES ;ii++) { if ( textFileNames[ii] == "NONE" ) throw cms::Exception("BadConfig") << fileType[ii] << " input file not found in configuration"; dataReader.readFile( textFileNames[ii], fileType[ii] ); } // Get info and map const MapType& mapIdToInfo = dataReader.detIdMap(); std::cout << "DATA HAS BEEN READ INTO THE MAP" << std::endl; std::cout << "DATA HAS BEEN CONVERTED IN ALIGNABLE COORDINATES" << std::endl; TrackerAlignment tr_align( iSetup ); if (applycoarseinfo) this->applyCoarseSurveyInfo(tr_align); if (applyfineinfo) this->applyFineSurveyInfo(tr_align, mapIdToInfo); if (adderrors) this->applyAPEs(tr_align); tr_align.saveToDB(); }
void SurveyDataConverter::applyAPEs | ( | TrackerAlignment & | tr_align | ) | [private] |
Definition at line 100 of file SurveyDataConverter.cc.
References AlignableModifier::addAlignmentPositionErrorLocal(), gather_cfg::cout, AlignableTracker::endcapGeomDets(), AlignableTracker::endcapLayers(), AlignableTracker::endcapPetals(), AlignableTracker::endCaps(), TrackerAlignment::getAlignableTracker(), edm::ParameterSet::getParameter(), AlignableTracker::innerBarrelGeomDets(), AlignableTracker::innerBarrelLayers(), AlignableTracker::innerHalfBarrels(), AlignableTracker::outerBarrelGeomDets(), AlignableTracker::outerBarrelRods(), AlignableTracker::outerHalfBarrels(), theParameterSet, AlignableTracker::TIDGeomDets(), AlignableTracker::TIDLayers(), AlignableTracker::TIDRings(), and AlignableTracker::TIDs().
Referenced by analyze().
{ std::cout << "Apply APEs: " << std::endl; // Neglect sensor-on-module mounting precision (10 um) // Irrelevant given other sizes .. std::vector<double> TIBerrors = theParameterSet.getParameter< std::vector<double> >("TIBerrors"); std::vector<double> TOBerrors = theParameterSet.getParameter< std::vector<double> >("TOBerrors"); std::vector<double> TIDerrors = theParameterSet.getParameter< std::vector<double> >("TIDerrors"); std::vector<double> TECerrors = theParameterSet.getParameter< std::vector<double> >("TECerrors"); if (TIBerrors.size() < 3 || TOBerrors.size() < 4 || TIDerrors.size() < 4 || TECerrors.size() < 4) { std::cout << "APE info not valid : please check test/run-converter.cfg" << std::endl; return; } AlignableModifier* theModifier = new AlignableModifier(); AlignableTracker* theAlignableTracker = tr_align.getAlignableTracker() ; align::Alignables::const_iterator iter; // TIB const align::Alignables& theTIBhb = theAlignableTracker->innerHalfBarrels(); for (iter = theTIBhb.begin(); iter != theTIBhb.end(); ++iter ) { theModifier->addAlignmentPositionErrorLocal( *iter, TIBerrors.at(0), TIBerrors.at(0), TIBerrors.at(0) ); } const align::Alignables& theTIBlayers = theAlignableTracker->innerBarrelLayers(); for (iter = theTIBlayers.begin(); iter != theTIBlayers.end(); ++iter) { theModifier->addAlignmentPositionErrorLocal( *iter, TIBerrors.at(1), TIBerrors.at(1), TIBerrors.at(1) ); } const align::Alignables& theTIBgd = theAlignableTracker->innerBarrelGeomDets(); for (iter = theTIBgd.begin(); iter != theTIBgd.end(); ++iter ) { theModifier->addAlignmentPositionErrorLocal( *iter, TIBerrors.at(2), TIBerrors.at(2), TIBerrors.at(2) ); } // TOB const align::Alignables& theTOBhb = theAlignableTracker->outerHalfBarrels(); for (iter = theTOBhb.begin(); iter != theTOBhb.end(); ++iter ) { theModifier->addAlignmentPositionErrorLocal( *iter, TOBerrors.at(0), TOBerrors.at(0), TOBerrors.at(1) ); } const align::Alignables& theTOBrods = theAlignableTracker->outerBarrelRods(); for (iter = theTOBrods.begin(); iter != theTOBrods.end(); ++iter ) { theModifier->addAlignmentPositionErrorLocal( *iter, TOBerrors.at(2), TOBerrors.at(2), TOBerrors.at(2) ); } const align::Alignables& theTOBgd = theAlignableTracker->outerBarrelGeomDets(); for (iter = theTOBgd.begin(); iter != theTOBgd.end(); ++iter ) { theModifier->addAlignmentPositionErrorLocal( *iter, TOBerrors.at(3), TOBerrors.at(3), TOBerrors.at(3) ); } // TID const align::Alignables& theTIDs = theAlignableTracker->TIDs(); for (iter = theTIDs.begin(); iter != theTIDs.end(); ++iter ) { theModifier->addAlignmentPositionErrorLocal( *iter, TIDerrors.at(0), TIDerrors.at(0), TIDerrors.at(0) ); } const align::Alignables& theTIDdiscs = theAlignableTracker->TIDLayers(); for (iter = theTIDdiscs.begin(); iter != theTIDdiscs.end(); ++iter ) { theModifier->addAlignmentPositionErrorLocal( *iter, TIDerrors.at(1), TIDerrors.at(1), TIDerrors.at(1) ); } const align::Alignables& theTIDrings = theAlignableTracker->TIDRings(); for (iter = theTIDrings.begin(); iter != theTIDrings.end(); ++iter ) { theModifier->addAlignmentPositionErrorLocal( *iter, TIDerrors.at(2), TIDerrors.at(2), TIDerrors.at(2) ); } const align::Alignables& theTIDgd = theAlignableTracker->TIDGeomDets(); for (iter = theTIDgd.begin(); iter != theTIDgd.end(); ++iter ) { theModifier->addAlignmentPositionErrorLocal( *iter, TIDerrors.at(3), TIDerrors.at(3), TIDerrors.at(3) ); } // TEC const align::Alignables& theTECs = theAlignableTracker->endCaps(); for (iter = theTECs.begin(); iter != theTECs.end(); ++iter ) { theModifier->addAlignmentPositionErrorLocal( *iter, TECerrors.at(0), TECerrors.at(0), TECerrors.at(0) ); } const align::Alignables& theTECdiscs = theAlignableTracker->endcapLayers(); for (iter = theTECdiscs.begin(); iter != theTECdiscs.end(); ++iter ) { theModifier->addAlignmentPositionErrorLocal( *iter, TECerrors.at(1), TECerrors.at(1), TECerrors.at(1) ); } const align::Alignables& theTECpetals = theAlignableTracker->endcapPetals(); for (iter = theTECpetals.begin(); iter != theTECpetals.end(); ++iter ) { theModifier->addAlignmentPositionErrorLocal( *iter, TECerrors.at(2), TECerrors.at(2), TECerrors.at(2) ); } const align::Alignables& theTECgd = theAlignableTracker->endcapGeomDets(); for (iter = theTECgd.begin(); iter != theTECgd.end(); ++iter ) { theModifier->addAlignmentPositionErrorLocal( *iter, TECerrors.at(3), TECerrors.at(3), TECerrors.at(3) ); } }
void SurveyDataConverter::applyCoarseSurveyInfo | ( | TrackerAlignment & | tr_align | ) | [private] |
Definition at line 88 of file SurveyDataConverter.cc.
References TrackerScenarioBuilder::applyScenario(), gather_cfg::cout, TrackerAlignment::getAlignableTracker(), edm::ParameterSet::getParameter(), MisalignScenario, and theParameterSet.
Referenced by analyze().
{ std::cout << "Apply coarse info: " << std::endl; MisalignScenario = theParameterSet.getParameter<edm::ParameterSet>( "MisalignmentScenario" ); TrackerScenarioBuilder scenarioBuilder( tr_align.getAlignableTracker() ); scenarioBuilder.applyScenario( MisalignScenario ); }
void SurveyDataConverter::applyFineSurveyInfo | ( | TrackerAlignment & | tr_align, |
const MapType & | map | ||
) | [private] |
Definition at line 60 of file SurveyDataConverter.cc.
References gather_cfg::cout, first, TrackerAlignment::moveAlignableTIBTIDs(), and edm::second().
Referenced by analyze().
{ std::cout << "Apply fine info: " << std::endl; for ( MapType::const_iterator it = map.begin(); it != map.end(); it++){ const align::Scalars& align_params = (it)->second; align::Scalars translations; translations.push_back(align_params[0]); translations.push_back(align_params[1]); translations.push_back(align_params[2]); align::RotationType bRotation(align_params[6], align_params[9], align_params[3], align_params[7], align_params[10], align_params[4], align_params[8], align_params[11], align_params[5]); align::RotationType fRotation(align_params[15], align_params[18], align_params[12], align_params[16], align_params[19], align_params[13], align_params[17], align_params[20], align_params[14]); // Use "false" for debugging only tr_align.moveAlignableTIBTIDs((it)->first, translations, bRotation, fRotation, true); } }
virtual void SurveyDataConverter::endJob | ( | void | ) | [inline, virtual] |
bool SurveyDataConverter::adderrors [private] |
Definition at line 57 of file SurveyDataConverter.h.
Referenced by analyze().
bool SurveyDataConverter::applycoarseinfo [private] |
Definition at line 57 of file SurveyDataConverter.h.
Referenced by analyze().
bool SurveyDataConverter::applyfineinfo [private] |
Definition at line 57 of file SurveyDataConverter.h.
Referenced by analyze().
Definition at line 53 of file SurveyDataConverter.h.
Referenced by applyCoarseSurveyInfo().
const int SurveyDataConverter::NFILES = 2 [static, private] |
Definition at line 41 of file SurveyDataConverter.h.
Referenced by analyze().
Definition at line 52 of file SurveyDataConverter.h.
Referenced by analyze(), applyAPEs(), and applyCoarseSurveyInfo().