CMS 3D CMS Logo

AlignPCLThresholds.cc
Go to the documentation of this file.
5 #include <iostream>
6 #include <iomanip> // std::setw
7 
8 //****************************************************************************//
10  m_thresholds[AlignableId]=Threshold;
11 }
12 
13 //****************************************************************************//
15  m_nrecords=Nrecords;
17 }
18 
19 //****************************************************************************//
20 void AlignPCLThresholds::setNRecords(const int &Nrecords){
21  m_nrecords=Nrecords;
22 }
23 
24 //****************************************************************************//
26  threshold_map::const_iterator it = m_thresholds.find(AlignableId);
27 
28  if (it != m_thresholds.end()){
29  return it->second;
30  } else {
31  throw cms::Exception("AlignPCLThresholds")<< "No Thresholds defined for Alignable id " << AlignableId << "\n";
32  }
33 }
34 
35 //****************************************************************************//
37  return m_thresholds[AlignableId];
38 }
39 
40 float AlignPCLThresholds::getSigCut(const std::string &AlignableId,const coordType &type) const {
42  switch(type){
43  case X:
44  return a.getSigXcut();
45  case Y:
46  return a.getSigYcut();
47  case Z:
48  return a.getSigZcut();
49  case theta_X:
50  return a.getSigThetaXcut();
51  case theta_Y:
52  return a.getSigThetaYcut();
53  case theta_Z:
54  return a.getSigThetaZcut();
55  default:
56  throw cms::Exception("AlignPCLThresholds")<<"Requested significance threshold for undefined coordinate"<< type << "\n";
57  }
58 }
59 
60 //****************************************************************************//
61 // overloaded method
62 //****************************************************************************//
63 std::array<float,6> AlignPCLThresholds::getSigCut(const std::string &AlignableId) const {
65  return {{a.getSigXcut(),a.getSigYcut(), a.getSigZcut(), a.getSigThetaXcut(), a.getSigThetaYcut(),a.getSigThetaZcut()}};
66 }
67 
68 float AlignPCLThresholds::getCut(const std::string &AlignableId,const coordType &type) const {
70  switch(type){
71  case X:
72  return a.getXcut();
73  case Y:
74  return a.getYcut();
75  case Z:
76  return a.getZcut();
77  case theta_X:
78  return a.getThetaXcut();
79  case theta_Y:
80  return a.getThetaYcut();
81  case theta_Z:
82  return a.getThetaZcut();
83  default:
84  throw cms::Exception("AlignPCLThresholds")<<"Requested significance threshold for undefined coordinate"<< type << "\n";
85  }
86 }
87 
88 //****************************************************************************//
89 // overloaded method
90 //****************************************************************************//
91 std::array<float,6> AlignPCLThresholds::getCut(const std::string &AlignableId) const {
93  return {{a.getXcut(),a.getYcut(), a.getZcut(), a.getThetaXcut(), a.getThetaYcut(),a.getThetaZcut()}};
94 }
95 
96 //****************************************************************************//
97 float AlignPCLThresholds::getMaxMoveCut(const std::string &AlignableId,const coordType &type) const {
99  switch(type){
100  case X:
101  return a.getMaxMoveXcut();
102  case Y:
103  return a.getMaxMoveYcut();
104  case Z:
105  return a.getMaxMoveZcut();
106  case theta_X:
107  return a.getMaxMoveThetaXcut();
108  case theta_Y:
109  return a.getMaxMoveThetaYcut();
110  case theta_Z:
111  return a.getMaxMoveThetaZcut();
112  default:
113  throw cms::Exception("AlignPCLThresholds")<<"Requested significance threshold for undefined coordinate"<< type << "\n";
114  }
115 }
116 
117 //****************************************************************************//
118 // overloaded method
119 //****************************************************************************//
120 std::array<float,6> AlignPCLThresholds::getMaxMoveCut(const std::string &AlignableId) const {
123 }
124 
125 //****************************************************************************//
126 float AlignPCLThresholds::getMaxErrorCut(const std::string &AlignableId,const coordType &type) const {
128  switch(type){
129  case X:
130  return a.getErrorXcut();
131  case Y:
132  return a.getErrorYcut();
133  case Z:
134  return a.getErrorZcut();
135  case theta_X:
136  return a.getErrorThetaXcut();
137  case theta_Y:
138  return a.getErrorThetaYcut();
139  case theta_Z:
140  return a.getErrorThetaZcut();
141  default:
142  throw cms::Exception("AlignPCLThresholds")<<"Requested significance threshold for undefined coordinate"<< type << "\n";
143  }
144 }
145 
146 //****************************************************************************//
147 // overloaded method
148 //****************************************************************************//
149 std::array<float,6> AlignPCLThresholds::getMaxErrorCut(const std::string &AlignableId) const {
152 }
153 
154 //****************************************************************************//
155 std::array<float,4> AlignPCLThresholds::getExtraDOFCutsForAlignable(const std::string &AlignableId,const unsigned int i) const {
157  return a.getExtraDOFCuts(i);
158 }
159 
160 //****************************************************************************//
161 std::string AlignPCLThresholds::getExtraDOFLabelForAlignable(const std::string &AlignableId,const unsigned int i) const
162 {
164  return a.getExtraDOFLabel(i);
165 }
166 
167 //****************************************************************************//
169 
170  edm::LogVerbatim("AlignPCLThresholds")<<"AlignPCLThresholds::printAll()";
171  edm::LogVerbatim("AlignPCLThresholds")<<" ===================================================================================================================";
172  edm::LogVerbatim("AlignPCLThresholds")<<"N records cut: "<<this->getNrecords();
173  for(auto it = m_thresholds.begin(); it != m_thresholds.end() ; ++it){
174  edm::LogVerbatim("AlignPCLThresholds")<<" ===================================================================================================================";
175  edm::LogVerbatim("AlignPCLThresholds")<<"key : " << it->first <<" \n"
176  <<"- Xcut : " <<std::setw(4)<< (it->second).getXcut() <<std::setw(5)<<" um"
177  <<"| sigXcut : " <<std::setw(4)<< (it->second).getSigXcut() <<std::setw(1)<<" "
178  <<"| maxMoveXcut : " <<std::setw(4)<< (it->second).getMaxMoveXcut() <<std::setw(5)<<" um"
179  <<"| ErrorXcut : " <<std::setw(4)<< (it->second).getErrorXcut() <<std::setw(5)<<" um\n"
180 
181  <<"- thetaXcut : " <<std::setw(4)<< (it->second).getThetaXcut() <<std::setw(5)<<" urad"
182  <<"| sigThetaXcut : " <<std::setw(4)<< (it->second).getSigThetaXcut() <<std::setw(1)<<" "
183  <<"| maxMoveThetaXcut : " <<std::setw(4)<< (it->second).getMaxMoveThetaXcut()<<std::setw(5)<<" urad"
184  <<"| ErrorThetaXcut : " <<std::setw(4)<< (it->second).getErrorThetaXcut() <<std::setw(5)<<" urad\n"
185 
186  <<"- Ycut : " <<std::setw(4)<< (it->second).getYcut() <<std::setw(5)<<" um"
187  <<"| sigYcut : " <<std::setw(4)<< (it->second).getSigXcut() <<std::setw(1)<<" "
188  <<"| maxMoveYcut : " <<std::setw(4)<< (it->second).getMaxMoveYcut() <<std::setw(5)<<" um"
189  <<"| ErrorYcut : " <<std::setw(4)<< (it->second).getErrorYcut() <<std::setw(5)<<" um\n"
190 
191  <<"- thetaYcut : " <<std::setw(4)<< (it->second).getThetaYcut() <<std::setw(5)<<" urad"
192  <<"| sigThetaYcut : " <<std::setw(4)<< (it->second).getSigThetaYcut() <<std::setw(1)<<" "
193  <<"| maxMoveThetaYcut : " <<std::setw(4)<< (it->second).getMaxMoveThetaYcut()<<std::setw(5)<<" urad"
194  <<"| ErrorThetaYcut : " <<std::setw(4)<< (it->second).getErrorThetaYcut() <<std::setw(5)<<" urad\n"
195 
196  <<"- Zcut : " <<std::setw(4)<< (it->second).getZcut() <<std::setw(5)<<" um"
197  <<"| sigZcut : " <<std::setw(4)<< (it->second).getSigZcut() <<std::setw(1)<<" "
198  <<"| maxMoveZcut : " <<std::setw(4)<< (it->second).getMaxMoveZcut() <<std::setw(5)<<" um"
199  <<"| ErrorZcut : " <<std::setw(4)<< (it->second).getErrorZcut() <<std::setw(5)<<" um\n"
200 
201  <<"- thetaZcut : " <<std::setw(4)<< (it->second).getThetaZcut() <<std::setw(5)<<" urad"
202  <<"| sigThetaZcut : " <<std::setw(4)<< (it->second).getSigThetaZcut() <<std::setw(1)<<" "
203  <<"| maxMoveThetaZcut : " <<std::setw(4)<< (it->second).getMaxMoveThetaZcut()<<std::setw(5)<<" urad"
204  <<"| ErrorThetaZcut : " <<std::setw(4)<< (it->second).getErrorThetaZcut() <<std::setw(5)<<" urad";
205 
206  if((it->second).hasExtraDOF()){
207  for (unsigned int j=0; j<(it->second).extraDOFSize(); j++){
208  std::array<float,4> extraDOFCuts = getExtraDOFCutsForAlignable(it->first,j);
209 
210  edm::LogVerbatim("AlignPCLThresholds")<<"Extra DOF " << j << " with label: "<<getExtraDOFLabelForAlignable(it->first,j);
211  edm::LogVerbatim("AlignPCLThresholds") <<"- cut : " << std::setw(4) << extraDOFCuts.at(0) << std::setw(5)<<" "
212  <<"| sigCut : " << std::setw(4) << extraDOFCuts.at(1) << std::setw(1)<<" "
213  <<"| maxMoveCut : " << std::setw(4) << extraDOFCuts.at(2) << std::setw(5)<<" "
214  <<"| maxErrorCut : " << std::setw(4) << extraDOFCuts.at(3) << std::setw(5)<<" ";
215  }
216  }
217  }
218 }
219 
220 //****************************************************************************//
221 std::vector<std::string> AlignPCLThresholds::getAlignableList() const {
222  std::vector<std::string> alignables_;
223  alignables_.reserve(m_thresholds.size());
224 
225  for(auto it = m_thresholds.begin(); it != m_thresholds.end() ; ++it){
226  alignables_.push_back(it->first);
227  }
228  return alignables_;
229 }
float getSigCut(const std::string &AlignableId, const coordType &type) const
float getErrorThetaZcut() const
type
Definition: HCALResponse.h:21
void setAlignPCLThreshold(const std::string &AlignableId, const AlignPCLThreshold &Threshold)
float getMaxMoveThetaXcut() const
float getMaxMoveZcut() const
float getSigXcut() const
float getMaxMoveThetaZcut() const
std::array< float, 4 > getExtraDOFCuts(const unsigned int i) const
float getXcut() const
float getSigThetaXcut() const
AlignPCLThreshold getAlignPCLThreshold(const std::string &AlignableId) const
std::vector< std::string > getAlignableList() const
float getMaxMoveYcut() const
float getErrorYcut() const
std::string getExtraDOFLabelForAlignable(const std::string &AlignableId, const unsigned int i) const
float getSigThetaYcut() const
float getCut(const std::string &AlignableId, const coordType &type) const
float getZcut() const
float getThetaXcut() const
float getMaxMoveCut(const std::string &AlignableId, const coordType &type) const
float getMaxErrorCut(const std::string &AlignableId, const coordType &type) const
std::array< float, 4 > getExtraDOFCutsForAlignable(const std::string &AlignableId, const unsigned int i) const
float getThetaZcut() const
float getErrorZcut() const
void setNRecords(const int &Nrecords)
threshold_map m_thresholds
float getSigYcut() const
void setAlignPCLThresholds(const int &Nrecords, const threshold_map &Thresholds)
float getYcut() const
std::map< std::string, AlignPCLThreshold > threshold_map
const int & getNrecords() const
float getErrorThetaYcut() const
double a
Definition: hdecay.h:121
float getSigZcut() const
float getMaxMoveThetaYcut() const
float getErrorXcut() const
float getMaxMoveXcut() const
float getErrorThetaXcut() const
float getThetaYcut() const
std::string getExtraDOFLabel(const unsigned int i) const
float getSigThetaZcut() const