CMS 3D CMS Logo

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