CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/Alignment/MuonAlignmentAlgorithms/plugins/MuonDTLocalMillepedeAlgorithm.cc

Go to the documentation of this file.
00001 #include <fstream>
00002 
00003 #include "TFile.h"
00004 #include "TTree.h"
00005 #include "TKey.h"
00006 #include "TMatrixD.h"
00007 #include "TF1.h"
00008 
00009 
00010 #include "Alignment/MuonAlignmentAlgorithms/src/DTMuonSLToSL.cc"
00011 #include "Alignment/MuonAlignmentAlgorithms/src/DTMuonMillepede.cc"
00012 
00013 
00014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00016 #include "FWCore/Framework/interface/MakerMacros.h"
00017 #include "FWCore/Framework/interface/Frameworkfwd.h"
00018 
00019 
00020 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00021 #include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
00022 
00023 #include "Alignment/CommonAlignment/interface/Alignable.h"  
00024 #include "Alignment/CommonAlignment/interface/AlignableDet.h"
00025 #include "Alignment/CommonAlignment/interface/AlignableNavigator.h"  
00026 #include "Alignment/CommonAlignment/interface/AlignableObjectId.h"  
00027 #include "Alignment/CommonAlignment/interface/AlignmentParameters.h"
00028 #include "Alignment/CommonAlignment/interface/SurveyResidual.h"
00029 #include "Alignment/CommonAlignmentAlgorithm/interface/AlignmentParameterStore.h"
00030 #include "Alignment/CommonAlignmentAlgorithm/interface/AlignmentParameterSelector.h"
00031 #include "Alignment/MuonAlignment/interface/AlignableMuon.h"
00032 #include <DataFormats/GeometrySurface/interface/LocalError.h> 
00033 #include "Alignment/TrackerAlignment/interface/AlignableTracker.h"
00034 #include "DataFormats/TrackReco/interface/Track.h"
00035 
00036 #include "Alignment/MuonAlignmentAlgorithms/plugins/MuonDTLocalMillepedeAlgorithm.h"
00037 #include "Alignment/CommonAlignmentAlgorithm/interface/AlignmentAlgorithmPluginFactory.h"
00038 
00039 // Constructor ----------------------------------------------------------------
00040 
00041 MuonDTLocalMillepedeAlgorithm::MuonDTLocalMillepedeAlgorithm(const edm::ParameterSet& cfg):
00042   AlignmentAlgorithmBase( cfg )
00043 {
00044   
00045   edm::LogInfo("Alignment") << "[MuonDTLocalMillepedeAlgorithm] constructed.";
00046 
00047   //Parse parameters. In the future this section should be completed with more options
00048   ntuplePath = cfg.getParameter<std::string>( "ntuplePath" );
00049   numberOfRootFiles = cfg.getParameter<int>( "numberOfRootFiles" ); 
00050   ptMax = cfg.getParameter<double>( "ptMax" );
00051   ptMin = cfg.getParameter<double>( "ptMin" );
00052   numberOfSigmasX = cfg.getParameter<double>( "numberOfSigmasX" );
00053   numberOfSigmasDXDZ = cfg.getParameter<double>( "numberOfSigmasDXDZ" );
00054   numberOfSigmasY = cfg.getParameter<double>( "numberOfSigmasY" );
00055   numberOfSigmasDYDZ = cfg.getParameter<double>( "numberOfSigmasDYDZ" );
00056   nPhihits = cfg.getParameter<double>( "nPhihits" );
00057   nThetahits = cfg.getParameter<double>( "nThetahits" );
00058   workingmode = cfg.getParameter<int>( "workingMode" );
00059   nMtxSection = cfg.getParameter<int>( "nMtxSection" );
00060 
00061  
00062   //The algorithm has three working modes:
00063   //0.- aligment information is extracted from the events and stored in root files
00064   //1.- The SLtoSl algorithm
00065   //2.- The local MuonMillepede algorithm
00066   if(workingmode == 0) {
00067     edm::LogInfo("Alignment") << "[MuonDTLocalMillepedeAlgorithm] Running on production mode."; 
00068     char nameOfFile[200];
00069     snprintf(nameOfFile, sizeof(nameOfFile), "%s/MyNtupleResidual.root", ntuplePath.c_str());
00070     f = new TFile(nameOfFile, "RECREATE");
00071     f->cd();
00072     setBranchTrees();
00073   } else if (workingmode == 1) {
00074     edm::LogInfo("Alignment") << "[MuonDTLocalMillepedeAlgorithm] Running SLToSL algorithm."; 
00075   } else {
00076     edm::LogInfo("Alignment") << "[MuonDTLocalMillepedeAlgorithm] Running Local Millepede algorithm."; 
00077   }
00078 
00079 }
00080 
00081 
00082 
00083 
00084 // Call at beginning of job ---------------------------------------------------
00085 void 
00086 MuonDTLocalMillepedeAlgorithm::initialize( const edm::EventSetup& setup, 
00087                                     AlignableTracker* tracker, AlignableMuon* muon, 
00088                                     AlignmentParameterStore* store )
00089 {
00090   
00091   // accessor Det->AlignableDet
00092   if ( !muon )
00093     theAlignableDetAccessor = new AlignableNavigator(tracker);
00094   else if ( !tracker )
00095     theAlignableDetAccessor = new AlignableNavigator(muon);
00096   else 
00097     theAlignableDetAccessor = new AlignableNavigator(tracker,muon);
00098   
00099   // set alignmentParameterStore
00100   theAlignmentParameterStore=store;
00101   
00102   // get alignables
00103   theAlignables = theAlignmentParameterStore->alignables();
00104   
00105 }
00106 
00107 
00108 
00109 // Call at end of job ---------------------------------------------------------
00110 void MuonDTLocalMillepedeAlgorithm::terminate(void)
00111 {
00112 
00113   //If workingmode equals 1 or 2, the algorithms are run before saving.   
00114   if(workingmode == 1) {
00115     edm::LogInfo("Alignment") << "[MuonDTLocalMillepedeAlgorithm] Starting SLToSL algorithm";
00116     DTMuonSLToSL mySLToSL(ntuplePath, numberOfRootFiles, ptMax, ptMin, f);
00117   } else if(workingmode >= 2) {
00118     edm::LogInfo("Alignment") << "[MuonDTLocalMillepedeAlgorithm] Starting local MuonMillepede algorithm";
00119     DTMuonMillepede myMillepede(ntuplePath, numberOfRootFiles, ptMax, ptMin, nPhihits, nThetahits, workingmode, nMtxSection);
00120   } 
00121 
00122   if (workingmode==0) {
00123     f->Write();
00124     f->Close();
00125   }
00126 
00127 }
00128 
00129 
00130 
00131 // Run the algorithm on trajectories and tracks -------------------------------
00132 void MuonDTLocalMillepedeAlgorithm::run(const edm::EventSetup& setup, const EventInfo &eventInfo)
00133 //void MuonDTLocalMillepedeAlgorithm::run( const edm::EventSetup& setup, const ConstTrajTrackPairCollection& tracks)
00134 {
00135 
00136   //Only important in the production mode
00137   if(workingmode != 0) return;
00138  
00139   const ConstTrajTrackPairCollection &tracks = eventInfo.trajTrackPairs_;
00140   for( ConstTrajTrackPairCollection::const_iterator it=tracks.begin();
00141        it!=tracks.end();it++) {
00142 
00143     const Trajectory *traj = (*it).first;
00144     const reco::Track *track = (*it).second;
00145     
00146     p    = track->p();
00147     pt    = track->pt();
00148     eta   = track->eta();
00149     phi   = track->phi();
00150     charge = track->charge();
00151  
00152 
00153     if(pt < ptMin || pt > ptMax) continue;
00154     
00155     vector<const TransientTrackingRecHit*> hitvec;
00156     vector<TrajectoryMeasurement> measurements = traj->measurements();
00157 
00158     //In this loop the measurements and hits are extracted and put into two vectors 
00159     int ch_muons = 0;
00160     for (vector<TrajectoryMeasurement>::iterator im=measurements.begin();
00161          im!=measurements.end(); im++) {
00162       TrajectoryMeasurement meas = *im;
00163       const TransientTrackingRecHit* hit = &(*meas.recHit());
00164       //We are not very strict at this point
00165       if (hit->isValid()) {
00166         if(hit->det()->geographicalId().det() == 2 && hit->det()->geographicalId().subdetId() == 1) {
00167           hitvec.push_back(hit);
00168           ch_muons++;
00169         } 
00170       }
00171     }
00172 
00173 
00174     vector<const TransientTrackingRecHit*>::const_iterator ihit=hitvec.begin();
00175     //Information is stored temporally in the myTrack1D object, which is analyzed
00176     //in the build4DSegments method in order to associate hits to segments. 
00177     int ch_counter = 0;
00178     while (ihit != hitvec.end()) 
00179       {
00180         const GeomDet* det=(*ihit)->det();
00181         if(det->geographicalId().det() == 2) {
00182           if(det->geographicalId().subdetId() == 1) {
00183             DTLayerId mLayer(det->geographicalId().rawId());
00184             DTChamberId mChamber(mLayer.wheel(), mLayer.station(), mLayer.sector());
00185             AlignableDet *aliDet = theAlignableDetAccessor->alignableDetFromDetId(mChamber);
00186             myTrack1D.wh[ch_counter] = mLayer.wheel();
00187             myTrack1D.st[ch_counter] = mLayer.station();
00188             myTrack1D.sr[ch_counter] = mLayer.sector();
00189             myTrack1D.sl[ch_counter] = mLayer.superlayer();
00190             myTrack1D.la[ch_counter] = mLayer.layer();
00191             myTrack1D.erx[ch_counter] = (*ihit)->localPositionError().xx();
00192             align::GlobalPoint globhit = det->toGlobal((*ihit)->localPosition());
00193             align::LocalPoint seghit = aliDet->surface().toLocal(globhit);
00194             myTrack1D.xc[ch_counter] = seghit.x();
00195             myTrack1D.yc[ch_counter] = seghit.y();
00196             myTrack1D.zc[ch_counter] = seghit.z();
00197             ch_counter++;
00198           }
00199         }
00200         ihit++;
00201       }
00202     myTrack1D.nhits = ch_counter;
00203     if(build4DSegments()) ttreeOutput->Fill();
00204   }
00205   
00206 }
00207 
00208 
00209 
00210 //This method separates the hits depending on the chamber and builds the 4D segments
00211 //through linear fits.
00212 //It returns true if 4Dsegments were correctly built and false otherwise.  
00213 //--------------------------------------------------------------------------------------------------------------
00214 bool MuonDTLocalMillepedeAlgorithm::build4DSegments() {
00215 
00216   bool saveThis = false;
00217   
00218   //Set to 0
00219   int id[20][5];
00220   int numlayer[20][12];
00221   for(int s = 0; s < 20; ++s) {
00222     for(int k = 0; k < 5; ++k) id[s][k] = 0;
00223     for(int k = 0; k < 12; ++k) numlayer[s][k] = 0;
00224   }
00225 
00226   
00227   int nChambers = 0;
00228   for(int counter = 0; counter < myTrack1D.nhits; ++counter) {
00229     bool isNew = true;
00230     for(int counterCham = 0; counterCham < nChambers; counterCham++) {
00231       if(myTrack1D.wh[counter] == id[counterCham][0] &&
00232          myTrack1D.st[counter] == id[counterCham][1] &&
00233          myTrack1D.sr[counter] == id[counterCham][2]) {
00234         if(myTrack1D.sl[counter] == 2) { 
00235           id[counterCham][4]++;
00236         } else {
00237           id[counterCham][3]++;
00238         }
00239         for (int ila = 1; ila<=4; ila++)
00240           if (myTrack1D.la[counter]==ila) {
00241             int jla = (myTrack1D.sl[counter]-1)*4 + ila -1;
00242             numlayer[counterCham][jla]++;
00243           }
00244         isNew = false;
00245       } 
00246     }
00247     if(isNew) {
00248       id[nChambers][0] = myTrack1D.wh[counter];
00249       id[nChambers][1] = myTrack1D.st[counter];
00250       id[nChambers][2] = myTrack1D.sr[counter];
00251       if(myTrack1D.sl[counter] == 2) {
00252         id[nChambers][4]++;
00253       } else {
00254         id[nChambers][3]++;
00255       }
00256       for (int ila = 1; ila<=4; ila++)
00257         if (myTrack1D.la[counter]==ila) {
00258           int jla = (myTrack1D.sl[counter]-1)*4 + ila -1;
00259           numlayer[nChambers][jla]++;
00260         }
00261       nChambers++;
00262     }
00263   }
00264   
00265   for (int iseg = 0; iseg<MAX_SEGMENT; iseg++) 
00266     for (int ihit = 0; ihit<MAX_HIT_CHAM; ihit++) {
00267       xc[iseg][ihit] = -250.;
00268       yc[iseg][ihit] = -250.;
00269       zc[iseg][ihit] = -250.;
00270       ex[iseg][ihit] = -250.;
00271       xcp[iseg][ihit] = -250.;
00272       ycp[iseg][ihit] = -250.;
00273       excp[iseg][ihit] = -250.;
00274       eycp[iseg][ihit] = -250.;
00275       sl[iseg][ihit] = 0;
00276       la[iseg][ihit] = 0;
00277     }
00278 
00279   nseg = 0;
00280   for(int counter = 0; counter < nChambers; ++counter) {
00281     
00282     bool GoodPhiChamber = true, GoodThetaChamber = true;
00283     for (int ila = 1; ila<=12; ila++) {
00284       if (numlayer[counter][ila-1]!=1 && (ila<5 || ila>8)) GoodPhiChamber = false;
00285       if (numlayer[counter][ila-1]!=1 && (ila<9 || ila>4) && id[counter][1]!=4) GoodThetaChamber = false;
00286     }
00287 
00288     if(id[counter][3] >= nPhihits && (id[counter][4] >= nThetahits || id[counter][1] == 4) &&
00289        GoodPhiChamber && GoodThetaChamber) {
00290 
00291       TMatrixD phiProjection(2,2);
00292       TMatrixD thetaProjection(2,2);
00293       TMatrixD bphiProjection(2,1);
00294       TMatrixD bthetaProjection(2,1);
00295       
00296       TMatrixD phiProjectionSL1(2,2);
00297       TMatrixD bphiProjectionSL1(2,1);
00298       TMatrixD phiProjectionSL3(2,2);
00299       TMatrixD bphiProjectionSL3(2,1);
00300      
00301       float SL1_z_ave = 0;
00302       float SL3_z_ave = 0;
00303 
00304       int numh1 = 0, numh2 = 0, numh3 = 0;
00305       for(int counterH = 0; counterH < myTrack1D.nhits; ++counterH) {
00306         if(myTrack1D.wh[counterH] == id[counter][0] && myTrack1D.st[counterH] == id[counter][1] &&
00307            myTrack1D.sr[counterH] == id[counter][2]) {
00308           if(myTrack1D.sl[counterH] == 2) {
00309             numh2++;
00310             thetaProjection(0,0) += 1.0/myTrack1D.erx[counterH];
00311             thetaProjection(0,1) += myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00312             thetaProjection(1,0) += myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00313             thetaProjection(1,1) += myTrack1D.zc[counterH]*myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00314             bthetaProjection(0,0) += myTrack1D.yc[counterH]/myTrack1D.erx[counterH];
00315             bthetaProjection(1,0) += myTrack1D.yc[counterH]*myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00316           } else {
00317             phiProjection(0,0) += 1.0/myTrack1D.erx[counterH];
00318             phiProjection(0,1) += myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00319             phiProjection(1,0) += myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00320             phiProjection(1,1) += myTrack1D.zc[counterH]*myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00321             bphiProjection(0,0) += myTrack1D.xc[counterH]/myTrack1D.erx[counterH];
00322             bphiProjection(1,0) += myTrack1D.xc[counterH]*myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00323             if(myTrack1D.sl[counterH] == 1) {
00324               numh1++;
00325               phiProjectionSL1(0,0) += 1.0/myTrack1D.erx[counterH];
00326               phiProjectionSL1(0,1) += myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00327               phiProjectionSL1(1,0) += myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00328               phiProjectionSL1(1,1) += myTrack1D.zc[counterH]*myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00329               bphiProjectionSL1(0,0) += myTrack1D.xc[counterH]/myTrack1D.erx[counterH];
00330               bphiProjectionSL1(1,0) += myTrack1D.xc[counterH]*myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00331               SL1_z_ave += myTrack1D.zc[counterH];
00332             } else {
00333               numh3++;
00334               phiProjectionSL3(0,0) += 1.0/myTrack1D.erx[counterH];
00335               phiProjectionSL3(0,1) += myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00336               phiProjectionSL3(1,0) += myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00337               phiProjectionSL3(1,1) += myTrack1D.zc[counterH]*myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00338               bphiProjectionSL3(0,0) += myTrack1D.xc[counterH]/myTrack1D.erx[counterH];
00339               bphiProjectionSL3(1,0) += myTrack1D.xc[counterH]*myTrack1D.zc[counterH]/myTrack1D.erx[counterH];
00340               SL3_z_ave += myTrack1D.zc[counterH];
00341             } 
00342           }
00343         }
00344       }
00345 
00346       SL1_z_ave /= 4.0;
00347       SL3_z_ave /= 4.0;
00348       
00349       if (phiProjection(0,0) != 0 && phiProjectionSL1(0,0) != 0 && phiProjectionSL3(0,0) != 0 && 
00350           (thetaProjection(0,0) != 0 || id[counter][1] == 4)) {
00351         
00352         wh[nseg] = id[counter][0];
00353         st[nseg] = id[counter][1];
00354         sr[nseg] = id[counter][2];
00355         
00356         if(thetaProjection(0,0) != 0 && id[counter][1] != 4) { // Already asked (almost)
00357           thetaProjection.Invert();
00358           TMatrixD solution = thetaProjection*bthetaProjection;
00359           ySl[nseg] = solution(0,0);
00360           dydzSl[nseg] = solution(1,0);
00361           eySl[nseg] = thetaProjection(0,0);
00362           edydzSl[nseg] = thetaProjection(1,1);
00363           eydydzSl[nseg] = thetaProjection(0,1);
00364         }
00365         phiProjection.Invert();
00366         phiProjectionSL1.Invert();
00367         phiProjectionSL3.Invert();
00368         TMatrixD solution = phiProjection*bphiProjection;
00369         TMatrixD solutionSL1 = phiProjectionSL1*bphiProjectionSL1;
00370         TMatrixD solutionSL3 = phiProjectionSL3*bphiProjectionSL3;
00371         xSl[nseg] = solution(0,0);
00372         dxdzSl[nseg] = solution(1,0);
00373         exSl[nseg] = phiProjection(0,0);
00374         edxdzSl[nseg] = phiProjection(1,1);
00375         exdxdzSl[nseg] = phiProjection(0,0);
00376         xSlSL1[nseg] = solutionSL1(0,0);
00377         dxdzSlSL1[nseg] = solutionSL1(1,0);
00378         exSlSL1[nseg] = phiProjectionSL1(0,0);
00379         edxdzSlSL1[nseg] = phiProjectionSL1(1,1);
00380         exdxdzSlSL1[nseg] = phiProjectionSL1(0,0);
00381         xSL1SL3[nseg] = solutionSL1(0,0) + SL3_z_ave * solutionSL1(1,0);
00382         xSlSL3[nseg] = solutionSL3(0,0);
00383         dxdzSlSL3[nseg] = solutionSL3(1,0);
00384         exSlSL3[nseg] = phiProjectionSL3(0,0);
00385         edxdzSlSL3[nseg] = phiProjectionSL3(1,1);
00386         exdxdzSlSL3[nseg] = phiProjectionSL3(0,0);
00387         xSL3SL1[nseg] = solutionSL3(0,0) + SL1_z_ave * solutionSL3(1,0);
00388         int hitcounter = 0;
00389         for(int counterH = 0; counterH < myTrack1D.nhits; ++counterH) {
00390           if(myTrack1D.wh[counterH] == wh[nseg] && myTrack1D.st[counterH] == st[nseg] &&
00391              myTrack1D.sr[counterH] == sr[nseg]) {
00392             xc[nseg][hitcounter] = myTrack1D.xc[counterH];
00393             yc[nseg][hitcounter] = myTrack1D.yc[counterH];
00394             zc[nseg][hitcounter] = myTrack1D.zc[counterH];
00395             ex[nseg][hitcounter] = myTrack1D.erx[counterH];
00396             xcp[nseg][hitcounter] = xSl[nseg]+dxdzSl[nseg]*myTrack1D.zc[counterH];
00397             ycp[nseg][hitcounter] = ySl[nseg]+dydzSl[nseg]*myTrack1D.zc[counterH];
00398             excp[nseg][hitcounter] = exSl[nseg]*exSl[nseg]+ (edxdzSl[nseg]*edxdzSl[nseg])*myTrack1D.zc[counterH];
00399             eycp[nseg][hitcounter] = eySl[nseg]*eySl[nseg]+ (edydzSl[nseg]*edydzSl[nseg])*myTrack1D.zc[counterH];
00400             sl[nseg][hitcounter] = myTrack1D.sl[counterH];
00401             la[nseg][hitcounter] = myTrack1D.la[counterH];
00402             saveThis = true;
00403             hitcounter++;
00404           }
00405         }
00406         nphihits[nseg] = id[counter][3];
00407         nthetahits[nseg] = id[counter][4];
00408         nhits[nseg] = hitcounter;
00409         nseg++;
00410       }
00411     }
00412   }
00413   return saveThis;
00414 }
00415 
00416 
00417 
00418 //The tree structure is defined and the variables associated ---------------------------------
00419 void MuonDTLocalMillepedeAlgorithm::setBranchTrees() {
00420 
00421   ttreeOutput = new TTree("InfoTuple", "InfoTuple");
00422   
00423   ttreeOutput->Branch("p", &p, "p/F");
00424   ttreeOutput->Branch("pt", &pt, "pt/F");
00425   ttreeOutput->Branch("eta", &eta, "eta/F");
00426   ttreeOutput->Branch("phi", &phi, "phi/F");
00427   ttreeOutput->Branch("charge", &charge, "charge/F");
00428   ttreeOutput->Branch("nseg", &nseg, "nseg/I");
00429   ttreeOutput->Branch("nphihits", nphihits, "nphihits[nseg]/I");
00430   ttreeOutput->Branch("nthetahits", nthetahits, "nthetahits[nseg]/I");
00431   ttreeOutput->Branch("nhits", nhits, "nhits[nseg]/I");
00432   ttreeOutput->Branch("xSl", xSl, "xSl[nseg]/F");
00433   ttreeOutput->Branch("dxdzSl", dxdzSl, "dxdzSl[nseg]/F");
00434   ttreeOutput->Branch("exSl", exSl, "exSl[nseg]/F");
00435   ttreeOutput->Branch("edxdzSl", edxdzSl, "edxdzSl[nseg]/F");
00436   ttreeOutput->Branch("exdxdzSl", edxdzSl, "exdxdzSl[nseg]/F");
00437   ttreeOutput->Branch("ySl", ySl, "ySl[nseg]/F");
00438   ttreeOutput->Branch("dydzSl", dydzSl, "dydzSl[nseg]/F");
00439   ttreeOutput->Branch("eySl", eySl, "eySl[nseg]/F");
00440   ttreeOutput->Branch("edydzSl", edydzSl, "edydzSl[nseg]/F");
00441   ttreeOutput->Branch("eydydzSl", eydydzSl, "eydydzSl[nseg]/F");
00442   ttreeOutput->Branch("xSlSL1", xSlSL1, "xSlSL1[nseg]/F");
00443   ttreeOutput->Branch("dxdzSlSL1", dxdzSlSL1, "dxdzSlSL1[nseg]/F");
00444   ttreeOutput->Branch("exSlSL1", exSlSL1, "exSlSL1[nseg]/F");
00445   ttreeOutput->Branch("edxdzSlSL1", edxdzSlSL1, "edxdzSlSL1[nseg]/F");
00446   ttreeOutput->Branch("xSL1SL3", xSL1SL3, "xSL1SL3[nseg]/F");
00447   ttreeOutput->Branch("xSlSL3", xSlSL3, "xSlSL3[nseg]/F");
00448   ttreeOutput->Branch("dxdzSlSL3", dxdzSlSL3, "dxdzSlSL3[nseg]/F");
00449   ttreeOutput->Branch("exSlSL3", exSlSL3, "exSlSL3[nseg]/F");
00450   ttreeOutput->Branch("edxdzSlSL3", edxdzSlSL3, "edxdzSlSL3[nseg]/F");
00451   ttreeOutput->Branch("xSL3SL1", xSL3SL1, "xSL3SL1[nseg]/F");
00452   ttreeOutput->Branch("xc", xc, "xc[nseg][14]/F");
00453   ttreeOutput->Branch("yc", yc, "yc[nseg][14]/F");
00454   ttreeOutput->Branch("zc", zc, "zc[nseg][14]/F");
00455   ttreeOutput->Branch("ex", ex, "ex[nseg][14]/F");
00456   ttreeOutput->Branch("xcp", xcp, "xcp[nseg][14]/F");
00457   ttreeOutput->Branch("ycp", ycp, "ycp[nseg][14]/F");
00458   ttreeOutput->Branch("excp", excp, "excp[nseg][14]/F");
00459   ttreeOutput->Branch("eycp", eycp, "eycp[nseg][14]/F");
00460   ttreeOutput->Branch("wh", wh, "wh[nseg]/I");
00461   ttreeOutput->Branch("st", st, "st[nseg]/I");
00462   ttreeOutput->Branch("sr", sr, "sr[nseg]/I");
00463   ttreeOutput->Branch("sl", sl, "sl[nseg][14]/I");
00464   ttreeOutput->Branch("la", la, "la[nseg][14]/I");
00465 
00466 }
00467 
00468 
00469 DEFINE_EDM_PLUGIN( AlignmentAlgorithmPluginFactory, MuonDTLocalMillepedeAlgorithm, "MuonDTLocalMillepedeAlgorithm" );
00470