Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:40:41

0001 #include "FWCore/Framework/interface/EventSetup.h"
0002 #include "FWCore/Framework/interface/ESHandle.h"
0003 #include "Geometry/Records/interface/PTrackerAdditionalParametersPerDetRcd.h"
0004 #include "Alignment/CommonAlignment/interface/SurveyDet.h"
0005 #include "Alignment/TrackerAlignment/interface/AlignableTracker.h"
0006 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeomBuilderFromGeometricDet.h"
0007 
0008 #include "FWCore/Framework/interface/EventSetup.h"
0009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 
0012 #include "DataFormats/DetId/interface/DetId.h"
0013 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0014 #include "Alignment/SurveyAnalysis/plugins/CreateSurveyRcds.h"
0015 #include "Geometry/CommonTopologies/interface/GeometryAligner.h"
0016 #include "CLHEP/Random/RandGauss.h"
0017 
0018 CreateSurveyRcds::CreateSurveyRcds(const edm::ParameterSet& cfg)
0019     : tTopoToken_(esConsumes()),
0020       geomDetToken_(esConsumes()),
0021       ptpToken_(esConsumes()),
0022       aliToken_(esConsumes()),
0023       aliErrToken_(esConsumes()) {
0024   m_inputGeom = cfg.getUntrackedParameter<std::string>("inputGeom");
0025   m_inputSimpleMis = cfg.getUntrackedParameter<double>("simpleMis");
0026   m_generatedRandom = cfg.getUntrackedParameter<bool>("generatedRandom");
0027   m_generatedSimple = cfg.getUntrackedParameter<bool>("generatedSimple");
0028 }
0029 
0030 void CreateSurveyRcds::analyze(const edm::Event& event, const edm::EventSetup& setup) {
0031   //Retrieve tracker topology from geometry
0032   const TrackerTopology* const tTopo = &setup.getData(tTopoToken_);
0033   const GeometricDet* geom = &setup.getData(geomDetToken_);
0034   const PTrackerParameters& ptp = setup.getData(ptpToken_);
0035   TrackerGeometry* tracker = TrackerGeomBuilderFromGeometricDet().build(geom, ptp, tTopo);
0036 
0037   //take geometry from DB or randomly generate geometry
0038   if (m_inputGeom == "sqlite") {
0039     //build the tracker
0040     const Alignments* alignments = &setup.getData(aliToken_);
0041     const AlignmentErrorsExtended* alignmentErrors = &setup.getData(aliErrToken_);
0042 
0043     //apply the latest alignments
0044     GeometryAligner aligner;
0045     aligner.applyAlignments<TrackerGeometry>(&(*tracker), alignments, alignmentErrors, AlignTransform());
0046   }
0047 
0048   addComponent(new AlignableTracker(tracker, tTopo));
0049 
0050   Alignable* ali = detector();
0051   if (m_inputGeom == "generated") {
0052     setGeometry(ali);
0053   }
0054 
0055   setSurveyErrors(ali);
0056 }
0057 
0058 void CreateSurveyRcds::setGeometry(Alignable* ali) {
0059   const align::Alignables& comp = ali->components();
0060   unsigned int nComp = comp.size();
0061   //move then do for lower level object
0062   //for issue of det vs detunit
0063   bool usecomps = true;
0064   if ((ali->alignableObjectId() == 2) && (nComp >= 1))
0065     usecomps = false;
0066   for (unsigned int i = 0; i < nComp; ++i) {
0067     if (usecomps)
0068       setGeometry(comp[i]);
0069   }
0070   DetId id(ali->id());
0071   int subdetlevel = id.subdetId();
0072   int level = ali->alignableObjectId();
0073 
0074   //for random misalignment
0075   if (m_generatedRandom) {
0076     if (subdetlevel > 0) {
0077       AlgebraicVector value = getStructureErrors(level, subdetlevel);
0078 
0079       double value0 = CLHEP::RandGauss::shoot(0, value[0]);
0080       double value1 = CLHEP::RandGauss::shoot(0, value[1]);
0081       double value2 = CLHEP::RandGauss::shoot(0, value[2]);
0082       double value3 = CLHEP::RandGauss::shoot(0, value[3]);
0083       double value4 = CLHEP::RandGauss::shoot(0, value[4]);
0084       double value5 = CLHEP::RandGauss::shoot(0, value[5]);
0085 
0086       //move/rotate the surface
0087       align::LocalVector diffR(value0, value1, value2);
0088       align::Scalar diffWx = value3;
0089       align::Scalar diffWy = value4;
0090       align::Scalar diffWz = value5;
0091       ali->move(ali->surface().toGlobal(diffR));
0092       ali->rotateAroundLocalX(diffWx);
0093       ali->rotateAroundLocalY(diffWy);
0094       ali->rotateAroundLocalZ(diffWz);
0095     }
0096   }
0097 
0098   // for simple misalignments
0099   if (m_generatedSimple) {
0100     if ((level == 2) || ((level == 1) && (ali->mother()->alignableObjectId() != 2))) {
0101       const double constMis = m_inputSimpleMis;
0102       const double dAngle = constMis / ali->surface().length();
0103       //std::cout << "Shift: " << constMis << ", Rot: " << dAngle << std::endl;
0104       double value0 = CLHEP::RandGauss::shoot(0, constMis);
0105       double value1 = CLHEP::RandGauss::shoot(0, constMis);
0106       double value2 = CLHEP::RandGauss::shoot(0, constMis);
0107       double value3 = CLHEP::RandGauss::shoot(0, dAngle);
0108       double value4 = CLHEP::RandGauss::shoot(0, dAngle);
0109       double value5 = CLHEP::RandGauss::shoot(0, dAngle);
0110 
0111       align::LocalVector diffR(value0, value1, value2);
0112       ali->move(ali->surface().toGlobal(diffR));
0113       align::Scalar diffWx = value3;
0114       align::Scalar diffWy = value4;
0115       align::Scalar diffWz = value5;
0116       ali->rotateAroundLocalX(diffWx);
0117       ali->rotateAroundLocalY(diffWy);
0118       ali->rotateAroundLocalZ(diffWz);
0119     }
0120   }
0121 }
0122 
0123 void CreateSurveyRcds::setSurveyErrors(Alignable* ali) {
0124   const align::Alignables& comp = ali->components();
0125   unsigned int nComp = comp.size();
0126   //move then do for lower level object
0127   //for issue of det vs detunit
0128   for (unsigned int i = 0; i < nComp; ++i) {
0129     setSurveyErrors(comp[i]);
0130   }
0131 
0132   DetId id(ali->id());
0133   int subdetlevel = id.subdetId();
0134   int level = ali->alignableObjectId();
0135 
0136   AlgebraicVector error = getStructureErrors(level, subdetlevel);
0137 
0138   double error0 = error[0];
0139   double error1 = error[1];
0140   double error2 = error[2];
0141   double error3 = error[3];
0142   double error4 = error[4];
0143   double error5 = error[5];
0144 
0145   // ----------------INFLATING ERRORS----------------------
0146   // inflating sensitive coordinates in each subdetector
0147   // tib
0148   if ((level <= 2) && (subdetlevel == 3)) {
0149     error0 = 0.01;
0150     error5 = 0.001;
0151     if ((level == 2) && (nComp == 2)) {
0152       error1 = 0.01;
0153     }
0154   }
0155   // tid
0156   if ((level <= 2) && (subdetlevel == 4)) {
0157     error0 = 0.01;
0158     error1 = 0.01;
0159     error2 = 0.01;
0160     error3 = 0.001;
0161     error4 = 0.001;
0162     error5 = 0.001;
0163     //error0=0.01; error1=0.002; error2=0.002; error3=0.0002; error4=0.0002; error5=0.001;
0164     //if ((level == 2)&&(nComp == 2)){
0165     //  error1 = 0.01;
0166     //}
0167   }
0168   if ((level == 23) && (subdetlevel == 4)) {  //Ring is a Disk
0169     error0 = 0.02;
0170     error1 = 0.02;
0171     error2 = 0.03;
0172     error3 = 0.0002;
0173     error4 = 0.0002;
0174     error5 = 0.0002;
0175   }
0176   if ((level == 22) && (subdetlevel == 4)) {  //Side of a Ring
0177     error0 = 0.01;
0178     error1 = 0.01;
0179     error2 = 0.01;
0180     error3 = 0.0002;
0181     error4 = 0.0002;
0182     error5 = 0.0002;
0183   }
0184   // tob
0185   if ((level <= 2) && (subdetlevel == 5)) {
0186     //error0 = 0.015; error1 = 0.015; error2 = 0.05; error3 = 0.001; error4 = 0.001; error5 = 0.001;
0187     error0 = 0.015;
0188     error1 = 0.003;
0189     error2 = 0.003;
0190     error3 = 0.0002;
0191     error4 = 0.0002;
0192     error5 = 0.001;
0193     if ((level == 2) && (nComp == 2)) {
0194       error1 = 0.015;
0195     }
0196   }
0197   if ((level == 27) && (subdetlevel == 5)) {  //Rod in a Layer
0198     error0 = 0.02;
0199     error1 = 0.02;
0200     error2 = 0.03;
0201     error3 = 0.001;
0202     error4 = 0.001;
0203     error5 = 0.001;
0204   }
0205   // tec
0206   if ((level <= 2) && (subdetlevel == 6)) {
0207     error0 = 0.02;
0208     error5 = 0.0005;
0209     if ((level == 2) && (nComp == 2)) {
0210       error1 = 0.02;
0211     }
0212   }
0213   if ((level == 34) && (subdetlevel == 6)) {  //Side on a Disk
0214     error0 = 0.01;
0215     error1 = 0.01;
0216     error2 = 0.02;
0217     error3 = 0.00005;
0218     error4 = 0.00005;
0219     error5 = 0.00005;
0220   }
0221   if ((level == 33) && (subdetlevel == 6)) {  //Petal on a Side of a Disk
0222     error0 = 0.01;
0223     error1 = 0.01;
0224     error2 = 0.02;
0225     error3 = 0.0001;
0226     error4 = 0.0001;
0227     error5 = 0.0001;
0228   }
0229   if ((level == 32) && (subdetlevel == 6)) {  //Ring on a Petal
0230     error0 = 0.007;
0231     error1 = 0.007;
0232     error2 = 0.015;
0233     error3 = 0.00015;
0234     error4 = 0.00015;
0235     error5 = 0.00015;
0236   }
0237   // ----------------INFLATING ERRORS----------------------
0238 
0239   //create the error matrix
0240   align::ErrorMatrix error_Matrix;
0241   double* errorData = error_Matrix.Array();
0242   errorData[0] = error0 * error0;
0243   errorData[2] = error1 * error1;
0244   errorData[5] = error2 * error2;
0245   errorData[9] = error3 * error3;
0246   errorData[14] = error4 * error4;
0247   errorData[20] = error5 * error5;
0248   errorData[1] = 0.0;
0249   errorData[3] = 0.0;
0250   errorData[4] = 0.0;
0251   errorData[6] = 0.0;
0252   errorData[7] = 0.0;
0253   errorData[8] = 0.0;
0254   errorData[10] = 0.0;
0255   errorData[11] = 0.0;
0256   errorData[12] = 0.0;
0257   errorData[13] = 0.0;
0258   errorData[15] = 0.0;
0259   errorData[16] = 0.0;
0260   errorData[17] = 0.0;
0261   errorData[18] = 0.0;
0262   errorData[19] = 0.0;
0263 
0264   ali->setSurvey(new SurveyDet(ali->surface(), error_Matrix));
0265 }
0266 
0267 //-------------------------------------------------------
0268 // DEFAULT VALUES FOR THE ASSEMBLY PRECISION
0269 //-------------------------------------------------------
0270 AlgebraicVector CreateSurveyRcds::getStructurePlacements(int level, int subdetlevel) {
0271   AlgebraicVector deltaRW(6);
0272   deltaRW(1) = 0.0;
0273   deltaRW(2) = 0.0;
0274   deltaRW(3) = 0.0;
0275   deltaRW(4) = 0.0;
0276   deltaRW(5) = 0.0;
0277   deltaRW(6) = 0.0;
0278   //PIXEL
0279   if ((level == 37) && (subdetlevel == 1)) {
0280     deltaRW(1) = 0.3;
0281     deltaRW(2) = 0.3;
0282     deltaRW(3) = 0.3;
0283     deltaRW(4) = 0.0017;
0284     deltaRW(5) = 0.0017;
0285     deltaRW(6) = 0.0017;
0286   }
0287   //STRIP
0288   if ((level == 38) && (subdetlevel == 3)) {
0289     deltaRW(1) = 0.3;
0290     deltaRW(2) = 0.3;
0291     deltaRW(3) = 0.3;
0292     deltaRW(4) = 0.0004;
0293     deltaRW(5) = 0.0004;
0294     deltaRW(6) = 0.0004;
0295   }
0296   //TRACKER
0297   if ((level == 39) && (subdetlevel == 1)) {
0298     deltaRW(1) = 0.0;
0299     deltaRW(2) = 0.0;
0300     deltaRW(3) = 0.0;
0301     deltaRW(4) = 0.0;
0302     deltaRW(5) = 0.0;
0303     deltaRW(6) = 0.0;
0304   }
0305   //TPB
0306   if ((level == 7) && (subdetlevel == 1)) {
0307     deltaRW(1) = 0.2;
0308     deltaRW(2) = 0.2;
0309     deltaRW(3) = 0.2;
0310     deltaRW(4) = 0.003;
0311     deltaRW(5) = 0.003;
0312     deltaRW(6) = 0.003;
0313   }
0314   if ((level == 6) && (subdetlevel == 1)) {
0315     deltaRW(1) = 0.05;
0316     deltaRW(2) = 0.05;
0317     deltaRW(3) = 0.05;
0318     deltaRW(4) = 0.0008;
0319     deltaRW(5) = 0.0008;
0320     deltaRW(6) = 0.0008;
0321   }
0322   if ((level == 5) && (subdetlevel == 1)) {
0323     deltaRW(1) = 0.02;
0324     deltaRW(2) = 0.02;
0325     deltaRW(3) = 0.02;
0326     deltaRW(4) = 0.0004;
0327     deltaRW(5) = 0.0004;
0328     deltaRW(6) = 0.0004;
0329   }
0330   if ((level == 4) && (subdetlevel == 1)) {
0331     deltaRW(1) = 0.01;
0332     deltaRW(2) = 0.01;
0333     deltaRW(3) = 0.005;
0334     deltaRW(4) = 0.0002;
0335     deltaRW(5) = 0.0002;
0336     deltaRW(6) = 0.0002;
0337   }
0338   if ((level == 2) && (subdetlevel == 1)) {
0339     deltaRW(1) = 0.005;
0340     deltaRW(2) = 0.005;
0341     deltaRW(3) = 0.003;
0342     deltaRW(4) = 0.001;
0343     deltaRW(5) = 0.001;
0344     deltaRW(6) = 0.001;
0345   }
0346   if ((level == 1) && (subdetlevel == 1)) {
0347     deltaRW(1) = 0.005;
0348     deltaRW(2) = 0.005;
0349     deltaRW(3) = 0.003;
0350     deltaRW(4) = 0.001;
0351     deltaRW(5) = 0.001;
0352     deltaRW(6) = 0.001;
0353   }
0354   //TPE
0355   if ((level == 13) && (subdetlevel == 2)) {
0356     deltaRW(1) = 0.2;
0357     deltaRW(2) = 0.2;
0358     deltaRW(3) = 0.2;
0359     deltaRW(4) = 0.0017;
0360     deltaRW(5) = 0.0017;
0361     deltaRW(6) = 0.0017;
0362   }
0363   if ((level == 12) && (subdetlevel == 2)) {
0364     deltaRW(1) = 0.05;
0365     deltaRW(2) = 0.05;
0366     deltaRW(3) = 0.05;
0367     deltaRW(4) = 0.0004;
0368     deltaRW(5) = 0.0004;
0369     deltaRW(6) = 0.0004;
0370   }
0371   if ((level == 11) && (subdetlevel == 2)) {
0372     deltaRW(1) = 0.02;
0373     deltaRW(2) = 0.02;
0374     deltaRW(3) = 0.02;
0375     deltaRW(4) = 0.001;
0376     deltaRW(5) = 0.001;
0377     deltaRW(6) = 0.001;
0378   }
0379   if ((level == 10) && (subdetlevel == 2)) {
0380     deltaRW(1) = 0.01;
0381     deltaRW(2) = 0.01;
0382     deltaRW(3) = 0.01;
0383     deltaRW(4) = 0.001;
0384     deltaRW(5) = 0.001;
0385     deltaRW(6) = 0.001;
0386   }
0387   if ((level == 9) && (subdetlevel == 2)) {
0388     deltaRW(1) = 0.01;
0389     deltaRW(2) = 0.01;
0390     deltaRW(3) = 0.005;
0391     deltaRW(4) = 0.002;
0392     deltaRW(5) = 0.002;
0393     deltaRW(6) = 0.002;
0394   }
0395   if ((level == 2) && (subdetlevel == 2)) {
0396     deltaRW(1) = 0.005;
0397     deltaRW(2) = 0.005;
0398     deltaRW(3) = 0.003;
0399     deltaRW(4) = 0.001;
0400     deltaRW(5) = 0.001;
0401     deltaRW(6) = 0.001;
0402   }
0403   if ((level == 1) && (subdetlevel == 2)) {
0404     deltaRW(1) = 0.005;
0405     deltaRW(2) = 0.005;
0406     deltaRW(3) = 0.003;
0407     deltaRW(4) = 0.001;
0408     deltaRW(5) = 0.001;
0409     deltaRW(6) = 0.001;
0410   }
0411   //TIB
0412   if ((level == 20) && (subdetlevel == 3)) {
0413     deltaRW(1) = 0.2;
0414     deltaRW(2) = 0.2;
0415     deltaRW(3) = 0.2;
0416     deltaRW(4) = 0.0017;
0417     deltaRW(5) = 0.0017;
0418     deltaRW(6) = 0.0017;
0419   }
0420   if ((level == 19) && (subdetlevel == 3)) {
0421     deltaRW(1) = 0.1;
0422     deltaRW(2) = 0.1;
0423     deltaRW(3) = 0.1;
0424     deltaRW(4) = 0.0008;
0425     deltaRW(5) = 0.0008;
0426     deltaRW(6) = 0.0008;
0427   }
0428   if ((level == 18) && (subdetlevel == 3)) {
0429     deltaRW(1) = 0.04;
0430     deltaRW(2) = 0.04;
0431     deltaRW(3) = 0.02;
0432     deltaRW(4) = 0.0006;
0433     deltaRW(5) = 0.0006;
0434     deltaRW(6) = 0.0006;
0435   }
0436   if ((level == 17) && (subdetlevel == 3)) {
0437     deltaRW(1) = 0.03;
0438     deltaRW(2) = 0.03;
0439     deltaRW(3) = 0.015;
0440     deltaRW(4) = 0.0004;
0441     deltaRW(5) = 0.0004;
0442     deltaRW(6) = 0.0004;
0443   }
0444   if ((level == 16) && (subdetlevel == 3)) {
0445     deltaRW(1) = 0.01;
0446     deltaRW(2) = 0.01;
0447     deltaRW(3) = 0.01;
0448     deltaRW(4) = 0.0004;
0449     deltaRW(5) = 0.0002;
0450     deltaRW(6) = 0.0002;
0451   }
0452   if ((level == 15) && (subdetlevel == 3)) {
0453     deltaRW(1) = 0.01;
0454     deltaRW(2) = 0.01;
0455     deltaRW(3) = 0.01;
0456     deltaRW(4) = 0.0004;
0457     deltaRW(5) = 0.0002;
0458     deltaRW(6) = 0.0002;
0459   }
0460   if ((level == 2) && (subdetlevel == 3)) {
0461     deltaRW(1) = 0.005;
0462     deltaRW(2) = 0.005;
0463     deltaRW(3) = 0.005;
0464     deltaRW(4) = 0.001;
0465     deltaRW(5) = 0.0005;
0466     deltaRW(6) = 0.0005;
0467   }
0468   if ((level == 1) && (subdetlevel == 3)) {
0469     deltaRW(1) = 0.005;
0470     deltaRW(2) = 0.005;
0471     deltaRW(3) = 0.005;
0472     deltaRW(4) = 0.001;
0473     deltaRW(5) = 0.0005;
0474     deltaRW(6) = 0.0005;
0475   }
0476   //TID
0477   if ((level == 25) && (subdetlevel == 4)) {
0478     deltaRW(1) = 0.2;
0479     deltaRW(2) = 0.2;
0480     deltaRW(3) = 0.2;
0481     deltaRW(4) = 0.0013;
0482     deltaRW(5) = 0.0013;
0483     deltaRW(6) = 0.0013;
0484   }
0485   if ((level == 24) && (subdetlevel == 4)) {
0486     deltaRW(1) = 0.05;
0487     deltaRW(2) = 0.05;
0488     deltaRW(3) = 0.05;
0489     deltaRW(4) = 0.0004;
0490     deltaRW(5) = 0.0004;
0491     deltaRW(6) = 0.0004;
0492   }
0493   if ((level == 23) && (subdetlevel == 4)) {
0494     deltaRW(1) = 0.01;
0495     deltaRW(2) = 0.01;
0496     deltaRW(3) = 0.01;
0497     deltaRW(4) = 0.0001;
0498     deltaRW(5) = 0.0001;
0499     deltaRW(6) = 0.0001;
0500   }
0501   if ((level == 22) && (subdetlevel == 4)) {
0502     deltaRW(1) = 0.005;
0503     deltaRW(2) = 0.005;
0504     deltaRW(3) = 0.005;
0505     deltaRW(4) = 0.0001;
0506     deltaRW(5) = 0.0001;
0507     deltaRW(6) = 0.0001;
0508   }
0509   if ((level == 2) && (subdetlevel == 4)) {
0510     deltaRW(1) = 0.005;
0511     deltaRW(2) = 0.005;
0512     deltaRW(3) = 0.005;
0513     deltaRW(4) = 0.0005;
0514     deltaRW(5) = 0.0005;
0515     deltaRW(6) = 0.0005;
0516   }
0517   if ((level == 1) && (subdetlevel == 4)) {
0518     deltaRW(1) = 0.005;
0519     deltaRW(2) = 0.005;
0520     deltaRW(3) = 0.005;
0521     deltaRW(4) = 0.0005;
0522     deltaRW(5) = 0.0005;
0523     deltaRW(6) = 0.0005;
0524   }
0525   //TOB
0526   if ((level == 30) && (subdetlevel == 5)) {
0527     deltaRW(1) = 0.2;
0528     deltaRW(2) = 0.2;
0529     deltaRW(3) = 0.2;
0530     deltaRW(4) = 0.0008;
0531     deltaRW(5) = 0.0008;
0532     deltaRW(6) = 0.0008;
0533   }
0534   if ((level == 29) && (subdetlevel == 5)) {
0535     deltaRW(1) = 0.014;
0536     deltaRW(2) = 0.014;
0537     deltaRW(3) = 0.05;
0538     deltaRW(4) = 0.0001;
0539     deltaRW(5) = 0.0001;
0540     deltaRW(6) = 0.0001;
0541   }
0542   if ((level == 28) && (subdetlevel == 5)) {
0543     deltaRW(1) = 0.02;
0544     deltaRW(2) = 0.02;
0545     deltaRW(3) = 0.02;
0546     deltaRW(4) = 0.0001;
0547     deltaRW(5) = 0.0001;
0548     deltaRW(6) = 0.0001;
0549   }
0550   if ((level == 27) && (subdetlevel == 5)) {
0551     deltaRW(1) = 0.01;
0552     deltaRW(2) = 0.01;
0553     deltaRW(3) = 0.02;
0554     deltaRW(4) = 0.0001;
0555     deltaRW(5) = 0.0001;
0556     deltaRW(6) = 0.0001;
0557   }
0558   if ((level == 2) && (subdetlevel == 5)) {
0559     deltaRW(1) = 0.003;
0560     deltaRW(2) = 0.003;
0561     deltaRW(3) = 0.01;
0562     deltaRW(4) = 0.0002;
0563     deltaRW(5) = 0.0002;
0564     deltaRW(6) = 0.0002;
0565   }
0566   if ((level == 1) && (subdetlevel == 5)) {
0567     deltaRW(1) = 0.003;
0568     deltaRW(2) = 0.003;
0569     deltaRW(3) = 0.01;
0570     deltaRW(4) = 0.0002;
0571     deltaRW(5) = 0.0002;
0572     deltaRW(6) = 0.0002;
0573   }
0574   //TEC
0575   if ((level == 36) && (subdetlevel == 6)) {
0576     deltaRW(1) = 0.2;
0577     deltaRW(2) = 0.2;
0578     deltaRW(3) = 0.2;
0579     deltaRW(4) = 0.0008;
0580     deltaRW(5) = 0.0008;
0581     deltaRW(6) = 0.0008;
0582   }
0583   if ((level == 35) && (subdetlevel == 6)) {
0584     deltaRW(1) = 0.05;
0585     deltaRW(2) = 0.05;
0586     deltaRW(3) = 0.05;
0587     deltaRW(4) = 0.0003;
0588     deltaRW(5) = 0.0003;
0589     deltaRW(6) = 0.0003;
0590   }
0591   if ((level == 34) && (subdetlevel == 6)) {
0592     deltaRW(1) = 0.01;
0593     deltaRW(2) = 0.01;
0594     deltaRW(3) = 0.02;
0595     deltaRW(4) = 0.00005;
0596     deltaRW(5) = 0.00005;
0597     deltaRW(6) = 0.00005;
0598   }
0599   if ((level == 33) && (subdetlevel == 6)) {
0600     deltaRW(1) = 0.01;
0601     deltaRW(2) = 0.01;
0602     deltaRW(3) = 0.02;
0603     deltaRW(4) = 0.0001;
0604     deltaRW(5) = 0.0001;
0605     deltaRW(6) = 0.0001;
0606   }
0607   if ((level == 32) && (subdetlevel == 6)) {
0608     deltaRW(1) = 0.007;
0609     deltaRW(2) = 0.007;
0610     deltaRW(3) = 0.015;
0611     deltaRW(4) = 0.00015;
0612     deltaRW(5) = 0.00015;
0613     deltaRW(6) = 0.00015;
0614   }
0615   if ((level == 2) && (subdetlevel == 6)) {
0616     deltaRW(1) = 0.002;
0617     deltaRW(2) = 0.002;
0618     deltaRW(3) = 0.005;
0619     deltaRW(4) = 0.0001;
0620     deltaRW(5) = 0.0001;
0621     deltaRW(6) = 0.0001;
0622   }
0623   if ((level == 1) && (subdetlevel == 6)) {
0624     deltaRW(1) = 0.002;
0625     deltaRW(2) = 0.002;
0626     deltaRW(3) = 0.005;
0627     deltaRW(4) = 0.0001;
0628     deltaRW(5) = 0.0001;
0629     deltaRW(6) = 0.0001;
0630   }
0631 
0632   return deltaRW;
0633 }
0634 //-------------------------------------------------------
0635 // DEFAULT VALUES FOR THE PRECISION OF THE SURVEY
0636 //-------------------------------------------------------
0637 AlgebraicVector CreateSurveyRcds::getStructureErrors(int level, int subdetlevel) {
0638   AlgebraicVector deltaRW(6);
0639   deltaRW(1) = 0.0;
0640   deltaRW(2) = 0.0;
0641   deltaRW(3) = 0.0;
0642   deltaRW(4) = 0.0;
0643   deltaRW(5) = 0.0;
0644   deltaRW(6) = 0.0;
0645   //PIXEL
0646   if ((level == 37) && (subdetlevel == 1)) {  //Pixel Detector in Tracker
0647     deltaRW(1) = 0.2;
0648     deltaRW(2) = 0.2;
0649     deltaRW(3) = 0.2;
0650     deltaRW(4) = 0.0017;
0651     deltaRW(5) = 0.0017;
0652     deltaRW(6) = 0.0017;
0653   }
0654   //STRIP
0655   if ((level == 38) && (subdetlevel == 3)) {  //Strip Tracker in Tracker
0656     deltaRW(1) = 0.2;
0657     deltaRW(2) = 0.2;
0658     deltaRW(3) = 0.2;
0659     deltaRW(4) = 0.0004;
0660     deltaRW(5) = 0.0004;
0661     deltaRW(6) = 0.0004;
0662   }
0663   //TRACKER
0664   if ((level == 39) && (subdetlevel == 1)) {  //Tracker
0665     deltaRW(1) = 0.0;
0666     deltaRW(2) = 0.0;
0667     deltaRW(3) = 0.0;
0668     deltaRW(4) = 0.0;
0669     deltaRW(5) = 0.0;
0670     deltaRW(6) = 0.0;
0671   }
0672   //TPB
0673   if ((level == 7) && (subdetlevel == 1)) {  //Barrel Pixel in Pixel
0674     deltaRW(1) = 0.05;
0675     deltaRW(2) = 0.05;
0676     deltaRW(3) = 0.1;
0677     deltaRW(4) = 0.0008;
0678     deltaRW(5) = 0.0008;
0679     deltaRW(6) = 0.0008;
0680   }
0681   if ((level == 6) && (subdetlevel == 1)) {  //HalfBarrel in Barrel Pixel
0682     deltaRW(1) = 0.015;
0683     deltaRW(2) = 0.015;
0684     deltaRW(3) = 0.03;
0685     deltaRW(4) = 0.0003;
0686     deltaRW(5) = 0.0003;
0687     deltaRW(6) = 0.0003;
0688   }
0689   if ((level == 5) && (subdetlevel == 1)) {  //HalfShell in HalfBarrel
0690     deltaRW(1) = 0.005;
0691     deltaRW(2) = 0.005;
0692     deltaRW(3) = 0.01;
0693     deltaRW(4) = 0.0001;
0694     deltaRW(5) = 0.0001;
0695     deltaRW(6) = 0.0001;
0696   }
0697   if ((level == 4) && (subdetlevel == 1)) {  //Ladder in HalfShell
0698     deltaRW(1) = 0.001;
0699     deltaRW(2) = 0.001;
0700     deltaRW(3) = 0.002;
0701     deltaRW(4) = 0.00005;
0702     deltaRW(5) = 0.00005;
0703     deltaRW(6) = 0.00005;
0704   }
0705   if ((level == 2) && (subdetlevel == 1)) {  //Det in Ladder
0706     deltaRW(1) = 0.0005;
0707     deltaRW(2) = 0.001;
0708     deltaRW(3) = 0.001;
0709     deltaRW(4) = 0.0001;
0710     deltaRW(5) = 0.0001;
0711     deltaRW(6) = 0.0001;
0712   }
0713   if ((level == 1) && (subdetlevel == 1)) {  //DetUnit in Ladder
0714     deltaRW(1) = 0.0005;
0715     deltaRW(2) = 0.001;
0716     deltaRW(3) = 0.001;
0717     deltaRW(4) = 0.0001;
0718     deltaRW(5) = 0.0001;
0719     deltaRW(6) = 0.0001;
0720   }
0721   //TPE
0722   if ((level == 13) && (subdetlevel == 2)) {  //Forward Pixel in Pixel
0723     deltaRW(1) = 0.05;
0724     deltaRW(2) = 0.05;
0725     deltaRW(3) = 0.1;
0726     deltaRW(4) = 0.0004;
0727     deltaRW(5) = 0.0004;
0728     deltaRW(6) = 0.0004;
0729   }
0730   if ((level == 12) && (subdetlevel == 2)) {  //HalfCylinder in Forward Pixel
0731     deltaRW(1) = 0.015;
0732     deltaRW(2) = 0.015;
0733     deltaRW(3) = 0.03;
0734     deltaRW(4) = 0.00015;
0735     deltaRW(5) = 0.00015;
0736     deltaRW(6) = 0.00015;
0737   }
0738   if ((level == 11) && (subdetlevel == 2)) {  //HalfDisk in HalfCylinder
0739     deltaRW(1) = 0.005;
0740     deltaRW(2) = 0.005;
0741     deltaRW(3) = 0.01;
0742     deltaRW(4) = 0.0001;
0743     deltaRW(5) = 0.0001;
0744     deltaRW(6) = 0.0001;
0745   }
0746   if ((level == 10) && (subdetlevel == 2)) {  //Blade in HalfDisk
0747     deltaRW(1) = 0.001;
0748     deltaRW(2) = 0.001;
0749     deltaRW(3) = 0.002;
0750     deltaRW(4) = 0.0001;
0751     deltaRW(5) = 0.0001;
0752     deltaRW(6) = 0.0001;
0753   }
0754   if ((level == 9) && (subdetlevel == 2)) {  //Panel in Blade
0755     deltaRW(1) = 0.001;
0756     deltaRW(2) = 0.0008;
0757     deltaRW(3) = 0.0006;
0758     deltaRW(4) = 0.0002;
0759     deltaRW(5) = 0.0002;
0760     deltaRW(6) = 0.0002;
0761   }
0762   if ((level == 2) && (subdetlevel == 2)) {  //Det in Panel
0763     deltaRW(1) = 0.0005;
0764     deltaRW(2) = 0.0004;
0765     deltaRW(3) = 0.0006;
0766     deltaRW(4) = 0.0001;
0767     deltaRW(5) = 0.0003;
0768     deltaRW(6) = 0.0001;
0769   }
0770   if ((level == 1) && (subdetlevel == 2)) {  //DetUnit in Panel
0771     deltaRW(1) = 0.0005;
0772     deltaRW(2) = 0.0004;
0773     deltaRW(3) = 0.0006;
0774     deltaRW(4) = 0.0001;
0775     deltaRW(5) = 0.0003;
0776     deltaRW(6) = 0.0001;
0777   }
0778   //TIB
0779   if ((level == 20) && (subdetlevel == 3)) {  //TIB in Strip Tracker
0780     deltaRW(1) = 0.08;
0781     deltaRW(2) = 0.08;
0782     deltaRW(3) = 0.04;
0783     deltaRW(4) = 0.0017;
0784     deltaRW(5) = 0.0017;
0785     deltaRW(6) = 0.0017;
0786   }
0787   if ((level == 19) && (subdetlevel == 3)) {  //HalfBarrel in TIB
0788     deltaRW(1) = 0.04;
0789     deltaRW(2) = 0.04;
0790     deltaRW(3) = 0.02;
0791     deltaRW(4) = 0.0003;
0792     deltaRW(5) = 0.0003;
0793     deltaRW(6) = 0.0003;
0794   }
0795   if ((level == 18) && (subdetlevel == 3)) {  //Layer in HalfBarrel
0796     deltaRW(1) = 0.02;
0797     deltaRW(2) = 0.02;
0798     deltaRW(3) = 0.01;
0799     deltaRW(4) = 0.0006;
0800     deltaRW(5) = 0.0006;
0801     deltaRW(6) = 0.0006;
0802   }
0803   if ((level == 17) && (subdetlevel == 3)) {  //HalfShell in Layer
0804     deltaRW(1) = 0.01;
0805     deltaRW(2) = 0.01;
0806     deltaRW(3) = 0.005;
0807     deltaRW(4) = 0.0002;
0808     deltaRW(5) = 0.0002;
0809     deltaRW(6) = 0.0002;
0810   }
0811   if ((level == 16) && (subdetlevel == 3)) {  //Surface in a HalfShell
0812     deltaRW(1) = 0.004;
0813     deltaRW(2) = 0.004;
0814     deltaRW(3) = 0.008;
0815     deltaRW(4) = 0.0002;
0816     deltaRW(5) = 0.0001;
0817     deltaRW(6) = 0.0001;
0818   }
0819   if ((level == 15) && (subdetlevel == 3)) {  //String in a Surface
0820     deltaRW(1) = 0.004;
0821     deltaRW(2) = 0.004;
0822     deltaRW(3) = 0.008;
0823     deltaRW(4) = 0.0002;
0824     deltaRW(5) = 0.0001;
0825     deltaRW(6) = 0.0001;
0826   }
0827   if ((level == 2) && (subdetlevel == 3)) {  //Det in a String
0828     deltaRW(1) = 0.002;
0829     deltaRW(2) = 0.002;
0830     deltaRW(3) = 0.004;
0831     deltaRW(4) = 0.0004;
0832     deltaRW(5) = 0.0002;
0833     deltaRW(6) = 0.0002;
0834   }
0835   if ((level == 1) && (subdetlevel == 3)) {  //DetUnit in a String
0836     deltaRW(1) = 0.002;
0837     deltaRW(2) = 0.002;
0838     deltaRW(3) = 0.004;
0839     deltaRW(4) = 0.0004;
0840     deltaRW(5) = 0.0002;
0841     deltaRW(6) = 0.0002;
0842   }
0843   //TID
0844   if ((level == 25) && (subdetlevel == 4)) {  //TID in Strip Tracker
0845     deltaRW(1) = 0.05;
0846     deltaRW(2) = 0.05;
0847     deltaRW(3) = 0.1;
0848     deltaRW(4) = 0.0003;
0849     deltaRW(5) = 0.0003;
0850     deltaRW(6) = 0.0003;
0851   }
0852   if ((level == 24) && (subdetlevel == 4)) {  //Disk in a TID
0853     deltaRW(1) = 0.01;
0854     deltaRW(2) = 0.01;
0855     deltaRW(3) = 0.02;
0856     deltaRW(4) = 0.0001;
0857     deltaRW(5) = 0.0001;
0858     deltaRW(6) = 0.0001;
0859   }
0860   if ((level == 23) && (subdetlevel == 4)) {  //Ring is a Disk
0861     deltaRW(1) = 0.004;
0862     deltaRW(2) = 0.004;
0863     deltaRW(3) = 0.005;
0864     deltaRW(4) = 0.00004;
0865     deltaRW(5) = 0.00004;
0866     deltaRW(6) = 0.00004;
0867   }
0868   if ((level == 22) && (subdetlevel == 4)) {  //Side of a Ring
0869     deltaRW(1) = 0.002;
0870     deltaRW(2) = 0.002;
0871     deltaRW(3) = 0.002;
0872     deltaRW(4) = 0.00004;
0873     deltaRW(5) = 0.00004;
0874     deltaRW(6) = 0.00004;
0875   }
0876   if ((level == 2) && (subdetlevel == 4)) {  //Det in a Side
0877     deltaRW(1) = 0.002;
0878     deltaRW(2) = 0.002;
0879     deltaRW(3) = 0.002;
0880     deltaRW(4) = 0.0002;
0881     deltaRW(5) = 0.0002;
0882     deltaRW(6) = 0.0002;
0883   }
0884   if ((level == 1) && (subdetlevel == 4)) {  //DetUnit is a Side
0885     deltaRW(1) = 0.002;
0886     deltaRW(2) = 0.002;
0887     deltaRW(3) = 0.002;
0888     deltaRW(4) = 0.0002;
0889     deltaRW(5) = 0.0002;
0890     deltaRW(6) = 0.0002;
0891   }
0892   //TOB
0893   if ((level == 30) && (subdetlevel == 5)) {  // TOB in Strip Tracker
0894     deltaRW(1) = 0.06;
0895     deltaRW(2) = 0.06;
0896     deltaRW(3) = 0.06;
0897     deltaRW(4) = 0.00025;
0898     deltaRW(5) = 0.00025;
0899     deltaRW(6) = 0.00025;
0900   }
0901   if ((level == 29) && (subdetlevel == 5)) {  //HalfBarrel in the TOB
0902     deltaRW(1) = 0.014;
0903     deltaRW(2) = 0.014;
0904     deltaRW(3) = 0.05;
0905     deltaRW(4) = 0.0001;
0906     deltaRW(5) = 0.0001;
0907     deltaRW(6) = 0.0001;
0908   }
0909   if ((level == 28) && (subdetlevel == 5)) {  //Layer in a HalfBarrel
0910     deltaRW(1) = 0.02;
0911     deltaRW(2) = 0.02;
0912     deltaRW(3) = 0.02;
0913     deltaRW(4) = 0.0001;
0914     deltaRW(5) = 0.0001;
0915     deltaRW(6) = 0.0001;
0916   }
0917   if ((level == 27) && (subdetlevel == 5)) {  //Rod in a Layer
0918     deltaRW(1) = 0.01;
0919     deltaRW(2) = 0.01;
0920     deltaRW(3) = 0.02;
0921     deltaRW(4) = 0.0001;
0922     deltaRW(5) = 0.0001;
0923     deltaRW(6) = 0.0001;
0924   }
0925   if ((level == 2) && (subdetlevel == 5)) {  //Det in a Rod
0926     deltaRW(1) = 0.003;
0927     deltaRW(2) = 0.003;
0928     deltaRW(3) = 0.01;
0929     deltaRW(4) = 0.0002;
0930     deltaRW(5) = 0.0002;
0931     deltaRW(6) = 0.0002;
0932   }
0933   if ((level == 1) && (subdetlevel == 5)) {  //DetUnit in a Rod
0934     deltaRW(1) = 0.003;
0935     deltaRW(2) = 0.003;
0936     deltaRW(3) = 0.01;
0937     deltaRW(4) = 0.0002;
0938     deltaRW(5) = 0.0002;
0939     deltaRW(6) = 0.0002;
0940   }
0941   //TEC
0942   if ((level == 36) && (subdetlevel == 6)) {  //TEC in the Strip Tracker
0943     deltaRW(1) = 0.06;
0944     deltaRW(2) = 0.06;
0945     deltaRW(3) = 0.1;
0946     deltaRW(4) = 0.0003;
0947     deltaRW(5) = 0.0003;
0948     deltaRW(6) = 0.0003;
0949   }
0950   if ((level == 35) && (subdetlevel == 6)) {  //Disk in the TEC
0951     deltaRW(1) = 0.015;
0952     deltaRW(2) = 0.015;
0953     deltaRW(3) = 0.03;
0954     deltaRW(4) = 0.0001;
0955     deltaRW(5) = 0.0001;
0956     deltaRW(6) = 0.0001;
0957   }
0958   if ((level == 34) && (subdetlevel == 6)) {  //Side on a Disk
0959     deltaRW(1) = 0.01;
0960     deltaRW(2) = 0.01;
0961     deltaRW(3) = 0.02;
0962     deltaRW(4) = 0.00005;
0963     deltaRW(5) = 0.00005;
0964     deltaRW(6) = 0.00005;
0965   }
0966   if ((level == 33) && (subdetlevel == 6)) {  //Petal on a Side of a Disk
0967     deltaRW(1) = 0.01;
0968     deltaRW(2) = 0.01;
0969     deltaRW(3) = 0.02;
0970     deltaRW(4) = 0.0001;
0971     deltaRW(5) = 0.0001;
0972     deltaRW(6) = 0.0001;
0973   }
0974   if ((level == 32) && (subdetlevel == 6)) {  //Ring on a Petal
0975     deltaRW(1) = 0.007;
0976     deltaRW(2) = 0.007;
0977     deltaRW(3) = 0.015;
0978     deltaRW(4) = 0.00015;
0979     deltaRW(5) = 0.00015;
0980     deltaRW(6) = 0.00015;
0981   }
0982   if ((level == 2) && (subdetlevel == 6)) {  //Det on a Ring
0983     deltaRW(1) = 0.002;
0984     deltaRW(2) = 0.002;
0985     deltaRW(3) = 0.005;
0986     deltaRW(4) = 0.0001;
0987     deltaRW(5) = 0.0001;
0988     deltaRW(6) = 0.0001;
0989   }
0990   if ((level == 1) && (subdetlevel == 6)) {  // DetUnit on a Ring
0991     deltaRW(1) = 0.002;
0992     deltaRW(2) = 0.002;
0993     deltaRW(3) = 0.005;
0994     deltaRW(4) = 0.0001;
0995     deltaRW(5) = 0.0001;
0996     deltaRW(6) = 0.0001;
0997   }
0998 
0999   return deltaRW;
1000 }
1001 
1002 // Plug in to framework
1003 
1004 #include "FWCore/Framework/interface/MakerMacros.h"
1005 
1006 DEFINE_FWK_MODULE(CreateSurveyRcds);