CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
DTMuonSLToSL Class Reference

#include <DTMuonSLToSL.h>

Inheritance diagram for DTMuonSLToSL:
DTMuonLocalAlignment

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 ()
 
- Public Member Functions inherited from DTMuonLocalAlignment
 DTMuonLocalAlignment ()
 
void initNTuples (int)
 
void setBranchAddressTree ()
 
 ~DTMuonLocalAlignment ()
 

Private Attributes

float cov [3][3]
 
float dx
 
float dz
 
float phiy
 
float ptMax
 
float ptMin
 
int srC
 
int stC
 
TTree * ttreeOutput
 
int whC
 

Additional Inherited Members

- Public Attributes inherited from DTMuonLocalAlignment
float charge
 
float dxdzSl [MAX_SEGMENT]
 
float dxdzSlSL1 [MAX_SEGMENT]
 
float dxdzSlSL3 [MAX_SEGMENT]
 
float dydzSl [MAX_SEGMENT]
 
float edxdzSl [MAX_SEGMENT]
 
float edxdzSlSL1 [MAX_SEGMENT]
 
float edxdzSlSL3 [MAX_SEGMENT]
 
float edydzSl [MAX_SEGMENT]
 
float eta
 
float ex [MAX_SEGMENT][MAX_HIT_CHAM]
 
float excp [MAX_SEGMENT][MAX_HIT_CHAM]
 
float exdxdzSl [MAX_SEGMENT]
 
float exdxdzSlSL1 [MAX_SEGMENT]
 
float exdxdzSlSL3 [MAX_SEGMENT]
 
float exSl [MAX_SEGMENT]
 
float exSlSL1 [MAX_SEGMENT]
 
float exSlSL3 [MAX_SEGMENT]
 
float eycp [MAX_SEGMENT][MAX_HIT_CHAM]
 
float eydydzSl [MAX_SEGMENT]
 
float eySl [MAX_SEGMENT]
 
TFile * f
 
int la [MAX_SEGMENT][MAX_HIT_CHAM]
 
int nhits [MAX_SEGMENT]
 
int nphihits [MAX_SEGMENT]
 
int nseg
 
int nthetahits [MAX_SEGMENT]
 
std::string ntuplePath
 
int numberOfRootFiles
 
float p
 
float phi
 
float pt
 
int sl [MAX_SEGMENT][MAX_HIT_CHAM]
 
int sr [MAX_SEGMENT]
 
int st [MAX_SEGMENT]
 
TChain * tali
 
int wh [MAX_SEGMENT]
 
float xc [MAX_SEGMENT][MAX_HIT_CHAM]
 
float xcp [MAX_SEGMENT][MAX_HIT_CHAM]
 
float xSl [MAX_SEGMENT]
 
float xSL1SL3 [MAX_SEGMENT]
 
float xSL3SL1 [MAX_SEGMENT]
 
float xSlSL1 [MAX_SEGMENT]
 
float xSlSL3 [MAX_SEGMENT]
 
float yc [MAX_SEGMENT][MAX_HIT_CHAM]
 
float ycp [MAX_SEGMENT][MAX_HIT_CHAM]
 
float ySl [MAX_SEGMENT]
 
float zc [MAX_SEGMENT][MAX_HIT_CHAM]
 

Detailed Description

$Date$

Revision:
1.3
Author
Luca Scodellaro Luca..nosp@m.Scod.nosp@m.ellar.nosp@m.o@ce.nosp@m.rn.ch

Definition at line 18 of file DTMuonSLToSL.h.

Constructor & Destructor Documentation

DTMuonSLToSL::DTMuonSLToSL ( std::string  path,
int  n_files,
float  MaxPt,
float  MinPt,
TFile *  f_ 
)
DTMuonSLToSL::~DTMuonSLToSL ( )

Definition at line 22 of file DTMuonSLToSL.cc.

22 {}

Member Function Documentation

void DTMuonSLToSL::calculationSLToSL ( )

Definition at line 26 of file DTMuonSLToSL.cc.

References trackerHits::c, cmsDriverOptions::counter, cov, dx, DTMuonLocalAlignment::dxdzSlSL1, DTMuonLocalAlignment::dxdzSlSL3, dz, DTMuonLocalAlignment::f, i, DTMuonLocalAlignment::nhits, DTMuonLocalAlignment::nseg, phiy, DTMuonLocalAlignment::pt, ptMax, ptMin, returnbSLMatrix(), returnCSLMatrix(), asciidump::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().

26  {
27 
28 
29  TMatrixD ****C13 = new TMatrixD ***[5];
30  TMatrixD ****b13 = new TMatrixD ***[5];
31  TMatrixD ****C31 = new TMatrixD ***[5];
32  TMatrixD ****b31 = new TMatrixD ***[5];
33 
34  for(int whI = -2; whI < 3; ++whI) {
35  C13[whI+2] = new TMatrixD **[4];
36  b13[whI+2] = new TMatrixD **[4];
37  C31[whI+2] = new TMatrixD **[4];
38  b31[whI+2] = new TMatrixD **[4];
39  for(int stI = 1; stI < 5; ++stI) {
40  C13[whI+2][stI-1] = new TMatrixD * [14];
41  b13[whI+2][stI-1] = new TMatrixD * [14];
42  C31[whI+2][stI-1] = new TMatrixD * [14];
43  b31[whI+2][stI-1] = new TMatrixD * [14];
44  for(int seI = 1; seI < 15; ++seI) {
45  if(seI > 12 && stI != 4) continue;
46  C13[whI+2][stI-1][seI-1] = new TMatrixD(3,3);
47  b13[whI+2][stI-1][seI-1] = new TMatrixD(3,1);
48  C31[whI+2][stI-1][seI-1] = new TMatrixD(3,3);
49  b31[whI+2][stI-1][seI-1] = new TMatrixD(3,1);
50  }
51  }
52  }
53 
54  //Run over the TTree
55  Int_t nentries = (Int_t)tali->GetEntries();
56  for (Int_t i=0;i<nentries;i++) {
57  tali->GetEntry(i);
58  //Basic cuts
59  if(pt > ptMax || pt < ptMin) continue;
60 
61  bool repeatedHits = false;
62  for(int counter = 0; counter < nseg; ++counter) {
63  //Make sure there are no repeated hits
64  for(int counterHi = 0; counterHi < nhits[counter]; counterHi++) {
65  for(int counterHj = 0; counterHj < nhits[counter]; counterHj++) {
66  if(counterHi == counterHj) continue;
67  if(zc[counter][counterHi] == zc[counter][counterHj]) {
68  repeatedHits = true;
69  }
70  }
71  }
72  if(repeatedHits == true) continue;
73 
74  float x_13 = xSlSL3[counter]; float xp_13 = xSL1SL3[counter];
75  float x_31 = xSlSL1[counter]; float xp_31 = xSL3SL1[counter];
76  //float tanphi = dxdzSl[counter];
77  float tanphi_13 = dxdzSlSL1[counter];
78  float tanphi_31 = dxdzSlSL3[counter];
79  int wheel = wh[counter];
80  int station = st[counter];
81  int sector = sr[counter];
82 
83  if(fabs(x_13-xp_13)< 3 && fabs(x_31-xp_31) && fabs(tanphi_13-tanphi_31)<0.06) {
84 
85  *(C13[wheel+2][station-1][sector-1]) += returnCSLMatrix(x_13, xp_13, tanphi_13);
86  *(b13[wheel+2][station-1][sector-1]) += returnbSLMatrix(x_13, xp_13, tanphi_13);
87 
88  *(C31[wheel+2][station-1][sector-1]) += returnCSLMatrix(x_31, xp_31, tanphi_31);
89  *(b31[wheel+2][station-1][sector-1]) += returnbSLMatrix(x_31, xp_31, tanphi_31);
90  }
91  }
92  }
93 
94  for(int wheel = -2; wheel < 3; ++wheel) {
95  for(int station = 1; station < 5; ++station) {
96  for(int sector = 1; sector < 15; ++sector) {
97  if(sector > 12 && station != 4) continue;
98  TMatrixD solution13(3,1);
99  TMatrixD solution31(3,1);
100  TMatrixD C31_copy = *(C31[wheel+2][station-1][sector-1]);
101  TMatrixD C13_copy = *(C13[wheel+2][station-1][sector-1]);
102  TMatrixD b31_copy = *(b31[wheel+2][station-1][sector-1]);
103  TMatrixD b13_copy = *(b13[wheel+2][station-1][sector-1]);
104 
105  C31_copy.Invert();
106  C13_copy.Invert();
107  solution13 = C13_copy * b13_copy;
108  solution31 = C31_copy * b31_copy;
109  whC = wheel; stC = station; srC = sector;
110  dx = solution13(0,0); dz = solution13(1,0);
111  phiy = solution13(2,0);
112  for(int c = 0; c < 3; ++c) {
113  for(int s = 0; s < 3; ++s) {
114  cov[c][s] = C13_copy(c, s);
115  }
116  }
117  ttreeOutput->Fill();
118  }
119  }
120  }
121  f->Write();
122 }
float dxdzSlSL3[MAX_SEGMENT]
int i
Definition: DBlmapReader.cc:9
float xSL3SL1[MAX_SEGMENT]
TTree * ttreeOutput
Definition: DTMuonSLToSL.h:45
float dxdzSlSL1[MAX_SEGMENT]
float xSL1SL3[MAX_SEGMENT]
float xSlSL1[MAX_SEGMENT]
float zc[MAX_SEGMENT][MAX_HIT_CHAM]
float xSlSL3[MAX_SEGMENT]
string s
Definition: asciidump.py:422
TMatrixD returnCSLMatrix(float, float, float)
float cov[3][3]
Definition: DTMuonSLToSL.h:40
TMatrixD returnbSLMatrix(float, float, float)
TMatrixD DTMuonSLToSL::returnbSLMatrix ( float  x,
float  xp,
float  tanphi 
)

Definition at line 145 of file DTMuonSLToSL.cc.

Referenced by calculationSLToSL().

145  {
146 
147  TMatrixD matrix(3,1);
148 
149  matrix(0,0) = -(x-xp);
150  matrix(1,0) = -(x-xp)*tanphi;
151  matrix(2,0) = -(x-xp)*tanphi*xp;
152 
153  return matrix;
154 
155 }
Definition: DDAxes.h:10
TMatrixD DTMuonSLToSL::returnCSLMatrix ( float  x,
float  xp,
float  tanphi 
)

Definition at line 126 of file DTMuonSLToSL.cc.

Referenced by calculationSLToSL().

126  {
127 
128  TMatrixD matrix(3,3);
129 
130  matrix(0,0) = 1.0;
131  matrix(1,0) = tanphi;
132  matrix(0,1) = tanphi;
133  matrix(1,1) = tanphi*tanphi;
134  matrix(0,2) = tanphi*xp;
135  matrix(2,0) = tanphi*xp;
136  matrix(2,2) = tanphi*tanphi*xp*xp;
137  matrix(2,1) = tanphi*tanphi*xp;
138  matrix(1,2) = tanphi*tanphi*xp;
139 
140  return matrix;
141 
142 }
void DTMuonSLToSL::setBranchTree ( )

Definition at line 158 of file DTMuonSLToSL.cc.

References cov, dx, dz, phiy, srC, stC, ttreeOutput, and whC.

Referenced by DTMuonSLToSL().

158  {
159 
160  ttreeOutput = new TTree("DTSLToSLResult", "DTSLToSLResult");
161 
162  ttreeOutput->Branch("wh", &whC, "wh/F");
163  ttreeOutput->Branch("st", &stC, "st/F");
164  ttreeOutput->Branch("sr", &srC, "sr/F");
165  ttreeOutput->Branch("dx", &dx, "dx/F");
166  ttreeOutput->Branch("dz", &dz, "dz/F");
167  ttreeOutput->Branch("phiy", &phiy, "phiy/F");
168  ttreeOutput->Branch("cov", cov, "cov[3][3]/F");
169 
170 }
TTree * ttreeOutput
Definition: DTMuonSLToSL.h:45
float cov[3][3]
Definition: DTMuonSLToSL.h:40

Member Data Documentation

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().