#include <DTMuonSLToSL.h>
Public Member Functions | |
void | calculationSLToSL () |
DTMuonSLToSL (std::string, int, float, float, TFile *) | |
TMatrixD | returnbSLMatrix (float, float, float) |
TMatrixD | returnCSLMatrix (float, float, float) |
void | setBranchTree () |
~DTMuonSLToSL () | |
Private Attributes | |
float | cov [3][3] |
float | dx |
float | dz |
float | phiy |
float | ptMax |
float | ptMin |
int | srC |
int | stC |
TTree * | ttreeOutput |
int | whC |
DTMuonSLToSL::DTMuonSLToSL | ( | std::string | path, |
int | n_files, | ||
float | MaxPt, | ||
float | MinPt, | ||
TFile * | f_ | ||
) |
Definition at line 4 of file DTMuonSLToSL.cc.
References calculationSLToSL(), DTMuonLocalAlignment::f, DTMuonLocalAlignment::initNTuples(), step1_ZMM_7Tev::MinPt, DTMuonLocalAlignment::ntuplePath, DTMuonLocalAlignment::numberOfRootFiles, scaleCards::path, ptMax, ptMin, and setBranchTree().
{ ntuplePath = path; numberOfRootFiles = n_files; f = f_; ptMax = MaxPt; ptMin = MinPt; setBranchTree(); initNTuples(0); calculationSLToSL(); }
DTMuonSLToSL::~DTMuonSLToSL | ( | ) |
Definition at line 22 of file DTMuonSLToSL.cc.
{}
void DTMuonSLToSL::calculationSLToSL | ( | ) |
Definition at line 26 of file DTMuonSLToSL.cc.
References trackerHits::c, cov, dx, DTMuonLocalAlignment::dxdzSlSL1, DTMuonLocalAlignment::dxdzSlSL3, dz, DTMuonLocalAlignment::f, i, DTMuonLocalAlignment::nhits, DTMuonLocalAlignment::nseg, phiy, DTMuonLocalAlignment::pt, ptMax, ptMin, returnbSLMatrix(), returnCSLMatrix(), alignCSCRings::s, DTMuonLocalAlignment::sr, srC, DTMuonLocalAlignment::st, relativeConstraints::station, stC, DTMuonLocalAlignment::tali, ttreeOutput, DTMuonLocalAlignment::wh, whC, DTMuonLocalAlignment::xSL1SL3, DTMuonLocalAlignment::xSL3SL1, DTMuonLocalAlignment::xSlSL1, DTMuonLocalAlignment::xSlSL3, and DTMuonLocalAlignment::zc.
Referenced by DTMuonSLToSL().
{ TMatrixD ****C13 = new TMatrixD ***[5]; TMatrixD ****b13 = new TMatrixD ***[5]; TMatrixD ****C31 = new TMatrixD ***[5]; TMatrixD ****b31 = new TMatrixD ***[5]; for(int whI = -2; whI < 3; ++whI) { C13[whI+2] = new TMatrixD **[4]; b13[whI+2] = new TMatrixD **[4]; C31[whI+2] = new TMatrixD **[4]; b31[whI+2] = new TMatrixD **[4]; for(int stI = 1; stI < 5; ++stI) { C13[whI+2][stI-1] = new TMatrixD * [14]; b13[whI+2][stI-1] = new TMatrixD * [14]; C31[whI+2][stI-1] = new TMatrixD * [14]; b31[whI+2][stI-1] = new TMatrixD * [14]; for(int seI = 1; seI < 15; ++seI) { if(seI > 12 && stI != 4) continue; C13[whI+2][stI-1][seI-1] = new TMatrixD(3,3); b13[whI+2][stI-1][seI-1] = new TMatrixD(3,1); C31[whI+2][stI-1][seI-1] = new TMatrixD(3,3); b31[whI+2][stI-1][seI-1] = new TMatrixD(3,1); } } } //Run over the TTree Int_t nentries = (Int_t)tali->GetEntries(); for (Int_t i=0;i<nentries;i++) { tali->GetEntry(i); //Basic cuts if(pt > ptMax || pt < ptMin) continue; bool repeatedHits = false; for(int counter = 0; counter < nseg; ++counter) { //Make sure there are no repeated hits for(int counterHi = 0; counterHi < nhits[counter]; counterHi++) { for(int counterHj = 0; counterHj < nhits[counter]; counterHj++) { if(counterHi == counterHj) continue; if(zc[counter][counterHi] == zc[counter][counterHj]) { repeatedHits = true; } } } if(repeatedHits == true) continue; float x_13 = xSlSL3[counter]; float xp_13 = xSL1SL3[counter]; float x_31 = xSlSL1[counter]; float xp_31 = xSL3SL1[counter]; //float tanphi = dxdzSl[counter]; float tanphi_13 = dxdzSlSL1[counter]; float tanphi_31 = dxdzSlSL3[counter]; int wheel = wh[counter]; int station = st[counter]; int sector = sr[counter]; if(fabs(x_13-xp_13)< 3 && fabs(x_31-xp_31) && fabs(tanphi_13-tanphi_31)<0.06) { *(C13[wheel+2][station-1][sector-1]) += returnCSLMatrix(x_13, xp_13, tanphi_13); *(b13[wheel+2][station-1][sector-1]) += returnbSLMatrix(x_13, xp_13, tanphi_13); *(C31[wheel+2][station-1][sector-1]) += returnCSLMatrix(x_31, xp_31, tanphi_31); *(b31[wheel+2][station-1][sector-1]) += returnbSLMatrix(x_31, xp_31, tanphi_31); } } } for(int wheel = -2; wheel < 3; ++wheel) { for(int station = 1; station < 5; ++station) { for(int sector = 1; sector < 15; ++sector) { if(sector > 12 && station != 4) continue; TMatrixD solution13(3,1); TMatrixD solution31(3,1); TMatrixD C31_copy = *(C31[wheel+2][station-1][sector-1]); TMatrixD C13_copy = *(C13[wheel+2][station-1][sector-1]); TMatrixD b31_copy = *(b31[wheel+2][station-1][sector-1]); TMatrixD b13_copy = *(b13[wheel+2][station-1][sector-1]); C31_copy.Invert(); C13_copy.Invert(); solution13 = C13_copy * b13_copy; solution31 = C31_copy * b31_copy; whC = wheel; stC = station; srC = sector; dx = solution13(0,0); dz = solution13(1,0); phiy = solution13(2,0); for(int c = 0; c < 3; ++c) { for(int s = 0; s < 3; ++s) { cov[c][s] = C13_copy(c, s); } } ttreeOutput->Fill(); } } } f->Write(); }
TMatrixD DTMuonSLToSL::returnbSLMatrix | ( | float | x, |
float | xp, | ||
float | tanphi | ||
) |
Definition at line 145 of file DTMuonSLToSL.cc.
References makeMuonMisalignmentScenario::matrix.
Referenced by calculationSLToSL().
TMatrixD DTMuonSLToSL::returnCSLMatrix | ( | float | x, |
float | xp, | ||
float | tanphi | ||
) |
Definition at line 126 of file DTMuonSLToSL.cc.
References makeMuonMisalignmentScenario::matrix.
Referenced by calculationSLToSL().
void DTMuonSLToSL::setBranchTree | ( | ) |
Definition at line 158 of file DTMuonSLToSL.cc.
References cov, dx, dz, phiy, srC, stC, ttreeOutput, and whC.
Referenced by DTMuonSLToSL().
{ ttreeOutput = new TTree("DTSLToSLResult", "DTSLToSLResult"); ttreeOutput->Branch("wh", &whC, "wh/F"); ttreeOutput->Branch("st", &stC, "st/F"); ttreeOutput->Branch("sr", &srC, "sr/F"); ttreeOutput->Branch("dx", &dx, "dx/F"); ttreeOutput->Branch("dz", &dz, "dz/F"); ttreeOutput->Branch("phiy", &phiy, "phiy/F"); ttreeOutput->Branch("cov", cov, "cov[3][3]/F"); }
float DTMuonSLToSL::cov[3][3] [private] |
Definition at line 40 of file DTMuonSLToSL.h.
Referenced by calculationSLToSL(), and setBranchTree().
float DTMuonSLToSL::dx [private] |
Definition at line 39 of file DTMuonSLToSL.h.
Referenced by calculationSLToSL(), and setBranchTree().
float DTMuonSLToSL::dz [private] |
Definition at line 39 of file DTMuonSLToSL.h.
Referenced by calculationSLToSL(), and setBranchTree().
float DTMuonSLToSL::phiy [private] |
Definition at line 39 of file DTMuonSLToSL.h.
Referenced by calculationSLToSL(), and setBranchTree().
float DTMuonSLToSL::ptMax [private] |
Definition at line 43 of file DTMuonSLToSL.h.
Referenced by calculationSLToSL(), and DTMuonSLToSL().
float DTMuonSLToSL::ptMin [private] |
Definition at line 43 of file DTMuonSLToSL.h.
Referenced by calculationSLToSL(), and DTMuonSLToSL().
int DTMuonSLToSL::srC [private] |
Definition at line 38 of file DTMuonSLToSL.h.
Referenced by calculationSLToSL(), and setBranchTree().
int DTMuonSLToSL::stC [private] |
Definition at line 38 of file DTMuonSLToSL.h.
Referenced by calculationSLToSL(), and setBranchTree().
TTree* DTMuonSLToSL::ttreeOutput [private] |
Definition at line 45 of file DTMuonSLToSL.h.
Referenced by calculationSLToSL(), and setBranchTree().
int DTMuonSLToSL::whC [private] |
Definition at line 38 of file DTMuonSLToSL.h.
Referenced by calculationSLToSL(), and setBranchTree().