CMS 3D CMS Logo

TrackAlignment Class Reference

Analyzer of the StandAlone muon tracks for alignment with Tracks. More...

#include <Alignment/MuonStandaloneAlgorithm/plugins/TrackAlignment.h>

Inheritance diagram for TrackAlignment:

edm::EDAnalyzer

List of all members.

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup)
virtual void beginJob (const edm::EventSetup &eventSetup)
virtual void endJob ()
 TrackAlignment (const edm::ParameterSet &pset)
 Constructor.
virtual ~TrackAlignment ()
 Destructor.

Private Attributes

std::vector< long > alignables
TMatrixD bMatrix
TMatrixD CMatrix
AlignmentDetectorCollection detectorCollection
TFile * f
long int numberOfTracks
TFile * theFile
std::string theMatrixFile
std::string theTrackForAlignment
std::vector< TH1F * > unitsDerXK
std::vector< TH1F * > unitsDerXL
std::vector< TH1F * > unitsDerXPhi
std::vector< TH1F * > unitsDerXT
std::vector< TH1F * > unitsDerXTheta
std::vector< TH1F * > unitsDerZK
std::vector< TH1F * > unitsDerZL
std::vector< TH1F * > unitsDerZPhi
std::vector< TH1F * > unitsDerZT
std::vector< TH1F * > unitsDerZTheta
std::vector< TH1F * > unitsPhi
std::vector< TH1F * > unitsRPhi
std::vector< TH1F * > unitsTheta
std::vector< TH1F * > unitsZ

Static Private Attributes

static const int NDOFAlign = 6
static const int NDOFChamber = 4
static const int NDOFCoor = 4
static const int NDOFTrack = 5


Detailed Description

Analyzer of the StandAlone muon tracks for alignment with Tracks.

Analyzer of the StandAlone muon tracks for alignment with tracks.

Date
2008/02/15 12:26:02
Revision
1.2
Author:
P. Martinez Ruiz del Arbol - IFCA (CSIC-UC) <Pablo.Martinez@cern.ch>
Date
2008/02/26 09:33:06
Revision
1.1
Author:
P. Martinez Ruiz del Arbol - IFCA (CSIC-UC) <Pablo.Martinez@cern.ch>
Date
2008/02/15 12:26:03
Revision
1.2
Author:
P. Martinez Ruiz del Arbol - IFCA (CSIC-UC) <Pablo.Martinez@cern.ch>

Definition at line 30 of file TrackAlignment.h.


Constructor & Destructor Documentation

TrackAlignment::TrackAlignment ( const edm::ParameterSet pset  ) 

Constructor.

Definition at line 45 of file TrackAlignment.cc.

References edm::ParameterSet::getParameter(), theMatrixFile, and theTrackForAlignment.

00045                                                       {
00046   
00047   theMatrixFile = pset.getParameter<string>("MatrixFile");
00048   theTrackForAlignment = pset.getParameter<string>("TrackAlignmentCollection");
00049   
00050 }

TrackAlignment::~TrackAlignment (  )  [virtual]

Destructor.

Definition at line 54 of file TrackAlignment.cc.

00054                                {
00055 }


Member Function Documentation

void TrackAlignment::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
) [virtual]

Implements edm::EDAnalyzer.

Definition at line 95 of file TrackAlignment.cc.

References bMatrix, CMatrix, detectorCollection, edm::EventSetup::get(), AlignmentDetector::index(), NDOFAlign, NDOFChamber, numberOfTracks, AlignmentDetector::setAlignmentDetector(), theTrackForAlignment, unitsDerXK, unitsDerXL, unitsDerXPhi, unitsDerXT, unitsDerXTheta, unitsDerZK, unitsDerZL, unitsDerZPhi, unitsDerZT, unitsDerZTheta, unitsPhi, unitsRPhi, unitsTheta, unitsZ, and muonGeometry::wheel.

00095                                                                              {
00096   
00097   //Take the collection of TrackForAlignment
00098   Handle<TrackForAlignmentCollection> aliTracks;
00099   event.getByLabel(theTrackForAlignment , aliTracks);
00100   
00101   //Take the TrackingGeometry
00102   ESHandle<GlobalTrackingGeometry> theTrackingGeometry;
00103   eventSetup.get<GlobalTrackingGeometryRecord>().get(theTrackingGeometry);
00104   
00105   edm::Service<TFileService> fs;
00106 
00107 
00108   //Loop over all the Tracks
00109   for(TrackForAlignmentCollection::const_iterator myTrack = aliTracks->begin(); myTrack != aliTracks->end(); myTrack++) {
00110     
00111     if((*myTrack).trackValid == true) {
00112       
00113       numberOfTracks++;
00114       //These are the matrizes and Residuals for track j
00115       TMatrixD Cp, bp;
00116       Cp.ResizeTo(NDOFChamber*NDOFAlign,NDOFChamber*NDOFAlign);
00117       bp.ResizeTo(NDOFChamber*NDOFAlign,1);
00118       Cp = (*myTrack).CpMatrix();
00119       bp = (*myTrack).BMatrix();
00120       
00121       std::vector<int> indexOfMatrix;
00122       
00123       //Get a vector of PointsForAlignment
00124       std::vector<PointForAlignment> myPoints = (*myTrack).pointsAlignment();
00125       
00126       //Loop over the points
00127       int countPoints = 0;
00128       for(std::vector<PointForAlignment>::iterator points = myPoints.begin(); points != myPoints.end(); points++) {
00129         
00130         int position = 0;
00131         bool newDetector = true;
00132         TMatrixD Residuals = (*points).residual();
00133         TMatrixD Jacobian = (*points).derivatives();
00134         //Loop over the AligmentDetectorCollection to see if the detector is new and requires a new entry
00135         for(AlignmentDetectorCollection::const_iterator myIds = detectorCollection.begin(); myIds != detectorCollection.end(); myIds++) {
00136           //If matches newDetector = false
00137           if((*myIds).rawId() == (*points).rawId()) {
00138             position = (*myIds).index();
00139             newDetector = false;
00140             break;
00141           }
00142         }
00143         DetId myDet((*points).rawId());
00144         DTChamberId myCham((*points).rawId());
00145         int wheel = myCham.wheel();
00146         //If the detector is new
00147         if(newDetector) {
00148           //Create an AlignmentDetector, fill it and push it into the collection 
00149           AlignmentDetector *det = new AlignmentDetector();
00150           det->setAlignmentDetector((*points).rawId(), detectorCollection.size());
00151           detectorCollection.push_back(*det);
00152           //This piece of code calculates the range of the residuals
00153           int station = 0;
00154           //If it's a DT
00155           if(myDet.subdetId() == 1) {
00156             DTChamberId myChamber((*points).rawId());
00157             station = myChamber.station();
00158           } else {
00159             CSCDetId myChamber((*points).rawId());
00160             station = myChamber.ring();
00161           }
00162           
00163           //Giving names to the Histograms
00164           char nameOfHistoRPhi[50];
00165           char nameOfHistoPhi[50];
00166           char nameOfHistoZ[50];
00167           char nameOfHistoTheta[50];
00168           // Derivatives
00169           char derivativeXK[50];
00170           char derivativeXPhi[50];
00171           char derivativeXT[50];
00172           char derivativeXTheta[50];
00173           char derivativeXL[50];
00174           char derivativeZK[50];
00175           char derivativeZPhi[50];
00176           char derivativeZT[50];
00177           char derivativeZTheta[50];
00178           char derivativeZL[50];
00179           
00180           sprintf(derivativeXK, "DerivativeXK%ld", (*det).rawId());
00181           sprintf(derivativeXPhi, "DerivativeXPhi%ld", (*det).rawId());
00182           sprintf(derivativeXT, "DerivativeXT%ld", (*det).rawId());
00183           sprintf(derivativeXTheta, "DerivativeXTheta%ld", (*det).rawId());
00184           sprintf(derivativeXL, "DerivativeXL%ld", (*det).rawId());
00185           sprintf(derivativeZK, "DerivativeZK%ld", (*det).rawId());
00186           sprintf(derivativeZPhi, "DerivativeZPhi%ld", (*det).rawId());
00187           sprintf(derivativeZT, "DerivativeZT%ld", (*det).rawId());
00188           sprintf(derivativeZTheta, "DerivativeZTheta%ld", (*det).rawId());
00189           sprintf(derivativeZL, "DerivativeZL%ld", (*det).rawId());
00190           
00191 
00192           sprintf(nameOfHistoPhi, "HistoPhi%ld", (*det).rawId());
00193           sprintf(nameOfHistoRPhi, "HistoRPhi%ld", (*det).rawId());
00194           sprintf(nameOfHistoZ, "HistoZ%ld", (*det).rawId());
00195           sprintf(nameOfHistoTheta, "HistoTheta%ld", (*det).rawId());
00196           double xrange = 0, zrange= 0;
00197           double rangeK = 0, rangePhi = 0, rangeT = 0, rangeTheta = 0, rangeL = 0; 
00198           if(station == 1) {
00199             rangeK = 500.0*500;
00200             rangePhi = 500;
00201             rangeT = 5.0;
00202             rangeTheta = 500.0;
00203             rangeL = 5.0;
00204             xrange = 0.1;
00205             zrange = 0.7;
00206           } else if(station == 2) {
00207             rangeK = 650*650;
00208             rangePhi = 650;
00209             rangeT = 5.0;
00210             rangeTheta = 650.0;
00211             rangeL = 5.0;
00212             xrange = 0.3;
00213             zrange = 0.7;
00214           } else if(station == 3) {
00215             rangeK = 750*750;
00216             rangePhi = 750.0;
00217             rangeT = 5.0;
00218             rangeTheta = 750.0;
00219             rangeL = 5.0;
00220             xrange = 3;
00221             zrange = 5;
00222           } else if(station == 4) {
00223             rangeK = 850.0*850.0;
00224             rangePhi = 850.0;
00225             rangeT = 5.0;
00226             rangeTheta = 850.0;
00227             rangeL = 5.0;
00228             xrange = 3;
00229             zrange = 3;
00230           }
00231  
00232           TH1F *DerivativeXK = fs->make<TH1F>(derivativeXK, derivativeXK, 100, -rangeK, rangeK);
00233           TH1F *DerivativeXPhi = fs->make<TH1F>(derivativeXPhi, derivativeXPhi, 100, -rangePhi, rangePhi);
00234           TH1F *DerivativeXT = fs->make<TH1F>(derivativeXT, derivativeXT, 100, -rangeT, rangeT);
00235           TH1F *DerivativeXTheta = fs->make<TH1F>(derivativeXTheta, derivativeXTheta, 100, -rangeTheta, rangeTheta);
00236           TH1F *DerivativeXL = fs->make<TH1F>(derivativeXL, derivativeXL, 100, -rangeL, rangeL);
00237           
00238           TH1F *DerivativeZK = fs->make<TH1F>(derivativeZK, derivativeZK, 100, -rangeK, rangeK);
00239           TH1F *DerivativeZPhi = fs->make<TH1F>(derivativeZPhi, derivativeZPhi, 100, -rangePhi, rangePhi);
00240           TH1F *DerivativeZT = fs->make<TH1F>(derivativeZT, derivativeZT, 100, -rangeT, rangeT);
00241           TH1F *DerivativeZTheta = fs->make<TH1F>(derivativeZTheta, derivativeZTheta, 100, -rangeTheta, rangeTheta);
00242           TH1F *DerivativeZL = fs->make<TH1F>(derivativeZL, derivativeZL, 100, -rangeL, rangeL);
00243  
00244           TH1F *histoPhi = fs->make<TH1F>(nameOfHistoPhi, nameOfHistoPhi, 500, -0.1, 0.1);
00245           TH1F *histoRPhi = fs->make<TH1F>(nameOfHistoRPhi, nameOfHistoRPhi, 100, -xrange, xrange);
00246           TH1F *histoZ = fs->make<TH1F>(nameOfHistoZ, nameOfHistoZ, 500, -zrange, zrange);
00247           TH1F *histoTheta = fs->make<TH1F>(nameOfHistoTheta, nameOfHistoTheta, 500, -0.1, 0.1);
00248           
00249           DerivativeXK->Fill(Jacobian(0,0));
00250           DerivativeXPhi->Fill(Jacobian(0,1));
00251           DerivativeXT->Fill(Jacobian(0,3));
00252           DerivativeXTheta->Fill(Jacobian(0,2));
00253           DerivativeXL->Fill(Jacobian(0,4));
00254           
00255           DerivativeZK->Fill(Jacobian(1,0));
00256           DerivativeZPhi->Fill(Jacobian(1,1));
00257           DerivativeZT->Fill(Jacobian(1,3));
00258           DerivativeZTheta->Fill(Jacobian(1,2));
00259           DerivativeZL->Fill(Jacobian(1,4));
00260           
00261           histoRPhi->Fill(Residuals(0,0));
00262           histoZ->Fill(Residuals(1,0));
00263           histoPhi->Fill(Residuals(2,0));
00264           histoTheta->Fill(Residuals(3,0));
00265           
00266           
00267           unitsDerXK.push_back(DerivativeXK);
00268           unitsDerXPhi.push_back(DerivativeXPhi);
00269           unitsDerXT.push_back(DerivativeXT);
00270           unitsDerXTheta.push_back(DerivativeXTheta);
00271           unitsDerXL.push_back(DerivativeXL);
00272           unitsDerZK.push_back(DerivativeZK);
00273           unitsDerZPhi.push_back(DerivativeZPhi);
00274           unitsDerZT.push_back(DerivativeZT);
00275           unitsDerZTheta.push_back(DerivativeZTheta);
00276           unitsDerZL.push_back(DerivativeZL);
00277           
00278 
00279           //Push them into their respective vectors
00280           unitsRPhi.push_back(histoRPhi);
00281           unitsZ.push_back(histoZ);
00282           unitsTheta.push_back(histoTheta);
00283           unitsPhi.push_back(histoPhi);
00284           
00285           //If the detector is new I update the dimension of the Matrix
00286           bMatrix.ResizeTo(bMatrix.GetNrows()+NDOFAlign,1);
00287           CMatrix.ResizeTo(CMatrix.GetNrows()+NDOFAlign,CMatrix.GetNcols()+NDOFAlign);
00288           for(int NDOFcounter = 0; NDOFcounter < NDOFAlign; ++NDOFcounter) {
00289             bMatrix(bMatrix.GetNrows()+NDOFcounter-NDOFAlign,0) = 0.0;
00290             for(int zeroCounter = 0; zeroCounter < CMatrix.GetNrows(); zeroCounter++) {
00291               CMatrix(CMatrix.GetNrows()+NDOFcounter-NDOFAlign, zeroCounter) = 0.0;
00292               CMatrix(zeroCounter, CMatrix.GetNrows()+NDOFcounter-NDOFAlign) = 0.0;
00293             }   
00294           }
00295           indexOfMatrix.push_back(det->index());
00296         } else {
00297           unitsDerXK.at(position)->Fill(Jacobian(0,0));
00298           unitsDerXPhi.at(position)->Fill(Jacobian(0,1));
00299           unitsDerXT.at(position)->Fill(Jacobian(0,3));
00300           unitsDerXTheta.at(position)->Fill(Jacobian(0,2));
00301           unitsDerXL.at(position)->Fill(Jacobian(0,4));
00302           unitsDerZK.at(position)->Fill(Jacobian(1,0));
00303           unitsDerZPhi.at(position)->Fill(Jacobian(1,1));
00304           unitsDerZT.at(position)->Fill(Jacobian(1,3));
00305           unitsDerZTheta.at(position)->Fill(Jacobian(1,2));
00306           unitsDerZL.at(position)->Fill(Jacobian(1,4));
00307           
00308           //If the detector was not new, just fill the histogram
00309           unitsRPhi.at(position)->Fill(Residuals(0,0));
00310           unitsZ.at(position)->Fill(Residuals(1,0));
00311           unitsPhi.at(position)->Fill(Residuals(2,0));
00312           unitsTheta.at(position)->Fill(Residuals(3,0));
00313           indexOfMatrix.push_back(position);
00314         }
00315         countPoints++;
00316       }
00317       int rowCounter = 0, columnCounter = 0;
00318       for(int NDOFcounter = 0; NDOFcounter < NDOFAlign; ++NDOFcounter) {
00319         for(int NDOFcounter2 = 0; NDOFcounter2 < NDOFAlign; ++NDOFcounter2) {
00320           rowCounter = 0;
00321           columnCounter = 0;
00322           for(std::vector<int>::iterator myIndex = indexOfMatrix.begin(); myIndex != indexOfMatrix.end(); ++myIndex){
00323             for(std::vector<int>::iterator myIndex2 = indexOfMatrix.begin(); myIndex2 != indexOfMatrix.end(); ++myIndex2){
00324               CMatrix(NDOFAlign*(*myIndex)+NDOFcounter, NDOFAlign*(*myIndex2)+NDOFcounter2) += Cp(NDOFAlign*rowCounter+NDOFcounter, NDOFAlign*columnCounter+NDOFcounter2);
00325               columnCounter++;
00326             }
00327             columnCounter = 0;
00328             bMatrix(NDOFAlign*(*myIndex)+NDOFcounter, 0) += bp(NDOFAlign*rowCounter+NDOFcounter, 0)/NDOFAlign;
00329             rowCounter++;
00330           }
00331         }
00332       }
00333     }
00334   }
00335 }

void TrackAlignment::beginJob ( const edm::EventSetup eventSetup  )  [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 59 of file TrackAlignment.cc.

References bMatrix, CMatrix, and numberOfTracks.

00059                                                          {
00060   
00061   //Set the size of Matrizes to 0
00062   CMatrix.ResizeTo(0,0);
00063   bMatrix.ResizeTo(0,0);
00064   numberOfTracks = 0; 
00065  
00066 }

void TrackAlignment::endJob ( void   )  [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 69 of file TrackAlignment.cc.

References bMatrix, CMatrix, counter(), detectorCollection, N, numberOfTracks, and theFile.

00069                            {
00070   
00071   //Output root file is opened for writting
00072   TMatrixD N(1,1);
00073   N(0,0) = numberOfTracks;  
00074 
00075   theFile = new TFile("SubMatrices.root", "RECREATE");
00076   theFile->cd();
00077   
00078   N.Write("numberOfTracks");
00079   CMatrix.Write("CMatrix");
00080   bMatrix.Write("bMatrix");
00081   //Save the index matrix
00082   TMatrixD realIndexMatrix(detectorCollection.size(), 1);
00083   int counter = 0;
00084   for(AlignmentDetectorCollection::const_iterator myIds = detectorCollection.begin(); myIds != detectorCollection.end(); myIds++) {
00085     realIndexMatrix(counter, 0) = (*myIds).rawId();
00086     counter++;
00087   }
00088   realIndexMatrix.Write("Index");
00089   theFile->Close();
00090   
00091 }


Member Data Documentation

std::vector<long> TrackAlignment::alignables [private]

Definition at line 70 of file TrackAlignment.h.

TMatrixD TrackAlignment::bMatrix [private]

Definition at line 72 of file TrackAlignment.h.

Referenced by analyze(), beginJob(), and endJob().

TMatrixD TrackAlignment::CMatrix [private]

Definition at line 71 of file TrackAlignment.h.

Referenced by analyze(), beginJob(), and endJob().

AlignmentDetectorCollection TrackAlignment::detectorCollection [private]

Definition at line 67 of file TrackAlignment.h.

Referenced by analyze(), and endJob().

TFile* TrackAlignment::f [private]

Definition at line 50 of file TrackAlignment.h.

const int TrackAlignment::NDOFAlign = 6 [static, private]

Definition at line 77 of file TrackAlignment.h.

Referenced by analyze().

const int TrackAlignment::NDOFChamber = 4 [static, private]

Definition at line 78 of file TrackAlignment.h.

Referenced by analyze().

const int TrackAlignment::NDOFCoor = 4 [static, private]

Definition at line 79 of file TrackAlignment.h.

const int TrackAlignment::NDOFTrack = 5 [static, private]

Definition at line 76 of file TrackAlignment.h.

long int TrackAlignment::numberOfTracks [private]

Definition at line 74 of file TrackAlignment.h.

Referenced by analyze(), beginJob(), and endJob().

TFile* TrackAlignment::theFile [private]

Definition at line 49 of file TrackAlignment.h.

Referenced by endJob().

std::string TrackAlignment::theMatrixFile [private]

Definition at line 68 of file TrackAlignment.h.

Referenced by TrackAlignment().

std::string TrackAlignment::theTrackForAlignment [private]

Definition at line 69 of file TrackAlignment.h.

Referenced by analyze(), and TrackAlignment().

std::vector<TH1F *> TrackAlignment::unitsDerXK [private]

Definition at line 51 of file TrackAlignment.h.

Referenced by analyze().

std::vector<TH1F *> TrackAlignment::unitsDerXL [private]

Definition at line 55 of file TrackAlignment.h.

Referenced by analyze().

std::vector<TH1F *> TrackAlignment::unitsDerXPhi [private]

Definition at line 52 of file TrackAlignment.h.

Referenced by analyze().

std::vector<TH1F *> TrackAlignment::unitsDerXT [private]

Definition at line 53 of file TrackAlignment.h.

Referenced by analyze().

std::vector<TH1F *> TrackAlignment::unitsDerXTheta [private]

Definition at line 54 of file TrackAlignment.h.

Referenced by analyze().

std::vector<TH1F *> TrackAlignment::unitsDerZK [private]

Definition at line 56 of file TrackAlignment.h.

Referenced by analyze().

std::vector<TH1F *> TrackAlignment::unitsDerZL [private]

Definition at line 60 of file TrackAlignment.h.

Referenced by analyze().

std::vector<TH1F *> TrackAlignment::unitsDerZPhi [private]

Definition at line 57 of file TrackAlignment.h.

Referenced by analyze().

std::vector<TH1F *> TrackAlignment::unitsDerZT [private]

Definition at line 58 of file TrackAlignment.h.

Referenced by analyze().

std::vector<TH1F *> TrackAlignment::unitsDerZTheta [private]

Definition at line 59 of file TrackAlignment.h.

Referenced by analyze().

std::vector<TH1F *> TrackAlignment::unitsPhi [private]

Definition at line 63 of file TrackAlignment.h.

Referenced by analyze().

std::vector<TH1F *> TrackAlignment::unitsRPhi [private]

Definition at line 62 of file TrackAlignment.h.

Referenced by analyze().

std::vector<TH1F *> TrackAlignment::unitsTheta [private]

Definition at line 65 of file TrackAlignment.h.

Referenced by analyze().

std::vector<TH1F *> TrackAlignment::unitsZ [private]

Definition at line 64 of file TrackAlignment.h.

Referenced by analyze().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:33:57 2009 for CMSSW by  doxygen 1.5.4