CMS 3D CMS Logo

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