CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SurveyDataConverter.cc
Go to the documentation of this file.
1 
2 // Framework
6 
11 
13 
16 
17 //__________________________________________________________________________________________________
19  theParameterSet( iConfig )
20 {
21 }
22 
23 //__________________________________________________________________________________________________
25 {
26  //Retrieve tracker topology from geometry
28  iSetup.get<IdealGeometryRecord>().get(tTopoHandle);
29  const TrackerTopology* const tTopo = tTopoHandle.product();
30 
31  edm::LogInfo("SurveyDataConverter") << "Analyzer called";
32  applyfineinfo = theParameterSet.getParameter<bool>("applyFineInfo");
33  applycoarseinfo = theParameterSet.getParameter<bool>("applyCoarseInfo");
34  adderrors = theParameterSet.getParameter<bool>("applyErrors");
35 
36  // Read in the information from the text files
37  edm::ParameterSet textFiles = theParameterSet.getParameter<edm::ParameterSet>( "textFileNames" );
38 
39  std::string textFileNames[NFILES];
40  std::string fileType[NFILES];
41  textFileNames[0] = textFiles.getUntrackedParameter<std::string>("forTIB","NONE");
42  fileType[0] = "TIB";
43  textFileNames[1] = textFiles.getUntrackedParameter<std::string>("forTID","NONE");
44  fileType[1] = "TID";
45 
46  SurveyDataReader dataReader;
47  for (int ii=0 ; ii<NFILES ;ii++) {
48  if ( textFileNames[ii] == "NONE" )
49  throw cms::Exception("BadConfig") << fileType[ii] << " input file not found in configuration";
50  dataReader.readFile( textFileNames[ii], fileType[ii], tTopo );
51  }
52 
53  // Get info and map
54  const MapType& mapIdToInfo = dataReader.detIdMap();
55  std::cout << "DATA HAS BEEN READ INTO THE MAP" << std::endl;
56  std::cout << "DATA HAS BEEN CONVERTED IN ALIGNABLE COORDINATES" << std::endl;
57 
58  TrackerAlignment tr_align( iSetup );
59  if (applycoarseinfo) this->applyCoarseSurveyInfo(tr_align);
60  if (applyfineinfo) this->applyFineSurveyInfo(tr_align, mapIdToInfo);
61  if (adderrors) this->applyAPEs(tr_align);
62  tr_align.saveToDB();
63 }
64 
65 //___________________________________
66 //
68 
69  std::cout << "Apply fine info: " << std::endl;
70 
71  for ( MapType::const_iterator it = map.begin(); it != map.end(); it++){
72 
73  const align::Scalars& align_params = (it)->second;
74 
75  align::Scalars translations;
76  translations.push_back(align_params[0]);
77  translations.push_back(align_params[1]);
78  translations.push_back(align_params[2]);
79 
80  align::RotationType bRotation(align_params[6], align_params[9], align_params[3],
81  align_params[7], align_params[10], align_params[4],
82  align_params[8], align_params[11], align_params[5]);
83 
84  align::RotationType fRotation(align_params[15], align_params[18], align_params[12],
85  align_params[16], align_params[19], align_params[13],
86  align_params[17], align_params[20], align_params[14]);
87 
88  // Use "false" for debugging only
89  tr_align.moveAlignableTIBTIDs((it)->first, translations, bRotation, fRotation, true);
90  }
91 }
92 
93 //___________________________________
94 //
96 
97  std::cout << "Apply coarse info: " << std::endl;
99 
100  TrackerScenarioBuilder scenarioBuilder( tr_align.getAlignableTracker() );
101  scenarioBuilder.applyScenario( MisalignScenario );
102 
103 }
104 
105 //___________________________________
106 //
108 
109  std::cout << "Apply APEs: " << std::endl;
110  // Neglect sensor-on-module mounting precision (10 um)
111  // Irrelevant given other sizes ..
112  std::vector<double> TIBerrors = theParameterSet.getParameter< std::vector<double> >("TIBerrors");
113  std::vector<double> TOBerrors = theParameterSet.getParameter< std::vector<double> >("TOBerrors");
114  std::vector<double> TIDerrors = theParameterSet.getParameter< std::vector<double> >("TIDerrors");
115  std::vector<double> TECerrors = theParameterSet.getParameter< std::vector<double> >("TECerrors");
116 
117  if (TIBerrors.size() < 3 || TOBerrors.size() < 4 || TIDerrors.size() < 4 || TECerrors.size() < 4) {
118  std::cout << "APE info not valid : please check test/run-converter.cfg" << std::endl;
119  return;
120  }
121 
122  AlignableModifier* theModifier = new AlignableModifier();
123  AlignableTracker* theAlignableTracker = tr_align.getAlignableTracker() ;
124  align::Alignables::const_iterator iter;
125 
126  // TIB
127  const align::Alignables& theTIBhb = theAlignableTracker->innerHalfBarrels();
128  for (iter = theTIBhb.begin(); iter != theTIBhb.end(); ++iter )
129  { theModifier->addAlignmentPositionErrorLocal( *iter, TIBerrors.at(0),
130  TIBerrors.at(0), TIBerrors.at(0) ); }
131  const align::Alignables& theTIBlayers = theAlignableTracker->innerBarrelLayers();
132  for (iter = theTIBlayers.begin(); iter != theTIBlayers.end(); ++iter)
133  { theModifier->addAlignmentPositionErrorLocal( *iter, TIBerrors.at(1),
134  TIBerrors.at(1), TIBerrors.at(1) ); }
135  const align::Alignables& theTIBgd = theAlignableTracker->innerBarrelGeomDets();
136  for (iter = theTIBgd.begin(); iter != theTIBgd.end(); ++iter )
137  { theModifier->addAlignmentPositionErrorLocal( *iter, TIBerrors.at(2),
138  TIBerrors.at(2), TIBerrors.at(2) ); }
139 
140  // TOB
141  const align::Alignables& theTOBhb = theAlignableTracker->outerHalfBarrels();
142  for (iter = theTOBhb.begin(); iter != theTOBhb.end(); ++iter )
143  { theModifier->addAlignmentPositionErrorLocal( *iter, TOBerrors.at(0),
144  TOBerrors.at(0), TOBerrors.at(1) ); }
145  const align::Alignables& theTOBrods = theAlignableTracker->outerBarrelRods();
146  for (iter = theTOBrods.begin(); iter != theTOBrods.end(); ++iter )
147  { theModifier->addAlignmentPositionErrorLocal( *iter, TOBerrors.at(2),
148  TOBerrors.at(2), TOBerrors.at(2) ); }
149  const align::Alignables& theTOBgd = theAlignableTracker->outerBarrelGeomDets();
150  for (iter = theTOBgd.begin(); iter != theTOBgd.end(); ++iter )
151  { theModifier->addAlignmentPositionErrorLocal( *iter, TOBerrors.at(3),
152  TOBerrors.at(3), TOBerrors.at(3) ); }
153 
154  // TID
155  const align::Alignables& theTIDs = theAlignableTracker->TIDs();
156  for (iter = theTIDs.begin(); iter != theTIDs.end(); ++iter )
157  { theModifier->addAlignmentPositionErrorLocal( *iter, TIDerrors.at(0),
158  TIDerrors.at(0), TIDerrors.at(0) ); }
159  const align::Alignables& theTIDdiscs = theAlignableTracker->TIDLayers();
160  for (iter = theTIDdiscs.begin(); iter != theTIDdiscs.end(); ++iter )
161  { theModifier->addAlignmentPositionErrorLocal( *iter, TIDerrors.at(1),
162  TIDerrors.at(1), TIDerrors.at(1) ); }
163  const align::Alignables& theTIDrings = theAlignableTracker->TIDRings();
164  for (iter = theTIDrings.begin(); iter != theTIDrings.end(); ++iter )
165  { theModifier->addAlignmentPositionErrorLocal( *iter, TIDerrors.at(2),
166  TIDerrors.at(2), TIDerrors.at(2) ); }
167  const align::Alignables& theTIDgd = theAlignableTracker->TIDGeomDets();
168  for (iter = theTIDgd.begin(); iter != theTIDgd.end(); ++iter )
169  { theModifier->addAlignmentPositionErrorLocal( *iter, TIDerrors.at(3),
170  TIDerrors.at(3), TIDerrors.at(3) ); }
171 
172  // TEC
173  const align::Alignables& theTECs = theAlignableTracker->endCaps();
174  for (iter = theTECs.begin(); iter != theTECs.end(); ++iter )
175  { theModifier->addAlignmentPositionErrorLocal( *iter, TECerrors.at(0),
176  TECerrors.at(0), TECerrors.at(0) ); }
177  const align::Alignables& theTECdiscs = theAlignableTracker->endcapLayers();
178  for (iter = theTECdiscs.begin(); iter != theTECdiscs.end(); ++iter )
179  { theModifier->addAlignmentPositionErrorLocal( *iter, TECerrors.at(1),
180  TECerrors.at(1), TECerrors.at(1) ); }
181  const align::Alignables& theTECpetals = theAlignableTracker->endcapPetals();
182  for (iter = theTECpetals.begin(); iter != theTECpetals.end(); ++iter )
183  { theModifier->addAlignmentPositionErrorLocal( *iter, TECerrors.at(2),
184  TECerrors.at(2), TECerrors.at(2) ); }
185  const align::Alignables& theTECgd = theAlignableTracker->endcapGeomDets();
186  for (iter = theTECgd.begin(); iter != theTECgd.end(); ++iter )
187  { theModifier->addAlignmentPositionErrorLocal( *iter, TECerrors.at(3),
188  TECerrors.at(3), TECerrors.at(3) ); }
189 }
190 
192 
T getParameter(std::string const &) const
SurveyDataReader::MapType MapType
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void applyFineSurveyInfo(TrackerAlignment &tr_align, const MapType &map)
int ii
Definition: cuy.py:588
void addAlignmentPositionErrorLocal(Alignable *alignable, float dx, float dy, float dz)
Add the AlignmentPositionError (in local frame) to Alignable.
dictionary map
Definition: Association.py:205
U second(std::pair< T, U > const &p)
void applyCoarseSurveyInfo(TrackerAlignment &tr_align)
int iEvent
Definition: GenABIO.cc:243
void applyAPEs(TrackerAlignment &tr_align)
edm::ParameterSet MisalignScenario
Alignables & innerBarrelGeomDets()
Return inner barrel GeomDets.
std::vector< Scalar > Scalars
Definition: Utilities.h:22
Alignables & endCaps()
Return TECs.
Alignables & TIDGeomDets()
Return TID GeomDets.
Alignables & outerBarrelRods()
Return outer barrel rods.
void applyScenario(const edm::ParameterSet &scenario)
Apply misalignment scenario to the tracker.
void moveAlignableTIBTIDs(int rawId, const align::Scalars &globalDisplacements, const align::RotationType &backwardRotation, const align::RotationType &forwardRotation, bool toAndFro)
bool first
Definition: L1TdeRCT.cc:94
edm::ParameterSet theParameterSet
Alignables & endcapGeomDets()
Return endcap GeomDets.
SurveyDataConverter(const edm::ParameterSet &iConfig)
Alignables & TIDs()
Return TIDs.
const MapType & detIdMap() const
Alignables & outerBarrelGeomDets()
Return outer barrel GeomDets.
Alignables & TIDLayers()
Return TID layers.
Alignables & TIDRings()
Return TID rings.
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
static const int NFILES
AlignableTracker * getAlignableTracker()
std::vector< Alignable * > Alignables
Definition: Utilities.h:28
Alignables & innerBarrelLayers()
Return inner barrel layers.
void readFile(const std::string &textFileName, const std::string &fileType, const TrackerTopology *tTopo)
Read given text file.
tuple cout
Definition: gather_cfg.py:121
Alignables & endcapPetals()
Return encap petals.
Alignables & innerHalfBarrels()
Return TIB half barrels.
Builds a scenario from configuration and applies it to the alignable tracker.
Alignables & endcapLayers()
Return endcap layers.
virtual void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup)
Alignables & outerHalfBarrels()
Return TOB half barrels.