00001 #include "Alignment/MuonAlignmentAlgorithms/interface/DTMuonMillepede.h"
00002
00003
00004 #include <iostream>
00005
00006 DTMuonMillepede::DTMuonMillepede(std::string path, int n_files,
00007 float MaxPt, float MinPt,
00008 int nPhihits, int nThetahits,
00009 int workingmode, int nMtxSection) {
00010
00011 ntuplePath = path;
00012 numberOfRootFiles = n_files;
00013
00014 ptMax = MaxPt;
00015 ptMin = MinPt;
00016
00017 nPhiHits = nPhihits;
00018 nThetaHits = nThetahits;
00019
00020 TDirectory * dirSave = gDirectory;
00021
00022
00023 myPG = new ReadPGInfo("./InternalData2009.root");
00024
00025 f = new TFile("./LocalMillepedeResults.root", "RECREATE");
00026 f->cd();
00027
00028 setBranchTree();
00029
00030 initNTuples(nMtxSection);
00031
00032 calculationMillepede(workingmode);
00033
00034 f->Write();
00035 f->Close();
00036
00037 dirSave->cd();
00038 }
00039
00040
00041 DTMuonMillepede::~DTMuonMillepede() {}
00042
00043
00044
00045 void DTMuonMillepede::calculationMillepede(int workingmode) {
00046
00047
00048 TMatrixD ****C = new TMatrixD ***[5];
00049 TMatrixD ****b = new TMatrixD ***[5];
00050
00051 for(int whI = -2; whI < 3; ++whI) {
00052 C[whI+2] = new TMatrixD **[4];
00053 b[whI+2] = new TMatrixD **[4];
00054 for(int stI = 1; stI < 5; ++stI) {
00055 C[whI+2][stI-1] = new TMatrixD * [14];
00056 b[whI+2][stI-1] = new TMatrixD * [14];
00057 for(int seI = 1; seI < 15; ++seI) {
00058 if(seI > 12 && stI != 4) continue;
00059 if(stI == 4) {
00060 C[whI+2][stI-1][seI-1] = new TMatrixD(24,24);
00061 b[whI+2][stI-1][seI-1] = new TMatrixD(24,1);
00062 } else {
00063 C[whI+2][stI-1][seI-1] = new TMatrixD(60,60);
00064 b[whI+2][stI-1][seI-1] = new TMatrixD(60,1);
00065 }
00066 }
00067 }
00068 }
00069
00070
00071 if (workingmode<=3) {
00072
00073 Int_t nentries = (Int_t)tali->GetEntries();
00074 for (Int_t i=0;i<nentries;i++) {
00075 tali->GetEntry(i);
00076
00077 if (i%100000==0) std::cout << "Analyzing track number " << i << std::endl;
00078
00079
00080 if(pt > ptMax || pt < ptMin) continue;
00081
00082 for(int counter = 0; counter < nseg; ++counter) {
00083
00084 if (nphihits[counter]<nPhiHits || (nthetahits[counter]<nThetaHits && st[counter]<4)) continue;
00085
00086 TMatrixD A(12,60);
00087 TMatrixD R(12,1);
00088 TMatrixD E(12,12);
00089 TMatrixD B(12, 4);
00090 TMatrixD I(12, 12);
00091 if(st[counter] == 4) {
00092 A.ResizeTo(8,24);
00093 R.ResizeTo(8,1);
00094 E.ResizeTo(8,8);
00095 B.ResizeTo(8,2);
00096 I.ResizeTo(8,8);
00097 }
00098
00099 for(int counterHi = 0; counterHi < nhits[counter]; counterHi++) {
00100
00101 int row = 0;
00102 if(sl[counter][counterHi] == 3) {
00103 row = 4 + (la[counter][counterHi]-1);
00104 } else if(sl[counter][counterHi] == 2) {
00105 row = 8 + (la[counter][counterHi]-1);
00106 } else {
00107 row = (la[counter][counterHi]-1);
00108 }
00109
00110 float x = xc[counter][counterHi];
00111 float y = yc[counter][counterHi];
00112 float xp = xcp[counter][counterHi];
00113 float yp = ycp[counter][counterHi];
00114 float zp = zc[counter][counterHi];
00115 float error = ex[counter][counterHi];
00116 float dxdz = dxdzSl[counter];
00117 float dydz = dydzSl[counter];
00118
00119 if(st[counter]!=4) {
00120 if(sl[counter][counterHi] == 2) {
00121 A(row, row*5) = -1.0;
00122 A(row, row*5+1) = dydz;
00123 A(row, row*5+2) = y*dydz;
00124 A(row, row*5+3) = -x*dydz;
00125 A(row, row*5+4) = -x;
00126 R(row, 0) = y-yp;
00127 } else {
00128 A(row, row*5+0) = -1.0;
00129 A(row, row*5+1) = dxdz;
00130 A(row, row*5+2) = y*dxdz;
00131 A(row, row*5+3) = -x*dxdz;
00132 A(row, row*5+4) = y;
00133 R(row, 0) = x-xp;
00134 }
00135 } else {
00136 if(sl[counter][counterHi]!=2) {
00137 A(row, row*3+0) = -1.0;
00138 A(row, row*3+1) = dxdz;
00139 A(row, row*3+2) = -x*dxdz;
00140 R(row, 0) = x-xp;
00141 }
00142 }
00143
00144 E(row, row) = 1.0/error;
00145 I(row, row) = 1.0;
00146
00147 if(sl[counter][counterHi] != 2) {
00148 B(row, 0) = 1.0;
00149 B(row, 1) = zp;
00150 } else {
00151 B(row, 2) = 1.0;
00152 B(row, 3) = zp;
00153 }
00154
00155 }
00156
00157 TMatrixD Bt = B; Bt.T();
00158 TMatrixD At = A; At.T();
00159 TMatrixD gamma = (Bt*E*B); gamma.Invert();
00160 *(C[wh[counter]+2][st[counter]-1][sr[counter]-1]) += At*E*(B*gamma*Bt*E*A-A);
00161 *(b[wh[counter]+2][st[counter]-1][sr[counter]-1]) += At*E*(B*gamma*Bt*E*I-I)*R;
00162
00163 }
00164
00165 }
00166
00167 }
00168
00169 if (workingmode==3)
00170 for(int wheel = -2; wheel < 3; ++wheel)
00171 for(int station = 1; station < 5; ++station)
00172 for(int sector = 1; sector < 15; ++sector) {
00173
00174 if(sector > 12 && station != 4) continue;
00175
00176 TMatrixD Ctr = *C[wheel+2][station-1][sector-1];
00177 TMatrixD btr = *b[wheel+2][station-1][sector-1];
00178
00179 TString CtrName = "Ctr_"; CtrName += wheel; CtrName += "_"; CtrName += station;
00180 CtrName += "_"; CtrName += sector;
00181 Ctr.Write(CtrName);
00182
00183 TString btrName = "btr_"; btrName += wheel; btrName += "_"; btrName += station;
00184 btrName += "_"; btrName += sector;
00185 btr.Write(btrName);
00186
00187 }
00188
00189
00190 if (workingmode%2==0) {
00191
00192 for(int wheel = -2; wheel < 3; ++wheel)
00193 for(int station = 1; station < 5; ++station)
00194 for(int sector = 1; sector < 15; ++sector) {
00195
00196 if(sector > 12 && station != 4) continue;
00197
00198 if (workingmode==4) {
00199
00200 for (int mf = -1; mf<=-1; mf++) {
00201
00202 TMatrixD Ch = getMatrixFromFile("Ctr_",wheel,station,sector,mf);
00203 TMatrixD bh = getMatrixFromFile("btr_",wheel,station,sector,mf);
00204
00205 *(C[wheel+2][station-1][sector-1]) += Ch;
00206 *(b[wheel+2][station-1][sector-1]) += bh;
00207
00208 }
00209
00210 }
00211
00212 TMatrixD Ctr = *C[wheel+2][station-1][sector-1];
00213 TMatrixD btr = *b[wheel+2][station-1][sector-1];
00214
00215 TMatrixD Ccs = getCcsMatrix(wheel, station, sector);
00216 TMatrixD bcs = getbcsMatrix(wheel, station, sector);
00217
00218 TMatrixD SC = Ctr + Ccs;
00219 TMatrixD Sb = btr + bcs;
00220
00221 SC.Invert();
00222
00223 TMatrixD solution = SC*Sb;
00224 for(int icounter = 0; icounter < SC.GetNrows(); icounter++) {
00225 for(int jcounter = 0; jcounter < SC.GetNrows(); jcounter++) {
00226 cov[icounter][jcounter] = SC(icounter, jcounter);
00227 }
00228 }
00229
00230 int nLayer = 12, nDeg = 5;
00231 if (station==4) { nLayer = 8; nDeg = 3; }
00232 for(int counterLayer = 0; counterLayer < nLayer; ++counterLayer)
00233 for(int counterDeg = 0; counterDeg < nDeg; ++counterDeg) {
00234
00235 if(counterLayer < 4) {
00236 slC[counterLayer] = 1;
00237 laC[counterLayer] = counterLayer+1;
00238 } else if(counterLayer > 3 && counterLayer < 8) {
00239 slC[counterLayer] = 3;
00240 laC[counterLayer] = counterLayer-3;
00241 } else {
00242 slC[counterLayer] = 2;
00243 laC[counterLayer] = counterLayer-7;
00244 }
00245
00246 if(counterLayer < 8) {
00247 dx[counterLayer] = solution(counterLayer*nDeg, 0);
00248 dy[counterLayer] = 0.0;
00249 } else {
00250 dx[counterLayer] = 0.0;
00251 dy[counterLayer] = solution(counterLayer*nDeg, 0);
00252 }
00253 dz[counterLayer] = solution(counterLayer*nDeg+1, 0);
00254
00255 if (station!=4) {
00256 phix[counterLayer] = solution(counterLayer*nDeg+2, 0);
00257 phiy[counterLayer] = solution(counterLayer*nDeg+3, 0);
00258 phiz[counterLayer] = solution(counterLayer*nDeg+4, 0);
00259 } else {
00260 phix[counterLayer] = 0.;
00261 phiy[counterLayer] = solution(counterLayer*nDeg+2, 0);
00262 phiz[counterLayer] = 0.;
00263 }
00264 }
00265
00266 whC = wheel; stC = station; srC = sector;
00267
00268 ttreeOutput->Fill();
00269
00270 }
00271
00272 }
00273
00274
00275
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325
00326
00327
00328
00329
00330
00331
00332
00333
00334
00335
00336
00337
00338
00339
00340 }
00341
00342
00343
00344 TMatrixD DTMuonMillepede::getCcsMatrix(int wh, int st, int se) {
00345
00346 int size = 60;
00347 if(st == 4) size = 24;
00348
00349 TMatrixD matrix(size, size);
00350
00351 float Error[3][6] = {{0.000001, 0.000001, 0.000001, 0.000001, 0.000001, 0.000001},
00352 {0.005, 0.005, 0.005, 0.00005, 0.00005, 0.00005},
00353 {0.005, 0.005, 0.005, 0.00005, 0.00005, 0.00005}};
00354
00355 float TollerancyPosition = 0.01;
00356 float TollerancyRotation = 0.0001;
00357
00358 TMatrixD ConsQC = myPG->giveQCCal(wh, st, se);
00359
00360 bool UseQC = true;
00361 if (ConsQC.GetNrows()<70) UseQC = false;
00362
00363 int nLayer = 12, nDeg = 5;
00364 if (st==4) { nLayer = 8; nDeg = 3; }
00365
00366 for (int la = 0; la<nLayer; la++)
00367 for (int dg = 0; dg<nDeg; dg++) {
00368
00369 int index = la*nDeg + dg;
00370
00371 int rdg = dg + 1;
00372 if (la<8 && rdg==1) rdg = 0;
00373 if (st==4 && rdg==3) rdg = 4;
00374
00375 double ThisError2 = Error[la/4][rdg]*Error[la/4][rdg];
00376 if (rdg<2) {
00377 if (UseQC) { ThisError2 += ConsQC(la, 1)*ConsQC(la, 1);
00378 } else { ThisError2 += TollerancyPosition*TollerancyPosition; }
00379 } else if (rdg==2) { ThisError2 += TollerancyPosition*TollerancyPosition;
00380 } else { ThisError2 += TollerancyRotation*TollerancyRotation; }
00381
00382 matrix(index, index) = 1./ThisError2;
00383
00384 }
00385
00386 return matrix;
00387
00388 }
00389
00390
00391
00392 TMatrixD DTMuonMillepede::getbcsMatrix(int wh, int st, int se) {
00393
00394 int size = 60;
00395 if(st == 4) size = 24;
00396
00397 TMatrixD matrix(size, 1);
00398
00399 float Error[3][6] = {{0.000001, 0.000001, 0.000001, 0.000001, 0.000001, 0.000001},
00400 {0.005, 0.005, 0.005, 0.00005, 0.00005, 0.00005},
00401 {0.005, 0.005, 0.005, 0.00005, 0.00005, 0.00005}};
00402
00403 float TollerancyPosition = 0.01;
00404 float TollerancyRotation = 0.0001;
00405
00406 TMatrixD ConsQC = myPG->giveQCCal(wh, st, se);
00407
00408 bool UseQC = true;
00409 if (ConsQC.GetNrows()<70) UseQC = false;
00410
00411 TMatrixD Survey = myPG->giveSurvey(wh, st, se);
00412
00413 int nLayer = 12, nDeg = 5;
00414 if (st==4) { nLayer = 8; nDeg = 3; }
00415
00416 for (int la = 0; la<nLayer; la++)
00417 for (int dg = 0; dg<nDeg; dg++) {
00418
00419 int index = la*nDeg + dg;
00420
00421 int rdg = dg + 1;
00422 if (la<8 && rdg==1) rdg = 0;
00423 if (st==4 && rdg==3) rdg = 4;
00424
00425 double ThisError2 = Error[la/4][rdg]*Error[la/4][rdg];
00426 if (rdg<2) {
00427 if (UseQC) { ThisError2 += ConsQC(la, 1)*ConsQC(la, 1);
00428 } else { ThisError2 += TollerancyPosition*TollerancyPosition; }
00429 } else if (rdg==2) { ThisError2 += TollerancyPosition*TollerancyPosition;
00430 } else { ThisError2 += TollerancyRotation*TollerancyRotation; }
00431
00432 float Constraint = 0.;
00433 if (la>3)
00434 Constraint += Survey(rdg, la/4-1);
00435 if (UseQC && rdg==0) Constraint += ConsQC(la, 0);;
00436 if (UseQC && rdg==1) Constraint -= ConsQC(la, 0);
00437
00438 matrix(index, 0) = Constraint/ThisError2;
00439
00440 }
00441
00442 return matrix;
00443
00444 }
00445
00446 TMatrixD DTMuonMillepede::getMatrixFromFile(TString Code, int wh, int st, int se, int mf) {
00447
00448 TString MtxFileName = "./LocalMillepedeMatrix_"; MtxFileName += mf; MtxFileName += ".root";
00449 if (mf==-1) MtxFileName = "./LocalMillepedeMatrix.root";
00450
00451 TDirectory *dirSave = gDirectory;
00452
00453 TFile *MatrixFile = new TFile(MtxFileName);
00454
00455 TString MtxName = Code; MtxName += wh; MtxName += "_"; MtxName += st; MtxName += "_"; MtxName += se;
00456 TMatrixD *ThisMtx = (TMatrixD *)MatrixFile->Get(MtxName);
00457
00458 MatrixFile->Close();
00459
00460 dirSave->cd();
00461
00462 return *ThisMtx;
00463
00464 }
00465
00466 TMatrixD DTMuonMillepede::getLagMatrix(int wh, int st, int se) {
00467
00468 TMatrixD matrix(60+6,60+6);
00469 if(st == 4) matrix.ResizeTo(40+6, 40+6);
00470
00471 for(int counterDeg = 0; counterDeg < 5; ++counterDeg) {
00472 for(int counterLayer = 0; counterLayer < 12; ++counterLayer) {
00473 if(st == 4) {
00474 matrix(40+counterDeg, counterLayer*5+counterDeg) = 10000.0;
00475 matrix(counterLayer*5+counterDeg, 40+counterDeg) = 10000.0;
00476 } else {
00477 int realCounter = counterDeg+1;
00478 if(counterDeg == 1 && counterLayer < 8) {
00479 realCounter = counterDeg-1;
00480 }
00481 matrix(counterLayer*5+counterDeg,40+realCounter) = 10000.0;
00482 if( (realCounter == 0 && counterLayer > 7) ||
00483 (realCounter == 1 && counterLayer < 7)) continue;
00484 matrix(60+realCounter,counterLayer*5+counterDeg) = 10000.0;
00485 }
00486 }
00487 }
00488 return matrix;
00489 }
00490
00491
00492 TMatrixD DTMuonMillepede::getCqcMatrix(int wh, int st, int se) {
00493
00494 TMatrixD surv = myPG->giveQCCal(wh, st, se);
00495 TMatrixD sigmaQC(5, 12);
00496
00497 TMatrixD matrix(60, 60);
00498 if(st == 4) matrix.ResizeTo(40, 40);
00499
00500
00501 if(surv.GetNrows() < 7) {
00502 for(int counterDeg = 0; counterDeg < 5; counterDeg++) {
00503 for(int counterLayer = 0; counterLayer < 12; ++counterLayer) {
00504 if(st != 4 && counterLayer > 7) continue;
00505 if(counterDeg == 0) {
00506 sigmaQC(counterDeg, counterLayer) = 0.05;
00507 } else if(counterDeg < 3) {
00508 sigmaQC(counterDeg, counterLayer) = 0.05;
00509 } else {
00510 sigmaQC(counterDeg, counterLayer) = 0.05;
00511 }
00512 }
00513 }
00514 } else {
00515 float meanvarSL1 = sqrt(surv(0,1)*surv(0,1)+surv(1,1)*surv(1,1)+surv(2,1)*surv(2,1)+surv(3,1)*surv(3,1))/10000.0;
00516 float meanvarSL2 = 0;
00517 if(surv.GetNrows() > 9) meanvarSL2 = sqrt(surv(8,1)*surv(8,1)+surv(9,1)*surv(9,1)+surv(10,1)*surv(10,1)+surv(11,1)*surv(11,1))/10000.0;
00518 float meanvarSL3 = sqrt(surv(4,1)*surv(4,1)+surv(5,1)*surv(5,1)+surv(6,1)*surv(6,1)+surv(7,1)*surv(7,1))/10000.0;
00519 for(int counterDeg = 0; counterDeg < 5; counterDeg++) {
00520 for(int counterLayer = 0; counterLayer < 12; ++counterLayer) {
00521 if(st != 4 && counterLayer > 7) continue;
00522 float meanerror = 0;
00523 if(counterLayer < 4) {
00524 meanerror = meanvarSL1;
00525 } else if(counterLayer > 3 && counterLayer < 8){
00526 meanerror = meanvarSL2;
00527 } else {
00528 meanerror = meanvarSL3;
00529 }
00530 if(counterDeg == 0) {
00531 sigmaQC(counterDeg, counterLayer) = sqrt((surv(counterLayer, 1)/10000.0)*(surv(counterLayer, 1)/10000.0)+meanerror*meanerror);
00532 } else if(counterDeg < 3) {
00533 sigmaQC(counterDeg, counterLayer) = 0.05;
00534 } else {
00535 sigmaQC(counterDeg, counterLayer) = 0.05;
00536 }
00537 }
00538 }
00539 double **Eta = new double *[12];
00540 for(int counterLayer = 0; counterLayer < 12; counterLayer++) {
00541
00542 if(counterLayer > 7 && st == 4) continue;
00543 Eta[counterLayer] = new double [5];
00544 for(int counterLayer2 = 0; counterLayer2 < 12; counterLayer2++) {
00545 if(counterLayer > 7 && st == 4) continue;
00546 if((counterLayer2 < 4 && counterLayer < 4) || (counterLayer2 > 3 && counterLayer > 3)) {
00547 if(counterLayer == counterLayer2) {
00548 Eta[counterLayer][counterLayer2] = 3.0/(4.0);
00549 } else {
00550 Eta[counterLayer][counterLayer2] = -1.0/(4.0);
00551 }
00552 } else {
00553 Eta[counterLayer][counterLayer2] = 0.0;
00554 }
00555 }
00556 }
00557
00558 for(int counterDeg = 0; counterDeg < 5; counterDeg++) {
00559 for(int counterLayer = 0; counterLayer < 12; counterLayer++) {
00560 if(counterLayer > 7 && st == 4) continue;
00561 for(int counterLayer2 = 0; counterLayer2 < 12; counterLayer2++) {
00562 if(counterLayer2 > 7 && st == 4) continue;
00563 for(int counterLayer3 = 0; counterLayer3 < 12; counterLayer3++) {
00564 if(counterLayer3 > 7 && st == 4) continue;
00565 matrix(5*counterLayer2+counterDeg, 5*counterLayer3+counterDeg) +=
00566 Eta[counterLayer][counterLayer2]*Eta[counterLayer][counterLayer3]/(sigmaQC(counterDeg,counterLayer)*sigmaQC(counterDeg,counterLayer));
00567 }
00568 }
00569 }
00570 }
00571 }
00572 return matrix;
00573 }
00574
00575
00576 TMatrixD DTMuonMillepede::getbqcMatrix(int wh, int st, int se) {
00577
00578 TMatrixD surv = myPG->giveQCCal(wh, st, se);
00579 TMatrixD ResQC(5, 12);
00580 TMatrixD sigmaQC(5, 12);
00581 TMatrixD matrix(60, 1);
00582 if(st == 4) matrix.ResizeTo(40, 1);
00583
00584 if(surv.GetNrows() < 7) {
00585 for(int counterDeg = 0; counterDeg < 5; counterDeg++) {
00586 for(int counterLayer = 0; counterLayer < 12; ++counterLayer) {
00587 if(st != 4 && counterLayer > 7) continue;
00588 if(counterDeg == 0) {
00589 sigmaQC(counterDeg, counterLayer) = 0.05;
00590 } else if(counterDeg < 3) {
00591 sigmaQC(counterDeg, counterLayer) = 0.05;
00592 } else {
00593 sigmaQC(counterDeg, counterLayer) = 0.05;
00594 }
00595 }
00596 }
00597 } else {
00598 for(int counterChamber = 0; counterChamber < 12; ++counterChamber) {
00599 ResQC(0, counterChamber) = -surv(counterChamber, 0)/10000.0;
00600 }
00601 float meanvarSL1 = sqrt(surv(0,1)*surv(0,1)+surv(1,1)*surv(1,1)+surv(2,1)*surv(2,1)+surv(3,1)*surv(3,1))/10000.0;
00602 float meanvarSL2 = 0;
00603 if(surv.GetNrows() > 9) {
00604 meanvarSL2 = sqrt(surv(8,1)*surv(8,1)+surv(9,1)*surv(9,1)+surv(10,1)*surv(10,1)+surv(11,1)*surv(11,1))/10000.0;
00605 }
00606 float meanvarSL3 = sqrt(surv(4,1)*surv(4,1)+surv(5,1)*surv(5,1)+surv(6,1)*surv(6,1)+surv(7,1)*surv(7,1))/10000.0;
00607 for(int counterDeg = 0; counterDeg < 5; counterDeg++) {
00608 for(int counterLayer = 0; counterLayer < 12; ++counterLayer) {
00609 if(st != 4 && counterLayer > 7) continue;
00610 float meanerror = 0;
00611 if(counterLayer < 4) {
00612 meanerror = meanvarSL1;
00613 } else if(counterLayer > 3 && counterLayer < 8){
00614 meanerror = meanvarSL2;
00615 } else {
00616 meanerror = meanvarSL3;
00617 }
00618 if(counterDeg == 0) {
00619 sigmaQC(counterDeg, counterLayer) = sqrt((surv(counterLayer, 1)/10000.0)*(surv(counterLayer, 1)/10000.0)+meanerror*meanerror);
00620 } else if(counterDeg < 3) {
00621 sigmaQC(counterDeg, counterLayer) = 0.05;
00622 } else {
00623 sigmaQC(counterDeg, counterLayer) = 0.05;
00624 }
00625 }
00626 }
00627 double **Eta = new double *[12];
00628 for(int counterLayer = 0; counterLayer < 12; counterLayer++) {
00629 if(counterLayer > 7 && st == 4) continue;
00630 Eta[counterLayer] = new double [5];
00631 for(int counterLayer2 = 0; counterLayer2 < 12; counterLayer2++) {
00632 if(counterLayer > 7 && st == 4) continue;
00633 if((counterLayer2 < 4 && counterLayer < 4) || (counterLayer2 > 3 && counterLayer > 3)) {
00634 if(counterLayer == counterLayer2) {
00635 Eta[counterLayer][counterLayer2] = 3.0/(4.0);
00636 } else {
00637 Eta[counterLayer][counterLayer2] = -1.0/(4.0);
00638 }
00639 } else {
00640 Eta[counterLayer][counterLayer2] = 0.0;
00641 }
00642 }
00643 }
00644
00645 for(int counterDeg = 0; counterDeg < 5; counterDeg++) {
00646 for(int counterLayer = 0; counterLayer < 12; counterLayer++) {
00647 if(counterLayer > 7 && st == 4) continue;
00648 for(int counterLayer2 = 0; counterLayer2 < 12; counterLayer2++) {
00649 if(counterLayer2 > 7 && st == 4) continue;
00650 float mean = 0;
00651 if(counterDeg != 0) {
00652 if(counterLayer < 4) {
00653 mean = (ResQC(counterDeg, 0)+ResQC(counterDeg, 1)+ResQC(counterDeg, 2)+ResQC(counterDeg, 3))/4.0;
00654 } else if(counterLayer > 3 && counterLayer < 8) {
00655 mean = (ResQC(counterDeg, 4)+ResQC(counterDeg, 5)+ResQC(counterDeg, 6)+ResQC(counterDeg, 7))/4.0;
00656 } else {
00657 mean = (ResQC(counterDeg, 8)+ResQC(counterDeg, 9)+ResQC(counterDeg, 10)+ResQC(counterDeg, 11))/4.0;
00658 }
00659 }
00660 matrix(5*counterLayer2+counterDeg, 0) += Eta[counterLayer][counterLayer2]*(ResQC(counterDeg,counterLayer)-mean)/(sigmaQC(counterDeg,counterLayer)*sigmaQC(counterDeg,counterLayer));
00661 }
00662 }
00663 }
00664 }
00665 return matrix;
00666 }
00667
00668
00669
00670
00671
00672
00673
00674
00675
00676
00677
00678
00679
00680
00681
00682
00683
00684
00685
00686
00687
00688
00689
00690
00691
00692
00693
00694 TMatrixD DTMuonMillepede::getCsurveyMatrix(int wh, int st, int se) {
00695
00696 int size = 60;
00697 if(st == 4) size = 40;
00698
00699 TMatrixD matrix(size+6, size+6);
00700
00701 float error[2][6] = {{0.000001, 0.000001, 0.000001, 0.000001, 0.000001, 0.000001},
00702 {0.05, 0.05, 0.05, 0.005, 0.005, 0.005}};
00703 for(int counterLayer = 0; counterLayer < size/5; counterLayer++) {
00704 for(int counterLayer2 = 0; counterLayer2 < size/5; counterLayer2++) {
00705 int sl1 = counterLayer/4;
00706 int sl2 = counterLayer2/4;
00707 if (sl1==sl2) {
00708 for(int counterDeg = 0; counterDeg < 5; counterDeg++) {
00709 int counterDegAux = counterDeg+1;
00710 if(counterLayer < 8 && counterDeg == 1) counterDegAux = 0;
00711 int sl = (sl1 + 1)/2;
00712 matrix(5*counterLayer+counterDeg, 5*counterLayer2+counterDeg) = 1.0/(16.0*error[sl][counterDegAux]*error[sl][counterDegAux]);
00713 }
00714 }
00715 }
00716 }
00717
00718 return matrix;
00719
00720 }
00721
00722
00723
00724
00725
00726
00727
00728
00729
00730
00731
00732
00733
00734
00735
00736
00737
00738
00739
00740
00741
00742
00743
00744
00745
00746
00747
00748
00749 TMatrixD DTMuonMillepede::getbsurveyMatrix(int wh, int st, int se) {
00750
00751 TMatrixD survey = myPG->giveSurvey(wh, st, se);
00752 int size = 60;
00753 if(st == 4) size = 40;
00754
00755 TMatrixD matrix(size+6, 1);
00756
00757 float error[2][6] = {{0.000001, 0.000001, 0.000001, 0.000001, 0.000001, 0.000001},
00758 {0.05, 0.05, 0.05, 0.005, 0.005, 0.005}};
00759 for(int counterLayer = 0; counterLayer < size/5; counterLayer++) {
00760 for(int counterDeg = 0; counterDeg < 5; counterDeg++) {
00761 int counterDegAux = counterDeg+1;
00762 if(counterLayer < 8 && counterDeg == 1) counterDegAux = 0;
00763 int superlayerAux = counterLayer/4;
00764 int sl = (superlayerAux + 1)/2;
00765 matrix(5*counterLayer+counterDeg, 0) = survey(counterDegAux, superlayerAux)/(4.0*error[sl][counterDegAux]*error[sl][counterDegAux]);
00766 }
00767 }
00768
00769 return matrix;
00770
00771 }
00772
00773
00774 TMatrixD DTMuonMillepede::prepareForLagrange(TMatrixD m) {
00775
00776 TMatrixD updatedMatrix = m;
00777 updatedMatrix.ResizeTo(m.GetNrows()+6, m.GetNcols()+6);
00778 return updatedMatrix;
00779
00780 }
00781
00782
00783 void DTMuonMillepede::setBranchTree() {
00784
00785 ttreeOutput = new TTree("DTLocalMillepede", "DTLocalMillepede");
00786
00787 ttreeOutput->Branch("wh", &whC, "wh/I");
00788 ttreeOutput->Branch("st", &stC, "st/I");
00789 ttreeOutput->Branch("sr", &srC, "sr/I");
00790 ttreeOutput->Branch("cov", cov, "cov[60][60]/F");
00791 ttreeOutput->Branch("sl", slC, "sl[12]/I");
00792 ttreeOutput->Branch("la", laC, "la[12]/I");
00793 ttreeOutput->Branch("dx", dx, "dx[12]/F");
00794 ttreeOutput->Branch("dy", dy, "dy[12]/F");
00795 ttreeOutput->Branch("dz", dz, "dz[12]/F");
00796 ttreeOutput->Branch("phix", phix, "phix[12]/F");
00797 ttreeOutput->Branch("phiy", phiy, "phiy[12]/F");
00798 ttreeOutput->Branch("phiz", phiz, "phiz[12]/F");
00799
00800 }
00801
00802
00803
00804