#include <Alignment/MuonStandaloneAlgorithm/plugins/AlignmentEstimator.h>
Public Member Functions | |
AlignmentEstimator (const edm::ParameterSet &pset) | |
Constructor. | |
void | analyze (const edm::Event &event, const edm::EventSetup &eventSetup) |
virtual void | beginJob (const edm::EventSetup &eventSetup) |
virtual void | endJob () |
virtual | ~AlignmentEstimator () |
Destructor. | |
Private Member Functions | |
int | getIndex (long) |
void | getMatrix () |
Private Attributes | |
TMatrixD | b |
TMatrixD | C |
TMatrixD | CMatrixInvert |
std::vector< TH1F * > | estimatorHist |
TFile * | f |
std::string | nameOfMatrix |
double | numberOfTracks |
TFile * | theFile |
TMatrixD | theIndex |
std::string | theTrackForAlignment |
bool | theWriteToDB |
Static Private Attributes | |
static const int | NDOFAlign = 6 |
static const int | NDOFChamber = 4 |
static const int | NDOFCoor = 4 |
static const int | NDOFTrack = 5 |
Definition at line 30 of file AlignmentEstimator.h.
AlignmentEstimator::AlignmentEstimator | ( | const edm::ParameterSet & | pset | ) |
Constructor.
Definition at line 45 of file AlignmentEstimator.cc.
References edm::ParameterSet::getParameter(), nameOfMatrix, theTrackForAlignment, and theWriteToDB.
00045 { 00046 00047 theTrackForAlignment = pset.getParameter<string>("TrackAlignmentCollection"); 00048 theWriteToDB = pset.getParameter<bool>("WriteToDB"); 00049 nameOfMatrix = pset.getParameter<string>("NameOfMatrix"); 00050 00051 }
AlignmentEstimator::~AlignmentEstimator | ( | ) | [virtual] |
void AlignmentEstimator::analyze | ( | const edm::Event & | event, | |
const edm::EventSetup & | eventSetup | |||
) | [virtual] |
Implements edm::EDAnalyzer.
Definition at line 140 of file AlignmentEstimator.cc.
References b, CMatrixInvert, counter(), HLT_VtxMuL3::estimator, estimatorHist, edm::EventSetup::get(), getIndex(), NDOFAlign, numberOfTracks, theIndex, and theTrackForAlignment.
00140 { 00141 00142 //Take the collection of TrackForAlignment 00143 Handle<TrackForAlignmentCollection> aliTracks; 00144 event.getByLabel(theTrackForAlignment , aliTracks); 00145 00146 //Take the TrackingGeometry 00147 ESHandle<GlobalTrackingGeometry> theTrackingGeometry; 00148 eventSetup.get<GlobalTrackingGeometryRecord>().get(theTrackingGeometry); 00149 00150 //Still not ready for admiting several degrees of freedom 00151 for(TrackForAlignmentCollection::const_iterator myTrack = aliTracks->begin(); myTrack != aliTracks->end(); myTrack++) { 00152 if((*myTrack).trackValid == true) { 00153 //Coger puntos 00154 //Calcular b y ponerlo en orden 00155 //Multiplicar 00156 TMatrixD bReduced = (*myTrack).BMatrix(); 00157 TMatrixD bComplete(b.GetNrows(), 1); 00158 00159 int counter = 0; 00160 std::vector<PointForAlignment> myPoints = (*myTrack).pointsAlignment(); 00161 for(std::vector<PointForAlignment>::iterator points = myPoints.begin(); points != myPoints.end(); ++points) { 00162 int indexPos = getIndex((*points).rawId()); 00163 for(int parCounter = 0; parCounter < NDOFAlign; parCounter++) { 00164 bComplete(indexPos*NDOFAlign+parCounter, 0) = bReduced(counter*NDOFAlign+parCounter, 0); 00165 } 00166 } 00167 00168 TMatrixD estimator = numberOfTracks*CMatrixInvert*bComplete; 00169 for(int countAlignable = 0; countAlignable < theIndex.GetNrows(); ++countAlignable) { 00170 for(int parCounter = 0; parCounter < NDOFAlign; ++parCounter) { 00171 estimatorHist.at(countAlignable*NDOFAlign+parCounter)->Fill(estimator(countAlignable*NDOFAlign+parCounter, 0)); 00172 } 00173 } 00174 } 00175 } 00176 }
void AlignmentEstimator::beginJob | ( | const edm::EventSetup & | eventSetup | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 86 of file AlignmentEstimator.cc.
References b, C, counter(), estimatorHist, getMatrix(), NDOFAlign, and theIndex.
00086 { 00087 00088 00089 //Set the size of Matrizes to 0 00090 C.ResizeTo(0,0); 00091 b.ResizeTo(0,0); 00092 theIndex.ResizeTo(0,0); 00093 00094 getMatrix(); // This method reads matrizes in C and alignables vector 00095 00096 edm::Service<TFileService> fs; 00097 00098 //Have to initialize the histograms 00099 for(int counter = 0; counter < theIndex.GetNrows(); ++counter) { 00100 for(int parCounter = 0; parCounter < NDOFAlign; ++parCounter) { 00101 char nameOfDegree[20]; 00102 switch(parCounter) { 00103 case 0: 00104 sprintf(nameOfDegree, "delta_x"); 00105 break; 00106 case 1: 00107 sprintf(nameOfDegree, "delta_r"); 00108 break; 00109 case 2: 00110 sprintf(nameOfDegree, "delta_z"); 00111 break; 00112 case 3: 00113 sprintf(nameOfDegree, "delta_alpha"); 00114 break; 00115 case 4: 00116 sprintf(nameOfDegree, "delta_beta"); 00117 break; 00118 case 5: 00119 sprintf(nameOfDegree, "delta_gamma"); 00120 break; 00121 default: 00122 break; 00123 } 00124 char nameOfEstimator[80]; 00125 sprintf(nameOfEstimator, "Estimator_%ld_%s", (long)theIndex(counter,0), nameOfDegree); 00126 TH1F *EstHist = fs->make<TH1F>(nameOfEstimator, nameOfEstimator, 100, -0.2, 0.2); 00127 estimatorHist.push_back(EstHist); 00128 } 00129 } 00130 }
int AlignmentEstimator::getIndex | ( | long | index | ) | [private] |
void AlignmentEstimator::getMatrix | ( | ) | [private] |
Definition at line 59 of file AlignmentEstimator.cc.
References b, C, CMatrixInvert, f, in, N, numberOfTracks, and theIndex.
Referenced by beginJob().
00059 { 00060 00061 f = TFile::Open("rfio:/dpm/ifca.es/home/cms/store/MuonStandaloneAlgorithm/Matrix.root", "READ"); 00062 TMatrixD *Cm = (TMatrixD *)f->Get("CMatrix"); 00063 TMatrixD *CmInvert = (TMatrixD *)f->Get("CMatrixInvert"); 00064 TMatrixD *bm = (TMatrixD *)f->Get("bMatrix"); 00065 TMatrixD *in = (TMatrixD *)f->Get("Index"); 00066 TMatrixD *N = (TMatrixD *)f->Get("numberOfTracks"); 00067 00068 numberOfTracks = (*N)(0,0); 00069 00070 C.ResizeTo(Cm->GetNrows(), Cm->GetNcols()); 00071 CMatrixInvert.ResizeTo(Cm->GetNrows(), Cm->GetNcols()); 00072 b.ResizeTo(bm->GetNrows(), 1); 00073 theIndex.ResizeTo(in->GetNrows(), 1); 00074 00075 C = *Cm; 00076 CMatrixInvert = *CmInvert; 00077 b = *bm; 00078 theIndex = *in; 00079 00080 f->Close(); 00081 delete f; 00082 00083 }
TMatrixD AlignmentEstimator::b [private] |
Definition at line 61 of file AlignmentEstimator.h.
Referenced by analyze(), beginJob(), and getMatrix().
TMatrixD AlignmentEstimator::C [private] |
TMatrixD AlignmentEstimator::CMatrixInvert [private] |
std::vector<TH1F *> AlignmentEstimator::estimatorHist [private] |
TFile* AlignmentEstimator::f [private] |
std::string AlignmentEstimator::nameOfMatrix [private] |
const int AlignmentEstimator::NDOFAlign = 6 [static, private] |
const int AlignmentEstimator::NDOFChamber = 4 [static, private] |
Definition at line 68 of file AlignmentEstimator.h.
const int AlignmentEstimator::NDOFCoor = 4 [static, private] |
Definition at line 69 of file AlignmentEstimator.h.
const int AlignmentEstimator::NDOFTrack = 5 [static, private] |
Definition at line 66 of file AlignmentEstimator.h.
double AlignmentEstimator::numberOfTracks [private] |
TFile* AlignmentEstimator::theFile [private] |
Definition at line 53 of file AlignmentEstimator.h.
TMatrixD AlignmentEstimator::theIndex [private] |
Definition at line 62 of file AlignmentEstimator.h.
Referenced by analyze(), beginJob(), getIndex(), and getMatrix().
std::string AlignmentEstimator::theTrackForAlignment [private] |
Definition at line 57 of file AlignmentEstimator.h.
Referenced by AlignmentEstimator(), and analyze().
bool AlignmentEstimator::theWriteToDB [private] |