Go to the documentation of this file.00001 #include "FWCore/Framework/interface/EventSetup.h"
00002 #include "FWCore/Framework/interface/ESHandle.h"
00003 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00004 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeomBuilderFromGeometricDet.h"
00005
00006 #include "Alignment/CommonAlignment/interface/SurveyDet.h"
00007 #include "Alignment/TrackerAlignment/interface/AlignableTracker.h"
00008
00009 #include "CondFormats/Alignment/interface/Alignments.h"
00010 #include "CondFormats/AlignmentRecord/interface/TrackerAlignmentRcd.h"
00011 #include "CondFormats/Alignment/interface/AlignmentErrors.h"
00012 #include "CondFormats/AlignmentRecord/interface/TrackerAlignmentErrorRcd.h"
00013
00014 #include "FWCore/Framework/interface/EventSetup.h"
00015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017
00018 #include "DataFormats/DetId/interface/DetId.h"
00019 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00020 #include "Alignment/SurveyAnalysis/plugins/CreateSurveyRcds.h"
00021 #include "Geometry/TrackingGeometryAligner/interface/GeometryAligner.h"
00022 #include "CLHEP/Random/RandGauss.h"
00023
00024
00025 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
00026 #include "FWCore/ServiceRegistry/interface/Service.h"
00027
00028 CreateSurveyRcds::CreateSurveyRcds(const edm::ParameterSet& cfg)
00029 : theParameterSet( cfg )
00030 {
00031 m_inputGeom = cfg.getUntrackedParameter< std::string > ("inputGeom");
00032 m_inputSimpleMis = cfg.getUntrackedParameter< double > ("simpleMis");
00033 m_generatedRandom = cfg.getUntrackedParameter< bool > ("generatedRandom");
00034 m_generatedSimple = cfg.getUntrackedParameter< bool > ("generatedSimple");
00035 }
00036
00037 void CreateSurveyRcds::analyze(const edm::Event& event, const edm::EventSetup& setup){
00038
00039 edm::ESHandle<GeometricDet> geom;
00040 setup.get<IdealGeometryRecord>().get(geom);
00041 TrackerGeometry* tracker = TrackerGeomBuilderFromGeometricDet().build(&*geom, theParameterSet);
00042
00043
00044 if (m_inputGeom == "sqlite"){
00045
00046 edm::ESHandle<Alignments> alignments;
00047 edm::ESHandle<AlignmentErrors> alignmentErrors;
00048
00049 setup.get<TrackerAlignmentRcd>().get(alignments);
00050 setup.get<TrackerAlignmentErrorRcd>().get(alignmentErrors);
00051
00052
00053 GeometryAligner aligner;
00054 aligner.applyAlignments<TrackerGeometry>( &(*tracker), &(*alignments), &(*alignmentErrors), AlignTransform() );
00055
00056 }
00057
00058
00059 addComponent(new AlignableTracker( tracker ) );
00060
00061 Alignable* ali = detector();
00062 if(m_inputGeom == "generated"){
00063 setGeometry( ali );
00064 }
00065
00066 setSurveyErrors( ali );
00067
00068 }
00069
00070 void CreateSurveyRcds::setGeometry(Alignable* ali)
00071 {
00072
00073 const align::Alignables& comp = ali->components();
00074 unsigned int nComp = comp.size();
00075
00076
00077 bool usecomps = true;
00078 if ((ali->alignableObjectId()==2)&&(nComp>=1)) usecomps = false;
00079 for (unsigned int i = 0; i < nComp; ++i){
00080 if (usecomps) setGeometry(comp[i]);
00081 }
00082 DetId id( ali->id() );
00083 int subdetlevel = id.subdetId();
00084 int level = ali->alignableObjectId();
00085
00086
00087 if (m_generatedRandom){
00088 if (subdetlevel > 0){
00089 AlgebraicVector value = getStructureErrors(level, subdetlevel);
00090
00091 double value0 = CLHEP::RandGauss::shoot(0,value[0]);
00092 double value1 = CLHEP::RandGauss::shoot(0,value[1]);
00093 double value2 = CLHEP::RandGauss::shoot(0,value[2]);
00094 double value3 = CLHEP::RandGauss::shoot(0,value[3]);
00095 double value4 = CLHEP::RandGauss::shoot(0,value[4]);
00096 double value5 = CLHEP::RandGauss::shoot(0,value[5]);
00097
00098
00099 align::LocalVector diffR(value0,value1,value2);
00100 align::Scalar diffWx = value3;
00101 align::Scalar diffWy = value4;
00102 align::Scalar diffWz = value5;
00103 ali->move(ali->surface().toGlobal(diffR));
00104 ali->rotateAroundLocalX(diffWx);
00105 ali->rotateAroundLocalY(diffWy);
00106 ali->rotateAroundLocalZ(diffWz);
00107 }
00108 }
00109
00110
00111 if (m_generatedSimple){
00112 if ((level == 2)||((level == 1)&&(ali->mother()->alignableObjectId() != 2))){
00113
00114 const double constMis = m_inputSimpleMis;
00115 const double dAngle = constMis/ali->surface().length();
00116
00117 double value0 = CLHEP::RandGauss::shoot(0, constMis);
00118 double value1 = CLHEP::RandGauss::shoot(0, constMis);
00119 double value2 = CLHEP::RandGauss::shoot(0, constMis);
00120 double value3 = CLHEP::RandGauss::shoot(0, dAngle);
00121 double value4 = CLHEP::RandGauss::shoot(0, dAngle);
00122 double value5 = CLHEP::RandGauss::shoot(0, dAngle);
00123
00124 align::LocalVector diffR(value0,value1,value2);
00125 ali->move(ali->surface().toGlobal(diffR));
00126 align::Scalar diffWx = value3;
00127 align::Scalar diffWy = value4;
00128 align::Scalar diffWz = value5;
00129 ali->rotateAroundLocalX(diffWx);
00130 ali->rotateAroundLocalY(diffWy);
00131 ali->rotateAroundLocalZ(diffWz);
00132 }
00133 }
00134
00135 }
00136
00137 void CreateSurveyRcds::setSurveyErrors( Alignable* ali ){
00138
00139 const align::Alignables& comp = ali->components();
00140 unsigned int nComp = comp.size();
00141
00142
00143 for (unsigned int i = 0; i < nComp; ++i){
00144 setSurveyErrors(comp[i]);
00145 }
00146
00147 DetId id( ali->id() );
00148 int subdetlevel = id.subdetId();
00149 int level = ali->alignableObjectId();
00150
00151 AlgebraicVector error = getStructureErrors(level, subdetlevel);
00152
00153 double error0 = error[0];
00154 double error1 = error[1];
00155 double error2 = error[2];
00156 double error3 = error[3];
00157 double error4 = error[4];
00158 double error5 = error[5];
00159
00160
00161
00162
00163 if ((level <= 2)&&(subdetlevel == 3)){
00164 error0 = 0.01; error5 = 0.001;
00165 if ((level==2)&&(nComp == 2)){
00166 error1 = 0.01;
00167 }
00168 }
00169
00170 if ((level <= 2)&&(subdetlevel == 4)){
00171 error0=0.01; error1=0.01; error2=0.01; error3=0.001; error4=0.001; error5=0.001;
00172
00173
00174
00175
00176 }
00177 if ((level == 23)&&(subdetlevel == 4)){
00178 error0=0.02; error1=0.02; error2=0.03; error3=0.0002; error4=0.0002; error5=0.0002;
00179 }
00180 if ((level == 22)&&(subdetlevel == 4)){
00181 error0=0.01; error1=0.01; error2=0.01; error3=0.0002; error4=0.0002; error5=0.0002;
00182 }
00183
00184 if ((level <= 2)&&(subdetlevel == 5)){
00185
00186 error0 = 0.015; error1 = 0.003; error2 = 0.003; error3 = 0.0002; error4 = 0.0002; error5 = 0.001;
00187 if ((level == 2)&&(nComp == 2)){
00188 error1 = 0.015;
00189 }
00190 }
00191 if ((level == 27)&&(subdetlevel == 5)){
00192 error0=0.02; error1=0.02; error2=0.03; error3=0.001; error4=0.001; error5=0.001;
00193 }
00194
00195 if ((level <= 2)&&(subdetlevel == 6)){
00196 error0 = 0.02; error5 = 0.0005;
00197 if ((level==2)&&(nComp == 2)){
00198 error1 = 0.02;
00199 }
00200 }
00201 if ((level == 34)&&(subdetlevel == 6)){
00202 error0=0.01; error1=0.01; error2=0.02; error3=0.00005; error4=0.00005; error5=0.00005;
00203 }
00204 if ((level == 33)&&(subdetlevel == 6)){
00205 error0=0.01; error1=0.01; error2=0.02; error3=0.0001; error4=0.0001; error5=0.0001;
00206 }
00207 if ((level == 32)&&(subdetlevel == 6)){
00208 error0=0.007; error1=0.007; error2=0.015; error3=0.00015; error4=0.00015; error5=0.00015;
00209 }
00210
00211
00212
00213 align::ErrorMatrix error_Matrix;
00214 double * errorData = error_Matrix.Array();
00215 errorData[0] = error0*error0; errorData[2] = error1*error1; errorData[5] = error2*error2;
00216 errorData[9] = error3*error3; errorData[14] = error4*error4; errorData[20] = error5*error5;
00217 errorData[1] = 0.0;
00218 errorData[3] = 0.0; errorData[4] = 0.0;
00219 errorData[6] = 0.0; errorData[7] = 0.0; errorData[8] = 0.0;
00220 errorData[10] = 0.0; errorData[11] = 0.0; errorData[12] = 0.0; errorData[13] = 0.0;
00221 errorData[15] = 0.0; errorData[16] = 0.0; errorData[17] = 0.0; errorData[18] = 0.0; errorData[19] = 0.0;
00222
00223
00224 ali->setSurvey( new SurveyDet(ali->surface(), error_Matrix ) );
00225
00226 }
00227
00228
00229
00230
00231 AlgebraicVector CreateSurveyRcds::getStructurePlacements(int level, int subdetlevel)
00232 {
00233 AlgebraicVector deltaRW(6);
00234 deltaRW(1)=0.0; deltaRW(2)=0.0; deltaRW(3)=0.0; deltaRW(4)=0.0; deltaRW(5)=0.0; deltaRW(6)=0.0;
00235
00236 if ((level == 37)&&(subdetlevel == 1)){
00237 deltaRW(1)=0.3; deltaRW(2)=0.3; deltaRW(3)=0.3; deltaRW(4)=0.0017; deltaRW(5)=0.0017; deltaRW(6)=0.0017;
00238 }
00239
00240 if ((level == 38)&&(subdetlevel == 3)){
00241 deltaRW(1)=0.3; deltaRW(2)=0.3; deltaRW(3)=0.3; deltaRW(4)=0.0004; deltaRW(5)=0.0004; deltaRW(6)=0.0004;
00242 }
00243
00244 if ((level == 39)&&(subdetlevel == 1)){
00245 deltaRW(1)=0.0; deltaRW(2)=0.0; deltaRW(3)=0.0; deltaRW(4)=0.0; deltaRW(5)=0.0; deltaRW(6)=0.0;
00246 }
00247
00248 if ((level == 7)&&(subdetlevel == 1)){
00249 deltaRW(1)=0.2; deltaRW(2)=0.2; deltaRW(3)=0.2; deltaRW(4)=0.003; deltaRW(5)=0.003; deltaRW(6)=0.003;
00250 }
00251 if ((level == 6)&&(subdetlevel == 1)){
00252 deltaRW(1)=0.05; deltaRW(2)=0.05; deltaRW(3)=0.05; deltaRW(4)=0.0008; deltaRW(5)=0.0008; deltaRW(6)=0.0008;
00253 }
00254 if ((level == 5)&&(subdetlevel == 1)){
00255 deltaRW(1)=0.02; deltaRW(2)=0.02; deltaRW(3)=0.02; deltaRW(4)=0.0004; deltaRW(5)=0.0004; deltaRW(6)=0.0004;
00256 }
00257 if ((level == 4)&&(subdetlevel == 1)){
00258 deltaRW(1)=0.01; deltaRW(2)=0.01; deltaRW(3)=0.005; deltaRW(4)=0.0002; deltaRW(5)=0.0002; deltaRW(6)=0.0002;
00259 }
00260 if ((level == 2)&&(subdetlevel == 1)){
00261 deltaRW(1)=0.005; deltaRW(2)=0.005; deltaRW(3)=0.003; deltaRW(4)=0.001; deltaRW(5)=0.001; deltaRW(6)=0.001;
00262 }
00263 if ((level == 1)&&(subdetlevel == 1)){
00264 deltaRW(1)=0.005; deltaRW(2)=0.005; deltaRW(3)=0.003; deltaRW(4)=0.001; deltaRW(5)=0.001; deltaRW(6)=0.001;
00265 }
00266
00267 if ((level == 13)&&(subdetlevel == 2)){
00268 deltaRW(1)=0.2; deltaRW(2)=0.2; deltaRW(3)=0.2; deltaRW(4)=0.0017; deltaRW(5)=0.0017; deltaRW(6)=0.0017;
00269 }
00270 if ((level == 12)&&(subdetlevel == 2)){
00271 deltaRW(1)=0.05; deltaRW(2)=0.05; deltaRW(3)=0.05; deltaRW(4)=0.0004; deltaRW(5)=0.0004; deltaRW(6)=0.0004;
00272 }
00273 if ((level == 11)&&(subdetlevel == 2)){
00274 deltaRW(1)=0.02; deltaRW(2)=0.02; deltaRW(3)=0.02; deltaRW(4)=0.001; deltaRW(5)=0.001; deltaRW(6)=0.001;
00275 }
00276 if ((level == 10)&&(subdetlevel == 2)){
00277 deltaRW(1)=0.01; deltaRW(2)=0.01; deltaRW(3)=0.01; deltaRW(4)=0.001; deltaRW(5)=0.001; deltaRW(6)=0.001;
00278 }
00279 if ((level == 9)&&(subdetlevel == 2)){
00280 deltaRW(1)=0.01; deltaRW(2)=0.01; deltaRW(3)=0.005; deltaRW(4)=0.002; deltaRW(5)=0.002; deltaRW(6)=0.002;
00281 }
00282 if ((level == 2)&&(subdetlevel == 2)){
00283 deltaRW(1)=0.005; deltaRW(2)=0.005; deltaRW(3)=0.003; deltaRW(4)=0.001; deltaRW(5)=0.001; deltaRW(6)=0.001;
00284 }
00285 if ((level == 1)&&(subdetlevel == 2)){
00286 deltaRW(1)=0.005; deltaRW(2)=0.005; deltaRW(3)=0.003; deltaRW(4)=0.001; deltaRW(5)=0.001; deltaRW(6)=0.001;
00287 }
00288
00289 if ((level == 20)&&(subdetlevel == 3)){
00290 deltaRW(1)=0.2; deltaRW(2)=0.2; deltaRW(3)=0.2; deltaRW(4)=0.0017; deltaRW(5)=0.0017; deltaRW(6)=0.0017;
00291 }
00292 if ((level == 19)&&(subdetlevel == 3)){
00293 deltaRW(1)=0.1; deltaRW(2)=0.1; deltaRW(3)=0.1; deltaRW(4)=0.0008; deltaRW(5)=0.0008; deltaRW(6)=0.0008;
00294 }
00295 if ((level == 18)&&(subdetlevel == 3)){
00296 deltaRW(1)=0.04; deltaRW(2)=0.04; deltaRW(3)=0.02; deltaRW(4)=0.0006; deltaRW(5)=0.0006; deltaRW(6)=0.0006;
00297 }
00298 if ((level == 17)&&(subdetlevel == 3)){
00299 deltaRW(1)=0.03; deltaRW(2)=0.03; deltaRW(3)=0.015; deltaRW(4)=0.0004; deltaRW(5)=0.0004; deltaRW(6)=0.0004;
00300 }
00301 if ((level == 16)&&(subdetlevel == 3)){
00302 deltaRW(1)=0.01; deltaRW(2)=0.01; deltaRW(3)=0.01; deltaRW(4)=0.0004; deltaRW(5)=0.0002; deltaRW(6)=0.0002;
00303 }
00304 if ((level == 15)&&(subdetlevel == 3)){
00305 deltaRW(1)=0.01; deltaRW(2)=0.01; deltaRW(3)=0.01; deltaRW(4)=0.0004; deltaRW(5)=0.0002; deltaRW(6)=0.0002;
00306 }
00307 if ((level == 2)&&(subdetlevel == 3)){
00308 deltaRW(1)=0.005; deltaRW(2)=0.005; deltaRW(3)=0.005; deltaRW(4)=0.001; deltaRW(5)=0.0005; deltaRW(6)=0.0005;
00309 }
00310 if ((level == 1)&&(subdetlevel == 3)){
00311 deltaRW(1)=0.005; deltaRW(2)=0.005; deltaRW(3)=0.005; deltaRW(4)=0.001; deltaRW(5)=0.0005; deltaRW(6)=0.0005;
00312 }
00313
00314 if ((level == 25)&&(subdetlevel == 4)){
00315 deltaRW(1)=0.2; deltaRW(2)=0.2; deltaRW(3)=0.2; deltaRW(4)=0.0013; deltaRW(5)=0.0013; deltaRW(6)=0.0013;
00316 }
00317 if ((level == 24)&&(subdetlevel == 4)){
00318 deltaRW(1)=0.05; deltaRW(2)=0.05; deltaRW(3)=0.05; deltaRW(4)=0.0004; deltaRW(5)=0.0004; deltaRW(6)=0.0004;
00319 }
00320 if ((level == 23)&&(subdetlevel == 4)){
00321 deltaRW(1)=0.01; deltaRW(2)=0.01; deltaRW(3)=0.01; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00322 }
00323 if ((level == 22)&&(subdetlevel == 4)){
00324 deltaRW(1)=0.005; deltaRW(2)=0.005; deltaRW(3)=0.005; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00325 }
00326 if ((level == 2)&&(subdetlevel == 4)){
00327 deltaRW(1)=0.005; deltaRW(2)=0.005; deltaRW(3)=0.005; deltaRW(4)=0.0005; deltaRW(5)=0.0005; deltaRW(6)=0.0005;
00328 }
00329 if ((level == 1)&&(subdetlevel == 4)){
00330 deltaRW(1)=0.005; deltaRW(2)=0.005; deltaRW(3)=0.005; deltaRW(4)=0.0005; deltaRW(5)=0.0005; deltaRW(6)=0.0005;
00331 }
00332
00333 if ((level == 30)&&(subdetlevel == 5)){
00334 deltaRW(1)=0.2; deltaRW(2)=0.2; deltaRW(3)=0.2; deltaRW(4)=0.0008; deltaRW(5)=0.0008; deltaRW(6)=0.0008;
00335 }
00336 if ((level == 29)&&(subdetlevel == 5)){
00337 deltaRW(1)=0.014; deltaRW(2)=0.014; deltaRW(3)=0.05; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00338 }
00339 if ((level == 28)&&(subdetlevel == 5)){
00340 deltaRW(1)=0.02; deltaRW(2)=0.02; deltaRW(3)=0.02; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00341 }
00342 if ((level == 27)&&(subdetlevel == 5)){
00343 deltaRW(1)=0.01; deltaRW(2)=0.01; deltaRW(3)=0.02; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00344 }
00345 if ((level == 2)&&(subdetlevel == 5)){
00346 deltaRW(1)=0.003; deltaRW(2)=0.003; deltaRW(3)=0.01; deltaRW(4)=0.0002; deltaRW(5)=0.0002; deltaRW(6)=0.0002;
00347 }
00348 if ((level == 1)&&(subdetlevel == 5)){
00349 deltaRW(1)=0.003; deltaRW(2)=0.003; deltaRW(3)=0.01; deltaRW(4)=0.0002; deltaRW(5)=0.0002; deltaRW(6)=0.0002;
00350 }
00351
00352 if ((level == 36)&&(subdetlevel == 6)){
00353 deltaRW(1)=0.2; deltaRW(2)=0.2; deltaRW(3)=0.2; deltaRW(4)=0.0008; deltaRW(5)=0.0008; deltaRW(6)=0.0008;
00354 }
00355 if ((level == 35)&&(subdetlevel == 6)){
00356 deltaRW(1)=0.05; deltaRW(2)=0.05; deltaRW(3)=0.05; deltaRW(4)=0.0003; deltaRW(5)=0.0003; deltaRW(6)=0.0003;
00357 }
00358 if ((level == 34)&&(subdetlevel == 6)){
00359 deltaRW(1)=0.01; deltaRW(2)=0.01; deltaRW(3)=0.02; deltaRW(4)=0.00005; deltaRW(5)=0.00005; deltaRW(6)=0.00005;
00360 }
00361 if ((level == 33)&&(subdetlevel == 6)){
00362 deltaRW(1)=0.01; deltaRW(2)=0.01; deltaRW(3)=0.02; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00363 }
00364 if ((level == 32)&&(subdetlevel == 6)){
00365 deltaRW(1)=0.007; deltaRW(2)=0.007; deltaRW(3)=0.015; deltaRW(4)=0.00015; deltaRW(5)=0.00015; deltaRW(6)=0.00015;
00366 }
00367 if ((level == 2)&&(subdetlevel == 6)){
00368 deltaRW(1)=0.002; deltaRW(2)=0.002; deltaRW(3)=0.005; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00369 }
00370 if ((level == 1)&&(subdetlevel == 6)){
00371 deltaRW(1)=0.002; deltaRW(2)=0.002; deltaRW(3)=0.005; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00372 }
00373
00374 return deltaRW;
00375 }
00376
00377
00378
00379 AlgebraicVector CreateSurveyRcds::getStructureErrors(int level, int subdetlevel)
00380 {
00381 AlgebraicVector deltaRW(6);
00382 deltaRW(1)=0.0; deltaRW(2)=0.0; deltaRW(3)=0.0; deltaRW(4)=0.0; deltaRW(5)=0.0; deltaRW(6)=0.0;
00383
00384 if ((level == 37)&&(subdetlevel == 1)){
00385 deltaRW(1)=0.2; deltaRW(2)=0.2; deltaRW(3)=0.2; deltaRW(4)=0.0017; deltaRW(5)=0.0017; deltaRW(6)=0.0017;
00386 }
00387
00388 if ((level == 38)&&(subdetlevel == 3)){
00389 deltaRW(1)=0.2; deltaRW(2)=0.2; deltaRW(3)=0.2; deltaRW(4)=0.0004; deltaRW(5)=0.0004; deltaRW(6)=0.0004;
00390 }
00391
00392 if ((level == 39)&&(subdetlevel == 1)){
00393 deltaRW(1)=0.0; deltaRW(2)=0.0; deltaRW(3)=0.0; deltaRW(4)=0.0; deltaRW(5)=0.0; deltaRW(6)=0.0;
00394 }
00395
00396 if ((level == 7)&&(subdetlevel == 1)){
00397 deltaRW(1)=0.05; deltaRW(2)=0.05; deltaRW(3)=0.1; deltaRW(4)=0.0008; deltaRW(5)=0.0008; deltaRW(6)=0.0008;
00398 }
00399 if ((level == 6)&&(subdetlevel == 1)){
00400 deltaRW(1)=0.015; deltaRW(2)=0.015; deltaRW(3)=0.03; deltaRW(4)=0.0003; deltaRW(5)=0.0003; deltaRW(6)=0.0003;
00401 }
00402 if ((level == 5)&&(subdetlevel == 1)){
00403 deltaRW(1)=0.005; deltaRW(2)=0.005; deltaRW(3)=0.01; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00404 }
00405 if ((level == 4)&&(subdetlevel == 1)){
00406 deltaRW(1)=0.001; deltaRW(2)=0.001; deltaRW(3)=0.002; deltaRW(4)=0.00005; deltaRW(5)=0.00005; deltaRW(6)=0.00005;
00407 }
00408 if ((level == 2)&&(subdetlevel == 1)){
00409 deltaRW(1)=0.0005; deltaRW(2)=0.001; deltaRW(3)=0.001; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00410 }
00411 if ((level == 1)&&(subdetlevel == 1)){
00412 deltaRW(1)=0.0005; deltaRW(2)=0.001; deltaRW(3)=0.001; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00413 }
00414
00415 if ((level == 13)&&(subdetlevel == 2)){
00416 deltaRW(1)=0.05; deltaRW(2)=0.05; deltaRW(3)=0.1; deltaRW(4)=0.0004; deltaRW(5)=0.0004; deltaRW(6)=0.0004;
00417 }
00418 if ((level == 12)&&(subdetlevel == 2)){
00419 deltaRW(1)=0.015; deltaRW(2)=0.015; deltaRW(3)=0.03; deltaRW(4)=0.00015; deltaRW(5)=0.00015; deltaRW(6)=0.00015;
00420 }
00421 if ((level == 11)&&(subdetlevel == 2)){
00422 deltaRW(1)=0.005; deltaRW(2)=0.005; deltaRW(3)=0.01; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00423 }
00424 if ((level == 10)&&(subdetlevel == 2)){
00425 deltaRW(1)=0.001; deltaRW(2)=0.001; deltaRW(3)=0.002; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00426 }
00427 if ((level == 9)&&(subdetlevel == 2)){
00428 deltaRW(1)=0.001; deltaRW(2)=0.0008; deltaRW(3)=0.0006; deltaRW(4)=0.0002; deltaRW(5)=0.0002; deltaRW(6)=0.0002;
00429 }
00430 if ((level == 2)&&(subdetlevel == 2)){
00431 deltaRW(1)=0.0005; deltaRW(2)=0.0004; deltaRW(3)=0.0006; deltaRW(4)=0.0001; deltaRW(5)=0.0003; deltaRW(6)=0.0001;
00432 }
00433 if ((level == 1)&&(subdetlevel == 2)){
00434 deltaRW(1)=0.0005; deltaRW(2)=0.0004; deltaRW(3)=0.0006; deltaRW(4)=0.0001; deltaRW(5)=0.0003; deltaRW(6)=0.0001;
00435 }
00436
00437 if ((level == 20)&&(subdetlevel == 3)){
00438 deltaRW(1)=0.08; deltaRW(2)=0.08; deltaRW(3)=0.04; deltaRW(4)=0.0017; deltaRW(5)=0.0017; deltaRW(6)=0.0017;
00439 }
00440 if ((level == 19)&&(subdetlevel == 3)){
00441 deltaRW(1)=0.04; deltaRW(2)=0.04; deltaRW(3)=0.02; deltaRW(4)=0.0003; deltaRW(5)=0.0003; deltaRW(6)=0.0003;
00442 }
00443 if ((level == 18)&&(subdetlevel == 3)){
00444 deltaRW(1)=0.02; deltaRW(2)=0.02; deltaRW(3)=0.01; deltaRW(4)=0.0006; deltaRW(5)=0.0006; deltaRW(6)=0.0006;
00445 }
00446 if ((level == 17)&&(subdetlevel == 3)){
00447 deltaRW(1)=0.01; deltaRW(2)=0.01; deltaRW(3)=0.005; deltaRW(4)=0.0002; deltaRW(5)=0.0002; deltaRW(6)=0.0002;
00448 }
00449 if ((level == 16)&&(subdetlevel == 3)){
00450 deltaRW(1)=0.004; deltaRW(2)=0.004; deltaRW(3)=0.008; deltaRW(4)=0.0002; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00451 }
00452 if ((level == 15)&&(subdetlevel == 3)){
00453 deltaRW(1)=0.004; deltaRW(2)=0.004; deltaRW(3)=0.008; deltaRW(4)=0.0002; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00454 }
00455 if ((level == 2)&&(subdetlevel == 3)){
00456 deltaRW(1)=0.002; deltaRW(2)=0.002; deltaRW(3)=0.004; deltaRW(4)=0.0004; deltaRW(5)=0.0002; deltaRW(6)=0.0002;
00457 }
00458 if ((level == 1)&&(subdetlevel == 3)){
00459 deltaRW(1)=0.002; deltaRW(2)=0.002; deltaRW(3)=0.004; deltaRW(4)=0.0004; deltaRW(5)=0.0002; deltaRW(6)=0.0002;
00460 }
00461
00462 if ((level == 25)&&(subdetlevel == 4)){
00463 deltaRW(1)=0.05; deltaRW(2)=0.05; deltaRW(3)=0.1; deltaRW(4)=0.0003; deltaRW(5)=0.0003; deltaRW(6)=0.0003;
00464 }
00465 if ((level == 24)&&(subdetlevel == 4)){
00466 deltaRW(1)=0.01; deltaRW(2)=0.01; deltaRW(3)=0.02; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00467 }
00468 if ((level == 23)&&(subdetlevel == 4)){
00469 deltaRW(1)=0.004; deltaRW(2)=0.004; deltaRW(3)=0.005; deltaRW(4)=0.00004; deltaRW(5)=0.00004; deltaRW(6)=0.00004;
00470 }
00471 if ((level == 22)&&(subdetlevel == 4)){
00472 deltaRW(1)=0.002; deltaRW(2)=0.002; deltaRW(3)=0.002; deltaRW(4)=0.00004; deltaRW(5)=0.00004; deltaRW(6)=0.00004;
00473 }
00474 if ((level == 2)&&(subdetlevel == 4)){
00475 deltaRW(1)=0.002; deltaRW(2)=0.002; deltaRW(3)=0.002; deltaRW(4)=0.0002; deltaRW(5)=0.0002; deltaRW(6)=0.0002;
00476 }
00477 if ((level == 1)&&(subdetlevel == 4)){
00478 deltaRW(1)=0.002; deltaRW(2)=0.002; deltaRW(3)=0.002; deltaRW(4)=0.0002; deltaRW(5)=0.0002; deltaRW(6)=0.0002;
00479 }
00480
00481 if ((level == 30)&&(subdetlevel == 5)){
00482 deltaRW(1)=0.06; deltaRW(2)=0.06; deltaRW(3)=0.06; deltaRW(4)=0.00025; deltaRW(5)=0.00025; deltaRW(6)=0.00025;
00483 }
00484 if ((level == 29)&&(subdetlevel == 5)){
00485 deltaRW(1)=0.014; deltaRW(2)=0.014; deltaRW(3)=0.05; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00486 }
00487 if ((level == 28)&&(subdetlevel == 5)){
00488 deltaRW(1)=0.02; deltaRW(2)=0.02; deltaRW(3)=0.02; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00489 }
00490 if ((level == 27)&&(subdetlevel == 5)){
00491 deltaRW(1)=0.01; deltaRW(2)=0.01; deltaRW(3)=0.02; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00492 }
00493 if ((level == 2)&&(subdetlevel == 5)){
00494 deltaRW(1)=0.003; deltaRW(2)=0.003; deltaRW(3)=0.01; deltaRW(4)=0.0002; deltaRW(5)=0.0002; deltaRW(6)=0.0002;
00495 }
00496 if ((level == 1)&&(subdetlevel == 5)){
00497 deltaRW(1)=0.003; deltaRW(2)=0.003; deltaRW(3)=0.01; deltaRW(4)=0.0002; deltaRW(5)=0.0002; deltaRW(6)=0.0002;
00498 }
00499
00500 if ((level == 36)&&(subdetlevel == 6)){
00501 deltaRW(1)=0.06; deltaRW(2)=0.06; deltaRW(3)=0.1; deltaRW(4)=0.0003; deltaRW(5)=0.0003; deltaRW(6)=0.0003;
00502 }
00503 if ((level == 35)&&(subdetlevel == 6)){
00504 deltaRW(1)=0.015; deltaRW(2)=0.015; deltaRW(3)=0.03; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00505 }
00506 if ((level == 34)&&(subdetlevel == 6)){
00507 deltaRW(1)=0.01; deltaRW(2)=0.01; deltaRW(3)=0.02; deltaRW(4)=0.00005; deltaRW(5)=0.00005; deltaRW(6)=0.00005;
00508 }
00509 if ((level == 33)&&(subdetlevel == 6)){
00510 deltaRW(1)=0.01; deltaRW(2)=0.01; deltaRW(3)=0.02; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00511 }
00512 if ((level == 32)&&(subdetlevel == 6)){
00513 deltaRW(1)=0.007; deltaRW(2)=0.007; deltaRW(3)=0.015; deltaRW(4)=0.00015; deltaRW(5)=0.00015; deltaRW(6)=0.00015;
00514 }
00515 if ((level == 2)&&(subdetlevel == 6)){
00516 deltaRW(1)=0.002; deltaRW(2)=0.002; deltaRW(3)=0.005; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00517 }
00518 if ((level == 1)&&(subdetlevel == 6)){
00519 deltaRW(1)=0.002; deltaRW(2)=0.002; deltaRW(3)=0.005; deltaRW(4)=0.0001; deltaRW(5)=0.0001; deltaRW(6)=0.0001;
00520 }
00521
00522 return deltaRW;
00523
00524
00525 }
00526
00527
00528
00529 #include "FWCore/Framework/interface/MakerMacros.h"
00530
00531 DEFINE_FWK_MODULE(CreateSurveyRcds);