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