CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/Alignment/SurveyAnalysis/plugins/CreateSurveyRcds.cc

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