ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/CMSSW/Alignment/CommonAlignmentAlgorithm/src/AlignmentParameterBuilder.cc
Revision: 1.8
Committed: Fri Nov 3 16:28:55 2006 UTC (18 years, 6 months ago) by flucke
Content type: text/plain
Branch: MAIN
Changes since 1.7: +125 -413 lines
Log Message:
- moving selection of Alignable into new class AlignableSelector
- changing level of messages: wrong config always gives an exception,
  info instead of warning if it is just information...
- little clean up

File Contents

# User Rev Content
1 flucke 1.5 /** \file AlignableParameterBuilder.cc
2     *
3 flucke 1.7 * $Date: 2006/11/03 11:00:55 $
4     * $Revision: 1.6 $
5 flucke 1.5 */
6    
7 fronga 1.1 #include "FWCore/MessageLogger/interface/MessageLogger.h"
8 flucke 1.8 #include "FWCore/ParameterSet/interface/ParameterSet.h"
9    
10 fronga 1.1 #include "Geometry/CommonDetAlgo/interface/AlgebraicObjects.h"
11     #include "Alignment/CommonAlignment/interface/Alignable.h"
12     #include "Alignment/CommonAlignment/interface/AlignableDet.h"
13    
14     #include "Alignment/CommonAlignmentParametrization/interface/RigidBodyAlignmentParameters.h"
15     #include "Alignment/CommonAlignmentParametrization/interface/CompositeRigidBodyAlignmentParameters.h"
16    
17 flucke 1.8 //#include "Alignment/TrackerAlignment/interface/AlignableTracker.h" not needed since only forwarded
18    
19     #include "Alignment/CommonAlignmentAlgorithm/interface/AlignableSelector.h"
20    
21 fronga 1.1 // This class's header
22    
23     #include "Alignment/CommonAlignmentAlgorithm/interface/AlignmentParameterBuilder.h"
24    
25    
26     //__________________________________________________________________________________________________
27 flucke 1.8 AlignmentParameterBuilder::AlignmentParameterBuilder(AlignableTracker* alignableTracker) :
28     theAlignables(), theAlignableTracker(alignableTracker)
29 fronga 1.1 {
30 flucke 1.8 }
31 fronga 1.1
32 flucke 1.8 //__________________________________________________________________________________________________
33     AlignmentParameterBuilder::AlignmentParameterBuilder(AlignableTracker* alignableTracker,
34     const edm::ParameterSet &pSet) :
35     theAlignables(), theAlignableTracker(alignableTracker)
36     {
37     this->addSelections(pSet);
38 fronga 1.1 }
39    
40     //__________________________________________________________________________________________________
41 flucke 1.8 unsigned int AlignmentParameterBuilder::addSelections(const edm::ParameterSet &pSet)
42     {
43    
44     const char *setName = "alignableParamSelector";
45     const std::vector<std::string> selections = pSet.getParameter<std::vector<std::string> >(setName);
46    
47     unsigned int addedSets = 0;
48     AlignableSelector selector(theAlignableTracker);
49     // loop via index instead of iterator due to possible enlargement inside loop
50     for (unsigned int iSel = 0; iSel < selections.size(); ++iSel) {
51    
52     std::vector<std::string> decompSel(this->decompose(selections[iSel], ','));
53     if (decompSel.empty()) continue; // edm::LogError or even cms::Exception??
54    
55     selector.clear();
56    
57     // // special scenarios mixing alignable and parameter selection
58     // const std::string geoSelSpecial(decompSel.size() > 1 ? "," + decompSel[1] : "");
59     // if (decompSel[0] == "ScenarioA") {
60     // selections.push_back(std::string("PixelHalfBarrelDets,111000") += geoSelSpecial);
61     // selections.push_back(std::string("BarrelDSRods,111000") += geoSelSpecial);
62     // selections.push_back(std::string("BarrelSSRods,101000") += geoSelSpecial);
63     // continue;
64     // } else if (decompSel[0] == "ScenarioB") {
65     // selections.push_back(std::string("PixelHalfBarrelLadders,111000") += geoSelSpecial);
66     // selections.push_back(std::string("BarrelDSLayers,111000") += geoSelSpecial);
67     // selections.push_back(std::string("BarrelSSLayers,101000") += geoSelSpecial);
68     // continue;
69     // } else if (decompSel[0] == "CustomStripLayers") {
70     // selections.push_back(std::string("BarrelDSLayers,111000") += geoSelSpecial);
71     // selections.push_back(std::string("BarrelSSLayers,110000") += geoSelSpecial);
72     // selections.push_back(std::string("TIDLayers,111000") += geoSelSpecial);
73     // selections.push_back(std::string("TECLayers,110000") += geoSelSpecial);
74     // continue;
75     // } else if (decompSel[0] == "CustomStripRods") {
76     // selections.push_back(std::string("BarrelDSRods,111000") += geoSelSpecial);
77     // selections.push_back(std::string("BarrelSSRods,101000") += geoSelSpecial);
78     // selections.push_back(std::string("TIDRings,111000") += geoSelSpecial);
79     // selections.push_back(std::string("TECPetals,110000") += geoSelSpecial);
80     // continue;
81     // } else if (decompSel[0] == "CSA06Selection") {
82     // selections.push_back(std::string("TOBDSRods,111111") += geoSelSpecial);
83     // selections.push_back(std::string("TOBSSRods15,100111") += geoSelSpecial);
84     // selections.push_back(std::string("TIBDSDets,111111") += geoSelSpecial);
85     // selections.push_back(std::string("TIBSSDets,100111") += geoSelSpecial);
86     // continue;
87     // }
88    
89 flucke 1.5 if (decompSel.size() < 2) {
90 flucke 1.8 throw cms::Exception("BadConfig") << "@SUB=AlignmentParameterBuilder::addSelections"
91     << selections[iSel]<<" from alignableParamSelector: "
92     << " should have at least 2 ','-separated parts";
93     } else if (decompSel.size() > 2) {
94     const edm::ParameterSet geoSel(pSet.getParameter<edm::ParameterSet>(decompSel[2].c_str()));
95     selector.addSelection(decompSel[0], geoSel);
96     } else {
97     selector.addSelection(decompSel[0]); // previous selection already cleared above
98 flucke 1.5 }
99    
100 flucke 1.8 this->add(selector.selectedAlignables(), this->decodeParamSel(decompSel[1]));
101    
102 flucke 1.5 ++addedSets;
103     }
104    
105 flucke 1.8 edm::LogInfo("Alignment") << "@SUB=AlignmentParameterBuilder::addSelections"
106     << " added " << addedSets << " sets of alignables"
107     << " from PSet " << setName;
108     return addedSets;
109 flucke 1.5 }
110    
111     //__________________________________________________________________________________________________
112     std::vector<bool> AlignmentParameterBuilder::decodeParamSel(const std::string &selString) const
113     {
114    
115 flucke 1.8 std::vector<bool> result(RigidBodyAlignmentParameters::N_PARAM, false);
116 flucke 1.5 if (selString.length() != RigidBodyAlignmentParameters::N_PARAM) {
117 flucke 1.8 throw cms::Exception("BadConfig") <<"@SUB=AlignmentParameterBuilder::decodeSelections"
118     << selString << " has wrong size != "
119     << RigidBodyAlignmentParameters::N_PARAM;
120 flucke 1.5 } else {
121     // shifts
122     if (selString.substr(0,1)=="1") result[RigidBodyAlignmentParameters::dx] = true;
123     if (selString.substr(1,1)=="1") result[RigidBodyAlignmentParameters::dy] = true;
124     if (selString.substr(2,1)=="1") result[RigidBodyAlignmentParameters::dz] = true;
125     // rotations
126     if (selString.substr(3,1)=="1") result[RigidBodyAlignmentParameters::dalpha] = true;
127     if (selString.substr(4,1)=="1") result[RigidBodyAlignmentParameters::dbeta] = true;
128     if (selString.substr(5,1)=="1") result[RigidBodyAlignmentParameters::dgamma] = true;
129     }
130    
131 flucke 1.8 return result;
132 flucke 1.5 }
133    
134    
135     //__________________________________________________________________________________________________
136     std::vector<std::string>
137     AlignmentParameterBuilder::decompose(const std::string &s, std::string::value_type delimiter) const
138     {
139    
140     std::vector<std::string> result;
141    
142     std::string::size_type previousPos = 0;
143     while (true) {
144     const std::string::size_type delimiterPos = s.find(delimiter, previousPos);
145     if (delimiterPos == std::string::npos) {
146     result.push_back(s.substr(previousPos)); // until end
147     break;
148     }
149     result.push_back(s.substr(previousPos, delimiterPos - previousPos));
150 flucke 1.8 previousPos = delimiterPos + 1; // +1: skip delimiter
151 flucke 1.5 }
152    
153     return result;
154     }
155    
156     //__________________________________________________________________________________________________
157 flucke 1.4 void AlignmentParameterBuilder::add(const std::vector<Alignable*> &alignables,
158     const std::vector<bool> &sel)
159 fronga 1.1 {
160    
161     int num_adu = 0;
162     int num_det = 0;
163     int num_hlo = 0;
164    
165     // loop on Alignable objects
166     for ( std::vector<Alignable*>::const_iterator ia=alignables.begin();
167     ia!=alignables.end(); ia++ ) {
168     Alignable* ali=(*ia);
169    
170 flucke 1.8 // // select on single/double sided barrel layers
171     // std::pair<int,int> tl=theTrackerAlignableId->typeAndLayerFromAlignable( ali );
172     // int type = tl.first;
173     // int layer = tl.second;
174     //
175     // bool keep=true;
176     // if (theOnlySS) // only single sided
177     // if ( (abs(type)==3 || abs(type)==5) && layer<=2 )
178     // keep=false;
179     //
180     // if (theOnlyDS) // only double sided
181     // if ( (abs(type)==3 || abs(type)==5) && layer>2 )
182     // keep=false;
183     //
184     // // reject layers
185     // if ( theSelLayers && (layer<theMinLayer || layer>theMaxLayer) )
186     // keep=false;
187     //
188     //
189     // if (keep) {
190     AlgebraicVector par(6,0);
191     AlgebraicSymMatrix cov(6,0);
192    
193     AlignableDet* alidet = dynamic_cast<AlignableDet*>(ali);
194     if (alidet !=0) { // alignable Det
195     RigidBodyAlignmentParameters* dap =
196     new RigidBodyAlignmentParameters(ali,par,cov,sel);
197     ali->setAlignmentParameters(dap);
198     num_det++;
199     } else { // higher level object
200     CompositeRigidBodyAlignmentParameters* dap =
201     new CompositeRigidBodyAlignmentParameters(ali,par,cov,sel);
202     ali->setAlignmentParameters(dap);
203     num_hlo++;
204 fronga 1.1 }
205 flucke 1.8
206     theAlignables.push_back(ali);
207     num_adu++;
208     // }
209 fronga 1.1 }
210    
211 flucke 1.8 edm::LogInfo("Alignment") << "@SUB=AlignmentParameterBuilder::add"
212     << "Added " << num_adu
213     << " Alignables, of which " << num_det << " are Dets and "
214     << num_hlo << " are higher level.";
215 fronga 1.1 }
216    
217    
218     //__________________________________________________________________________________________________
219     void AlignmentParameterBuilder::fixAlignables(int n)
220     {
221    
222 fpschill 1.3 if (n<1 || n>3) {
223 flucke 1.5 edm::LogError("BadArgument") << " n = " << n << " is not in [1,3]";
224 fronga 1.1 return;
225     }
226    
227     std::vector<Alignable*> theNewAlignables;
228     int i=0;
229     int imax = theAlignables.size();
230     for ( std::vector<Alignable*>::const_iterator ia=theAlignables.begin();
231     ia!=theAlignables.end(); ia++ )
232     {
233     i++;
234     if ( n==1 && i>1 )
235     theNewAlignables.push_back(*ia);
236     else if ( n==2 && i>1 && i<imax )
237     theNewAlignables.push_back(*ia);
238     else if ( n==3 && i>2 && i<imax)
239     theNewAlignables.push_back(*ia);
240     }
241    
242     theAlignables = theNewAlignables;
243    
244 fpschill 1.3 edm::LogWarning("Alignment") << "removing " << n
245 flucke 1.5 << " alignables, so that " << theAlignables.size()
246     << " alignables left";
247 fronga 1.1
248     }
249