File indexing completed on 2024-04-06 11:56:32
0001
0002
0003
0004
0005 #include "PlotMillePede.h"
0006 #include <TH1.h>
0007 #include <TH2.h>
0008 #include <TProfile.h>
0009 #include <TArrow.h>
0010 #include <TEllipse.h>
0011 #include <TF1.h>
0012 #include <TMath.h>
0013 #include <TTree.h>
0014 #include <TPaveText.h>
0015 #include <TObjString.h>
0016 #include <TError.h>
0017 #include <TROOT.h>
0018 #include <TCanvas.h>
0019
0020 #include <iostream>
0021
0022 #include "GFUtils/GFHistManager.h"
0023
0024
0025 PlotMillePede::PlotMillePede(const char *fileName, Int_t iov, Int_t hieraLevel, bool useDiff)
0026 : MillePedeTrees(fileName, iov),
0027 fHistManager(new GFHistManager),
0028 fHieraLevel(hieraLevel),
0029 fUseDiff(useDiff),
0030 fSubDetIds(),
0031 fAlignableTypeId(-1),
0032 fMaxDevUp(500.),
0033 fMaxDevDown(-500.),
0034 fNbins(101) {
0035 fHistManager->SetLegendX1Y1X2Y2(0.14, 0.7, 0.45, 0.9);
0036 }
0037
0038 PlotMillePede::PlotMillePede(const char *fileName, Int_t iov, Int_t hieraLevel, const char *treeNameAdd)
0039 : MillePedeTrees(fileName, iov, treeNameAdd),
0040 fHistManager(new GFHistManager),
0041 fHieraLevel(hieraLevel),
0042 fUseDiff(false),
0043 fSubDetIds(),
0044 fAlignableTypeId(-1),
0045 fMaxDevUp(500.),
0046 fMaxDevDown(-500.),
0047 fNbins(101) {
0048 fHistManager->SetLegendX1Y1X2Y2(0.14, 0.7, 0.45, 0.9);
0049 }
0050
0051
0052 PlotMillePede::~PlotMillePede() { delete fHistManager; }
0053
0054
0055 void PlotMillePede::DrawAll(Option_t *opt) {
0056 const TString o(opt);
0057 bool wasBatch = fHistManager->SetBatch();
0058 fHistManager->Clear();
0059
0060
0061 if (o.Contains("r", TString::kIgnoreCase))
0062 this->DrawParamResult("add");
0063 if (o.Contains("o", TString::kIgnoreCase))
0064 this->DrawOrigParam(true);
0065 if (o.Contains("g", TString::kIgnoreCase))
0066 this->DrawGlobCorr(true);
0067 if (o.Contains("p", TString::kIgnoreCase))
0068 this->DrawPull("add");
0069 if (o.Contains("m", TString::kIgnoreCase))
0070 this->DrawMisVsLocation(true);
0071 if (o.Contains("e", TString::kIgnoreCase))
0072 this->DrawErrorVsHit(true);
0073 if (o.Contains("h", TString::kIgnoreCase))
0074 this->DrawHitMaps(true);
0075
0076 fHistManager->SetBatch(wasBatch);
0077 if (!wasBatch)
0078 fHistManager->Draw();
0079 }
0080
0081
0082 void PlotMillePede::DrawParam(bool addPlots, const TString &sel) {
0083 const Int_t layer = this->PrepareAdd(addPlots);
0084
0085 const TString titleAdd = this->TitleAdd();
0086 UInt_t nPlot = 0;
0087 for (UInt_t iPar = 0; iPar < kNpar; ++iPar) {
0088 TString parSel(sel.Length() ? sel : Fixed(iPar, false));
0089 this->AddBasicSelection(parSel);
0090 const TString toMum(this->ToMumMuRad(iPar));
0091
0092
0093 const TString hNameA(this->Unique(Form("after%d", iPar)));
0094 const TString hNameB(this->Unique(Form("before%d", iPar)));
0095 TH1 *hAfter = this->CreateHist(FinalMisAlignment(iPar) += toMum, parSel, hNameA);
0096 TH1 *hBefore = this->CreateHist(Parenth(MisParT() += Par(iPar)) += toMum, parSel, hNameB);
0097
0098 if (0. == hAfter->GetEntries())
0099 continue;
0100 hBefore->SetTitle(DelName(iPar) += titleAdd + ";" + DelNameU(iPar) += ";#parameters");
0101 hAfter->SetTitle(DelName(iPar) += titleAdd + ";" + DelNameU(iPar) += ";#parameters");
0102 fHistManager->AddHist(hBefore, layer, "before");
0103 fHistManager->AddHistSame(hAfter, layer, nPlot, "after");
0104 fHistManager->AddHist(static_cast<TH1 *>(hAfter->Clone()), layer + 1, "after");
0105 fHistManager->AddHistSame(static_cast<TH1 *>(hBefore->Clone()), layer + 1, nPlot, "before");
0106 ++nPlot;
0107 }
0108
0109 fHistManager->Draw();
0110 }
0111
0112
0113 void PlotMillePede::DrawPedeParam(Option_t *option, unsigned int nNonRigidParam) {
0114 const bool addParVsPar = TString(option).Contains("vs", TString::kIgnoreCase);
0115 const Int_t layer = this->PrepareAdd(TString(option).Contains("add", TString::kIgnoreCase));
0116 const unsigned int nPar = kNpar + nNonRigidParam;
0117 const TString titleAdd = this->TitleAdd();
0118 UInt_t nPlot = 0;
0119 TObjArray primitivesParVsPar;
0120 for (UInt_t iPar = 0; iPar < nPar; ++iPar) {
0121 TString parSel(Fixed(iPar, false) += AndL() += Valid(iPar));
0122 this->AddBasicSelection(parSel);
0123 const TString toMum(this->ToMumMuRadPede(iPar));
0124 const TString hName(this->Unique(Form("pedePar%d", iPar)));
0125 TH1 *h = this->CreateHist(Parenth(MpT() += Par(iPar)) += toMum, parSel, hName);
0126
0127 if (0. == h->GetEntries())
0128 continue;
0129
0130 TObjArray histsParVsPar;
0131 if (addParVsPar) {
0132 for (UInt_t iPar2 = iPar + 1; iPar2 < nPar; ++iPar2) {
0133 const TString hNameVs(this->Unique(Form("pedePar%d_%d", iPar, iPar2)));
0134 const TString toMum2(this->ToMumMuRadPede(iPar2));
0135 TH2 *hVs = this->CreateHist2D(Parenth(MpT() += Par(iPar)) += toMum,
0136 Parenth(MpT() += Par(iPar2)) += toMum2,
0137 Parenth(parSel) += AndL() += Fixed(iPar2, false),
0138 hNameVs,
0139 "BOX");
0140 if (0. == hVs->GetEntries())
0141 continue;
0142 else {
0143 hVs->SetTitle("pede: " + NamePede(iPar2) += " vs. " + NamePede(iPar) += titleAdd + ";" + NamePede(iPar) +=
0144 UnitPede(iPar) += ";" + NamePede(iPar2) += UnitPede(iPar2));
0145 histsParVsPar.Add(hVs);
0146
0147 TPaveText *pave = new TPaveText(.15, .15, .5, .25, "NDC");
0148 pave->SetBorderSize(1);
0149 pave->AddText(Form("#rho = %.3f", hVs->GetCorrelationFactor()));
0150 primitivesParVsPar.Add(pave);
0151 }
0152 }
0153 }
0154
0155
0156 const TString hNameBySi(this->Unique(Form("pedeParBySi%d", iPar)));
0157 TH1 *hBySi = this->CreateHist(
0158 Parenth(MpT() += Par(iPar)) += Div() += ParSi(iPar), parSel + AndL() += ParSiOk(iPar), hNameBySi);
0159 TH1 *hBySiInv = 0;
0160 TH2 *hSiVsPar = 0;
0161 if (hBySi->GetEntries() == 0.) {
0162 delete hBySi;
0163 hBySi = 0;
0164 } else {
0165 const TString hNameBySiInv(this->Unique(Form("pedeParBySiInv%d", iPar)) += "(100,-20,20)");
0166 hBySiInv = this->CreateHist(
0167 ParSi(iPar) += Div() += Parenth(MpT() += Par(iPar)), parSel + AndL() += ParSiOk(iPar), hNameBySiInv);
0168 const TString hNameSiVsPar(this->Unique(Form("pedeParVsSi%d", iPar)));
0169 hSiVsPar = this->CreateHist2D(Parenth(MpT() += Par(iPar)) += toMum,
0170 ParSi(iPar) += toMum,
0171 parSel + AndL() += ParSiOk(iPar),
0172 hNameSiVsPar,
0173 "BOX");
0174 }
0175
0176
0177 const TString hNameH(this->Unique(Form("pedeParVsHits%d", iPar)));
0178 TH2 *hHits = this->CreateHist2D(HitsX(), Parenth(MpT() += Par(iPar)) += toMum, parSel, hNameH, "BOX");
0179
0180
0181 const TString hNameG(this->Unique(Form("pedeParVsGlob%d", iPar)));
0182 TH2 *hGlobCor = this->CreateHist2D(
0183 Cor(iPar), Parenth(MpT() += Par(iPar)) += toMum, parSel + AndL() += Cor(iPar) += ">-0.1999", hNameG, "BOX");
0184 if (!hGlobCor->GetEntries()) {
0185 delete hGlobCor;
0186 hGlobCor = 0;
0187 }
0188
0189 h->SetTitle("determined pede " + NamePede(iPar) += titleAdd + ";" + NamePede(iPar) += UnitPede(iPar) +=
0190 ";#alignables");
0191 hHits->SetTitle("determined pede " + NamePede(iPar) += titleAdd + " vs #n(hit_{x});N_{hit,x};" + NamePede(iPar) +=
0192 UnitPede(iPar));
0193 if (hGlobCor)
0194 hGlobCor->SetTitle("determined pede " + NamePede(iPar) +=
0195 titleAdd + " vs glob. corr;Global Correlation;" + NamePede(iPar) += UnitPede(iPar));
0196 fHistManager->AddHist(h, layer);
0197 fHistManager->AddHist(hHits, layer + 1);
0198 if (addParVsPar)
0199 fHistManager->AddHists(&histsParVsPar, layer + 2);
0200 if (hGlobCor)
0201 fHistManager->AddHist(hGlobCor, layer + (addParVsPar ? 3 : 2));
0202
0203 if (hBySi) {
0204 const TString namI(NamePede(iPar));
0205 hBySi->SetTitle("pede: " + namI + "/#sigma_{" + namI + "}" + titleAdd + ";" + namI + Div() +=
0206 Fun("#sigma", namI) += ";#alignables");
0207 fHistManager->AddHist(hBySi, layer + 2 + (hGlobCor ? 1 : 0) + addParVsPar);
0208 hBySiInv->SetTitle("pede: #sigma_{" + namI + "}/" + namI + titleAdd + ";" + Fun("#sigma", namI) += Div() +=
0209 namI + ";#alignables");
0210 fHistManager->AddHist(hBySiInv, layer + 3 + (hGlobCor ? 1 : 0) + addParVsPar);
0211 hSiVsPar->SetTitle("pede: #sigma_{" + namI + " } vs " + namI + titleAdd + ";" + namI + UnitPede(iPar) + ";" +
0212 Fun("#sigma", namI) += UnitPede(iPar));
0213 fHistManager->AddHist(hSiVsPar, layer + 4 + (hGlobCor ? 1 : 0) + addParVsPar);
0214 }
0215 ++nPlot;
0216 }
0217
0218 if (addParVsPar) {
0219 for (Int_t i = 0; i < primitivesParVsPar.GetEntries(); ++i) {
0220 fHistManager->AddObject(primitivesParVsPar[i], layer + 2, i);
0221 }
0222 }
0223
0224 fHistManager->Draw();
0225 }
0226
0227
0228 void PlotMillePede::DrawPedeParamVsLocation(Option_t *option, unsigned int nNonRigidParam) {
0229 const Int_t layer = this->PrepareAdd(TString(option).Contains("add", TString::kIgnoreCase));
0230 const TString titleAdd = this->TitleAdd();
0231
0232 const unsigned int nPar = kNpar + nNonRigidParam;
0233 UInt_t nPlot = 0;
0234 for (UInt_t iPar = 0; iPar < nPar; ++iPar) {
0235 TString parSel(Fixed(iPar, false) += AndL() += Valid(iPar));
0236 this->AddBasicSelection(parSel);
0237 const TString toMum(this->ToMumMuRadPede(iPar));
0238 const TString pedePar(Parenth(MpT() += Par(iPar)) += toMum);
0239
0240 const TString nDpz(this->Unique(Form("pedeParZ%d", iPar)));
0241 TH2 *hPedeParZ = this->CreateHist2D(OrgPosT() += ZPos(), pedePar, parSel, nDpz, "BOX");
0242
0243 if (0. == hPedeParZ->GetEntries())
0244 continue;
0245
0246 const TString nDpr(this->Unique(Form("pedeParR%d", iPar)));
0247 TH2 *hPedeParR = this->CreateHist2D(RPos(OrgPosT()), pedePar, parSel, nDpr, "BOX");
0248
0249 const TString nDpp(this->Unique(Form("pedeParPhi%d", iPar)));
0250 TH2 *hPedeParPhi = this->CreateHist2D(Phi(OrgPosT()), pedePar, parSel, nDpp, "BOX");
0251
0252
0253
0254
0255
0256 const TString nDpy(this->Unique(Form("pedeParY%d", iPar)));
0257 TH2 *hPedeParY = this->CreateHist2D(OrgPosT() += YPos(), pedePar, parSel, nDpy, "BOX");
0258
0259 const TString title("determined pede " + NamePede(iPar) += " vs. %s" + titleAdd + ";%s;" + NamePede(iPar) +=
0260 UnitPede(iPar));
0261 hPedeParZ->SetTitle(Form(title.Data(), "z", "z [cm]"));
0262 hPedeParR->SetTitle(Form(title.Data(), "r", "r [cm]"));
0263 hPedeParPhi->SetTitle(Form(title.Data(), "#phi", "#phi"));
0264
0265 hPedeParY->SetTitle(Form(title.Data(), "y", "y [cm]"));
0266
0267 fHistManager->AddHist(hPedeParR, layer + nPlot);
0268 fHistManager->AddHist(hPedeParZ, layer + nPlot);
0269 fHistManager->AddHist(hPedeParPhi, layer + nPlot);
0270
0271 fHistManager->AddHist(hPedeParY, layer + nPlot);
0272
0273
0274 ++nPlot;
0275 }
0276
0277 fHistManager->Draw();
0278 }
0279
0280
0281
0282
0283
0284
0285
0286
0287
0288
0289
0290
0291
0292
0293
0294
0295
0296
0297
0298
0299
0300
0301
0302
0303
0304
0305
0306
0307
0308
0309
0310
0311
0312
0313
0314
0315
0316
0317
0318
0319
0320
0321
0322
0323
0324
0325
0326
0327
0328
0329
0330
0331
0332
0333
0334
0335
0336
0337
0338
0339
0340
0341
0342
0343
0344
0345
0346
0347
0348
0349
0350
0351
0352
0353
0354
0355
0356
0357
0358
0359
0360
0361
0362
0363
0364
0365
0366
0367
0368
0369
0370
0371
0372
0373
0374
0375
0376
0377
0378
0379
0380
0381
0382
0383
0384
0385
0386
0387
0388
0389
0390
0391
0392
0393
0394
0395
0396
0397
0398
0399
0400
0401
0402
0403
0404
0405
0406
0407
0408
0409
0410
0411
0412
0413
0414 void PlotMillePede::DrawSurfaceDeformations(const TString &whichOne,
0415 Option_t *option,
0416 unsigned int maxNumPars,
0417 unsigned int firstPar) {
0418 const Int_t layer = this->PrepareAdd(TString(option).Contains("add", TString::kIgnoreCase));
0419 const TString titleAdd = this->TitleAdd();
0420 const bool noLimit = TString(option).Contains("nolimit", TString::kIgnoreCase);
0421
0422 TString parSel(Valid(0) += AndL() += Fixed(0, false));
0423 if (TString(option).Contains("all", TString::kIgnoreCase))
0424 parSel = "";
0425 this->AddBasicSelection(parSel);
0426
0427 TObjArray whichOnes;
0428 if (whichOne.Contains("result", TString::kIgnoreCase))
0429 whichOnes.Add(new TObjString("result"));
0430 if (whichOne.Contains("start", TString::kIgnoreCase))
0431 whichOnes.Add(new TObjString("start"));
0432 if (whichOne.Contains("diff", TString::kIgnoreCase))
0433 whichOnes.Add(new TObjString("diff"));
0434
0435 for (Int_t wi = 0; wi < whichOnes.GetEntriesFast(); ++wi) {
0436 unsigned int nPlot = 0;
0437 for (unsigned int i = firstPar; i < maxNumPars; ++i) {
0438 TString hName(this->Unique(Form("hSurf%s%u", whichOnes[wi]->GetName(), i)));
0439 if (!noLimit)
0440 hName += Form("(%d,%f,%f)", fNbins, fMaxDevDown, fMaxDevUp);
0441 TH1 *h = this->CreateHist(DeformValue(i, whichOnes[wi]->GetName()) += this->ToMumMuRadSurfDef(i),
0442 parSel + AndL() += Parenth(NumDeformValues(whichOnes[wi]->GetName())),
0443 hName);
0444
0445 if (!h || 0. == h->GetEntries())
0446 continue;
0447 h->SetTitle(Form("SurfaceDeformation %s ", whichOnes[wi]->GetName()) + NameSurfDef(i) +=
0448 titleAdd + ";" + NameSurfDef(i) += UnitSurfDef(i));
0449 fHistManager->AddHistSame(h, layer, nPlot, whichOnes[wi]->GetName());
0450
0451 ++nPlot;
0452 }
0453 }
0454
0455 whichOnes.Delete();
0456 const bool old = fHistManager->SameWithStats(true);
0457 fHistManager->Draw();
0458 fHistManager->SameWithStats(old);
0459 }
0460
0461
0462 void PlotMillePede::DrawSurfaceDeformationsVsLocation(const TString &whichOne,
0463 Option_t *option,
0464 unsigned int maxNumPar,
0465 unsigned int firstPar) {
0466 const Int_t layer = this->PrepareAdd(TString(option).Contains("add", TString::kIgnoreCase));
0467 const TString titleAdd = this->TitleAdd();
0468
0469
0470
0471
0472
0473
0474 UInt_t nPlot = 0;
0475 for (UInt_t iPar = firstPar; iPar <= maxNumPar; ++iPar) {
0476 TString parSel(Valid(0) += AndL() += Fixed(0, false));
0477 if (TString(option).Contains("all", TString::kIgnoreCase))
0478 parSel = "";
0479 this->AddBasicSelection(parSel);
0480
0481
0482 TString hNameR(this->Unique(Form("hSurfR%u", iPar)));
0483
0484 TH1 *hR = this->CreateHist2D(RPos(OrgPosT()),
0485 DeformValue(iPar, whichOne) += this->ToMumMuRadSurfDef(iPar),
0486 parSel + AndL() += Parenth(NumDeformValues(whichOne)),
0487 hNameR,
0488 "BOX");
0489
0490 if (!hR || 0. == hR->GetEntries())
0491 continue;
0492
0493 TString hNameZ(this->Unique(Form("hSurfZ%u", iPar)));
0494
0495 TH1 *hZ = this->CreateHist2D(OrgPosT() += ZPos(),
0496 DeformValue(iPar, whichOne) += this->ToMumMuRadSurfDef(iPar),
0497 parSel + AndL() += Parenth(NumDeformValues(whichOne)),
0498 hNameZ,
0499 "BOX");
0500
0501 const TString title("Surface deformation " + NameSurfDef(iPar) +=
0502 " vs. %s" + titleAdd + ";%s;" + NameSurfDef(iPar) += UnitSurfDef(iPar));
0503 hR->SetTitle(Form(title.Data(), "r", "r [cm]"));
0504 hZ->SetTitle(Form(title.Data(), "z", "z [cm]"));
0505
0506 fHistManager->AddHist(hR, layer + nPlot);
0507 fHistManager->AddHist(hZ, layer + nPlot);
0508 ++nPlot;
0509 }
0510
0511 fHistManager->Draw();
0512 }
0513
0514
0515 void PlotMillePede::DrawSurfaceDeformationsLayer(Option_t *option,
0516 const unsigned int firstDetLayer,
0517 const unsigned int lastDetLayer,
0518 const TString &whichOne,
0519 unsigned int maxNumPars) {
0520 const TString opt(option);
0521 const Int_t firstLayer = this->PrepareAdd(opt.Contains("add", TString::kIgnoreCase));
0522 const bool noLimit = opt.Contains("nolimit", TString::kIgnoreCase);
0523 const bool spread = opt.Contains("spread", TString::kIgnoreCase);
0524 const bool verbose = opt.Contains("verbose", TString::kIgnoreCase);
0525 const bool verbose2 = opt.Contains("verboseByParam", TString::kIgnoreCase);
0526
0527 this->SetDetLayerCuts(0, false);
0528
0529 unsigned int iParUsed = 0;
0530 for (unsigned int iPar = 0; iPar < maxNumPars; ++iPar) {
0531
0532 const unsigned int numDetLayers = lastDetLayer - firstDetLayer + 1;
0533 TH1 *layerHist = new TH1F(
0534 this->Unique("hSurfAll" + whichOne += iPar),
0535 "Average deformations " + NameSurfDef(iPar) += ";;#LT" + (NameSurfDef(iPar) += "#GT") += UnitSurfDef(iPar),
0536 numDetLayers,
0537 0,
0538 numDetLayers);
0539 TH1 *layerHistWithSpread =
0540 (spread ? static_cast<TH1 *>(layerHist->Clone(Form("%s_spread", layerHist->GetName()))) : 0);
0541 if (spread)
0542 layerHistWithSpread->SetTitle(Form("%s (err is spread)", layerHist->GetTitle()));
0543
0544 TH1 *layerHistRms = new TH1F(
0545 this->Unique("hSurfAllRms" + whichOne += iPar),
0546 "RMS of deformations " + NameSurfDef(iPar) += ";;RMS(" + (NameSurfDef(iPar) += ")") += UnitSurfDef(iPar),
0547 numDetLayers,
0548 0,
0549 numDetLayers);
0550
0551
0552 unsigned int iDetLayerUsed = 0;
0553 for (unsigned int iDetLayer = firstDetLayer; iDetLayer <= lastDetLayer; ++iDetLayer) {
0554 if (!this->SetDetLayerCuts(iDetLayer, true))
0555 continue;
0556 TString sel;
0557 this->AddBasicSelection(sel);
0558
0559 const TString hName(this->Unique(Form("hSurf%s%u_%u", whichOne.Data(), iPar, iDetLayer)) +=
0560 (noLimit ? "" : Form("(%d,%f,%f)", fNbins, fMaxDevDown, fMaxDevUp)));
0561
0562 TH1 *h = this->CreateHist(DeformValue(iPar, whichOne) += this->ToMumMuRadSurfDef(iPar),
0563 (sel + AndL()) += Parenth(DeformValue(iPar, whichOne) += "!= 0."),
0564 hName);
0565
0566 if (!h || 0. == h->GetEntries())
0567 continue;
0568 ++iDetLayerUsed;
0569 layerHist->GetXaxis()->SetBinLabel(iDetLayerUsed, this->DetLayerLabel(iDetLayer));
0570 layerHist->SetBinContent(iDetLayerUsed, h->GetMean());
0571 layerHist->SetBinError(iDetLayerUsed, h->GetMeanError());
0572 layerHistRms->GetXaxis()->SetBinLabel(iDetLayerUsed, this->DetLayerLabel(iDetLayer));
0573 layerHistRms->SetBinContent(iDetLayerUsed, h->GetRMS());
0574 if (spread)
0575 layerHistWithSpread->SetBinContent(iDetLayerUsed, h->GetMean());
0576 if (spread)
0577 layerHistWithSpread->SetBinError(iDetLayerUsed, h->GetRMS());
0578 if (verbose) {
0579 h->SetTitle(("SurfaceDeformation " + whichOne) += " " + NameSurfDef(iPar) +=
0580 this->TitleAdd() + ";" + NameSurfDef(iPar) += UnitSurfDef(iPar));
0581 if (verbose2)
0582 fHistManager->AddHist(h, firstLayer + 2 + iPar);
0583 else
0584 fHistManager->AddHist(h, firstLayer + 2 + iDetLayerUsed);
0585 } else
0586 delete h;
0587 }
0588 layerHist->LabelsDeflate();
0589 if (spread)
0590 layerHistWithSpread->LabelsDeflate();
0591 layerHistRms->LabelsDeflate();
0592 if (layerHist->GetEntries()) {
0593 fHistManager->AddHistSame(layerHist, firstLayer, iParUsed, (spread ? "error: #sigma(mean)" : ""));
0594 if (spread)
0595 fHistManager->AddHistSame(layerHistWithSpread, firstLayer, iParUsed, "error: RMS");
0596 fHistManager->AddHist(layerHistRms, firstLayer + 1);
0597
0598 ++iParUsed;
0599 } else {
0600 delete layerHist;
0601 delete layerHistRms;
0602 }
0603 }
0604
0605 fHistManager->Draw();
0606
0607 this->SetSubDetId(-1);
0608 this->ClearAdditionalSel();
0609 }
0610
0611
0612 bool PlotMillePede::SetDetLayerCuts(unsigned int detLayer, bool silent) {
0613 if (!silent) {
0614
0615 if (fAdditionalSelTitle.Length()) {
0616 ::Warning("PlotMillePede::SetDetLayerCuts", "Clearing selection '%s'!", fAdditionalSelTitle.Data());
0617 }
0618 if (fSubDetIds.GetSize()) {
0619 ::Warning("PlotMillePede::SetDetLayerCuts", "Possibly changing subdet selection!");
0620 }
0621 }
0622 this->ClearAdditionalSel();
0623
0624 switch (detLayer) {
0625 case 0:
0626 this->SetSubDetId(1);
0627 this->AddAdditionalSel("r", 0., 5.5);
0628 return true;
0629 case 1:
0630 this->SetSubDetId(1);
0631 this->AddAdditionalSel("r", 5.5, 8.5);
0632 return true;
0633 case 2:
0634 this->SetSubDetId(1);
0635 this->AddAdditionalSel("r", 8.5, 12.);
0636 return true;
0637
0638 case 3:
0639 case 4:
0640 this->SetSubDetId(2);
0641 return false;
0642 break;
0643
0644 case 5:
0645 this->SetSubDetId(3);
0646 this->AddAdditionalSel("r", 20., 30.);
0647 this->AddAdditionalSel("StripRphi");
0648 return true;
0649 case 6:
0650 this->SetSubDetId(3);
0651 this->AddAdditionalSel("r", 20., 30.);
0652 this->AddAdditionalSel("StripStereo");
0653 return true;
0654 case 7:
0655 this->SetSubDetId(3);
0656 this->AddAdditionalSel("r", 30., 38.);
0657 this->AddAdditionalSel("StripRphi");
0658 return true;
0659 case 8:
0660 this->SetSubDetId(3);
0661 this->AddAdditionalSel("r", 30., 38.);
0662 this->AddAdditionalSel("StripStereo");
0663 return true;
0664 case 9:
0665 this->SetSubDetId(3);
0666 this->AddAdditionalSel("r", 38., 46.);
0667 return true;
0668 case 10:
0669 this->SetSubDetId(3);
0670 this->AddAdditionalSel("r", 46., 55.);
0671 return true;
0672
0673 case 11:
0674 this->SetSubDetId(4);
0675 this->AddAdditionalSel("r", 20., 33.);
0676 this->AddAdditionalSel("StripRphi");
0677 return true;
0678 case 12:
0679 this->SetSubDetId(4);
0680 this->AddAdditionalSel("r", 20., 33.);
0681 this->AddAdditionalSel("StripStereo");
0682 return true;
0683 case 13:
0684 this->SetSubDetId(4);
0685 this->AddAdditionalSel("r", 33., 41.);
0686 this->AddAdditionalSel("StripRphi");
0687 return true;
0688 case 14:
0689 this->SetSubDetId(4);
0690 this->AddAdditionalSel("r", 33., 41.);
0691 this->AddAdditionalSel("StripStereo");
0692 return true;
0693 case 15:
0694 this->SetSubDetId(4);
0695 this->AddAdditionalSel("r", 41., 50.);
0696 return true;
0697
0698
0699
0700
0701 case 16:
0702 this->SetSubDetId(6);
0703 this->AddAdditionalSel("r", 20., 33.);
0704 this->AddAdditionalSel("StripRphi");
0705 return true;
0706 case 17:
0707 this->SetSubDetId(6);
0708 this->AddAdditionalSel("r", 20., 33.);
0709 this->AddAdditionalSel("StripStereo");
0710 return true;
0711 case 18:
0712 this->SetSubDetId(6);
0713 this->AddAdditionalSel("r", 33., 41.);
0714 this->AddAdditionalSel("StripRphi");
0715 return true;
0716 case 19:
0717 this->SetSubDetId(6);
0718 this->AddAdditionalSel("r", 33., 41.);
0719 this->AddAdditionalSel("StripStereo");
0720 return true;
0721 case 20:
0722 this->SetSubDetId(6);
0723 this->AddAdditionalSel("r", 41., 50.);
0724 return true;
0725 case 21:
0726 this->SetSubDetId(6);
0727 this->AddAdditionalSel("r", 50., 60.);
0728 return true;
0729 case 22:
0730 this->SetSubDetId(6);
0731 this->AddAdditionalSel("r", 60., 75.);
0732 this->AddAdditionalSel("StripRphi");
0733 return true;
0734 case 23:
0735 this->SetSubDetId(6);
0736 this->AddAdditionalSel("r", 60., 75.);
0737 this->AddAdditionalSel("StripStereo");
0738 return true;
0739 case 24:
0740 this->SetSubDetId(6);
0741 this->AddAdditionalSel("r", 75., 90.);
0742 return true;
0743 case 25:
0744 this->SetSubDetId(6);
0745 this->AddAdditionalSel("r", 90., 120.);
0746 return true;
0747
0748 case 26:
0749 this->SetSubDetId(5);
0750 this->AddAdditionalSel("r", 50., 65.);
0751 this->AddAdditionalSel("StripRphi");
0752 return true;
0753 case 27:
0754 this->SetSubDetId(5);
0755 this->AddAdditionalSel("r", 50., 65.);
0756 this->AddAdditionalSel("StripStereo");
0757 return true;
0758 case 28:
0759 this->SetSubDetId(5);
0760 this->AddAdditionalSel("r", 65., 73.);
0761 this->AddAdditionalSel("StripRphi");
0762 return true;
0763 case 29:
0764 this->SetSubDetId(5);
0765 this->AddAdditionalSel("r", 65., 73.);
0766 this->AddAdditionalSel("StripStereo");
0767 return true;
0768 case 30:
0769 this->SetSubDetId(5);
0770 this->AddAdditionalSel("r", 73., 82.5);
0771 return true;
0772 case 31:
0773 this->SetSubDetId(5);
0774 this->AddAdditionalSel("r", 82.5, 92.);
0775 return true;
0776 case 32:
0777 this->SetSubDetId(5);
0778 this->AddAdditionalSel("r", 92., 102.);
0779 return true;
0780 case 33:
0781 this->SetSubDetId(5);
0782 this->AddAdditionalSel("r", 102., 120.);
0783 return true;
0784 }
0785
0786 return false;
0787 }
0788
0789
0790 TString PlotMillePede::DetLayerLabel(unsigned int detLayer) const {
0791 switch (detLayer) {
0792 case 0:
0793 return "BPIX L1";
0794 case 1:
0795 return "BPIX L2";
0796 case 2:
0797 return "BPIX L3";
0798
0799
0800
0801 case 5:
0802 return "TIB L1#phi";
0803 case 6:
0804 return "TIB L1s";
0805 case 7:
0806 return "TIB L2#phi";
0807 case 8:
0808 return "TIB L2s";
0809 case 9:
0810 return "TIB L3";
0811 case 10:
0812 return "TIB L4";
0813
0814 case 11:
0815 return "TID R1#phi";
0816 case 12:
0817 return "TID R1s";
0818 case 13:
0819 return "TID R2#phi";
0820 case 14:
0821 return "TID R2s";
0822 case 15:
0823 return "TID R3";
0824
0825 case 16:
0826 return "TEC R1#phi";
0827 case 17:
0828 return "TEC R1s";
0829 case 18:
0830 return "TEC R2#phi";
0831 case 19:
0832 return "TEC R2s";
0833 case 20:
0834 return "TEC R3";
0835 case 21:
0836 return "TEC R4";
0837 case 22:
0838 return "TEC R5#phi";
0839 case 23:
0840 return "TEC R5s";
0841 case 24:
0842 return "TEC R6";
0843 case 25:
0844 return "TEC R7";
0845
0846 case 26:
0847 return "TOB L1#phi";
0848 case 27:
0849 return "TOB L1s";
0850 case 28:
0851 return "TOB L2#phi";
0852 case 29:
0853 return "TOB L2s";
0854 case 30:
0855 return "TOB L3";
0856 case 31:
0857 return "TOB L4";
0858 case 32:
0859 return "TOB L5";
0860 case 33:
0861 return "TOB L6";
0862 }
0863
0864 return Form("unknown DetLayer %u", detLayer);
0865 }
0866
0867
0868 void PlotMillePede::DrawOrigParam(bool addPlots, const TString &sel) {
0869
0870 const Int_t layer = this->PrepareAdd(addPlots);
0871 const TString titleAdd = this->TitleAdd();
0872 TString aSel(sel);
0873 this->AddBasicSelection(aSel);
0874
0875 for (UInt_t iPar = 0; iPar < kNpar; ++iPar) {
0876 const TString toMum(this->ToMumMuRad(iPar));
0877
0878 const TString hNameB(this->Unique(Form("beforeO%d", iPar)));
0879 TH1 *hBefore = this->CreateHist(Parenth(MisParT() += Par(iPar)) += toMum, aSel, hNameB);
0880 if (0. == hBefore->GetEntries())
0881 continue;
0882 hBefore->SetTitle(DelName(iPar) += ": misplacement" + titleAdd + ";" + DelNameU(iPar) += ";#parameters");
0883 fHistManager->AddHist(hBefore, layer);
0884 }
0885
0886 fHistManager->Draw();
0887 }
0888
0889
0890 void PlotMillePede::DrawOrigPos(bool addPlots, const TString &sel) {
0891 const Int_t layer = this->PrepareAdd(addPlots);
0892 const TString titleAdd = this->TitleAdd();
0893 TString aSel(sel);
0894 this->AddBasicSelection(aSel);
0895
0896 const TString posNames[] = {"phi", "r", "x", "y", "z"};
0897 for (UInt_t iPos = 0; iPos < sizeof(posNames) / sizeof(posNames[0]); ++iPos) {
0898 const TString &pos = posNames[iPos];
0899 TH1 *h = this->CreateHist(OrgPos(pos), aSel, this->Unique("org_" + pos));
0900 h->SetTitle("original position " + Name(pos) + titleAdd + ";" + Name(pos) + ";#alignables");
0901 fHistManager->AddHist(h, layer);
0902 }
0903
0904 fHistManager->Draw();
0905 }
0906
0907
0908 void PlotMillePede::DrawParamResult(Option_t *option) {
0909 const TString opt(option);
0910 const Int_t layer = this->PrepareAdd(opt.Contains("add", TString::kIgnoreCase));
0911
0912 const TString titleAdd = this->TitleAdd();
0913 UInt_t nPlot = 0;
0914 for (UInt_t iPar = 0; iPar < kNpar; ++iPar) {
0915 TString sel(opt.Contains("withfixed", TString::kIgnoreCase) ? "" : Fixed(iPar, false).Data());
0916 this->AddBasicSelection(sel);
0917
0918 const TString toMu(this->ToMumMuRad(iPar));
0919 const TString finalMis(this->FinalMisAlignment(iPar) += toMu);
0920 const TString startMis(Parenth(MisParT() += Par(iPar)) += toMu);
0921
0922 const TString hNameB(this->Unique(Form("before%d", iPar)) += Form("(%d,%f,%f)", fNbins, fMaxDevDown, fMaxDevUp));
0923 TH1 *hBef = this->CreateHist(startMis, sel, hNameB);
0924 const TString hNameD(
0925 this->Unique(Form("end%d", iPar)) +=
0926 Form("(%d,%f,%f)", hBef->GetNbinsX(), hBef->GetXaxis()->GetXmin(), hBef->GetXaxis()->GetXmax()));
0927 TH1 *hEnd = this->CreateHist(finalMis, sel, hNameD);
0928 const TString hName2D(this->Unique(Form("vs%d", iPar)) += Form("(30,%f,%f,30,-500,500)", fMaxDevDown, fMaxDevUp));
0929 TH1 *hVs = this->CreateHist(startMis + ":" + finalMis, sel, hName2D, "BOX");
0930 if (0. == hEnd->GetEntries())
0931 continue;
0932 hEnd->SetTitle(DelName(iPar) += titleAdd + ";" + DelNameU(iPar) += ";#parameters");
0933 hVs->SetTitle(DelName(iPar) += titleAdd + ";" + DelNameU(iPar) += "(end);" + DelNameU(iPar) += "(start)");
0934 if (this->GetTitle().Length() != 0) {
0935 fHistManager->AddHist(hEnd, layer, this->GetTitle());
0936 } else {
0937 fHistManager->AddHist(hEnd, layer, "remaining misal.");
0938 }
0939 fHistManager->AddHistSame(hBef, layer, nPlot, "misaligned");
0940 fHistManager->AddHist(hVs, layer + 1);
0941
0942 ++nPlot;
0943 }
0944
0945 fHistManager->Draw();
0946 }
0947
0948
0949 void PlotMillePede::DrawPosResult(bool addPlots, const TString &selection) {
0950 const Int_t layer = this->PrepareAdd(addPlots);
0951 TString sel(selection);
0952 this->AddBasicSelection(sel);
0953
0954 const TString posNames[] = {"rphi", "r", "z", "phi", "x", "y"};
0955
0956 const TString titleAdd = this->TitleAdd();
0957 UInt_t nPlot = 0;
0958 for (UInt_t iPos = 0; iPos < sizeof(posNames) / sizeof(posNames[0]); ++iPos) {
0959 const TString &posName = posNames[iPos];
0960
0961 const TString toMu(this->ToMumMuRad(posName));
0962 const TString hNameB(this->Unique(posName + "Before") += Form("(%d,%f,%f)", fNbins, fMaxDevDown, fMaxDevUp));
0963 const TString misPos(Parenth(DeltaPos(posName, MisPosT())) += toMu);
0964 TH1 *hBef = this->CreateHist(misPos, sel, hNameB);
0965 if (0. == hBef->GetEntries()) {
0966 delete hBef;
0967 continue;
0968 }
0969 TString hNameD(this->Unique(posName + "End"));
0970 this->CopyAddBinning(hNameD, hBef);
0971 const TString endPos(Parenth(DeltaPos(posName, PosT())) += toMu);
0972 TH1 *hEnd = this->CreateHist(endPos, sel, hNameD);
0973 const TString hName2D(this->Unique(posName + "Vs(30,-100,100,30,-500,500)"));
0974 TH1 *hVs = this->CreateHist2D(endPos, misPos, sel, hName2D, "BOX");
0975
0976 hEnd->SetTitle(DelName(posName) += titleAdd + ";" + DelNameU(posName) += ";#alignables");
0977 hVs->SetTitle(DelName(posName) += titleAdd + ";" + DelNameU(posName) += "(end);" + DelNameU(posName) += "(start)");
0978
0979 fHistManager->AddHist(hVs, layer);
0980 if (this->GetTitle().Length() != 0) {
0981 fHistManager->AddHist(hEnd, layer + 1, this->GetTitle());
0982 } else {
0983 fHistManager->AddHist(hEnd, layer + 1, "remaining misal.");
0984 }
0985
0986 fHistManager->AddHistSame(hBef, layer + 1, nPlot, "misaligned");
0987
0988 ++nPlot;
0989 }
0990
0991 fHistManager->Draw();
0992 }
0993
0994
0995 void PlotMillePede::DrawMisVsLocation(bool addPlots, const TString &sel, Option_t *option) {
0996 const Int_t layer = this->PrepareAdd(addPlots);
0997
0998 TString opt(option);
0999 opt.ToLower();
1000 int vsEuler = -1;
1001 if (opt.Contains("vse0"))
1002 vsEuler = 0;
1003 else if (opt.Contains("vse1"))
1004 vsEuler = 1;
1005
1006 const bool addStartMis = (opt.Contains("mis") ? true : false);
1007 const bool addFixed = opt.Contains("withfixed");
1008
1009 const TString titleAdd = this->TitleAdd();
1010 UInt_t nPlot = 0;
1011 for (UInt_t iPar = 0; iPar < kNpar; ++iPar) {
1012 TString fixSel(sel);
1013 if (!addFixed) {
1014 if (fixSel.Length())
1015 fixSel.Prepend(Fixed(iPar, false) += AndL());
1016 else
1017 fixSel = Fixed(iPar, false);
1018 }
1019 this->AddBasicSelection(fixSel);
1020
1021 const TString misPar(this->FinalMisAlignment(iPar) += ToMumMuRad(iPar));
1022 const TString startMisPar(MisParT() += Par(iPar) += ToMumMuRad(iPar));
1023
1024 const TString nDpr(this->Unique(Form("misParR%d", iPar)));
1025 TH2 *hMisParR = this->CreateHist2D(RPos(OrgPosT()), misPar, fixSel, nDpr, "BOX");
1026 TProfile *hProfParR = this->CreateHistProf(RPos(OrgPosT()), misPar, fixSel, "prof" + nDpr);
1027 if (0. == hMisParR->GetEntries())
1028 continue;
1029
1030
1031 TProfile *hProfParStartR =
1032 !addStartMis ? 0 : this->CreateHistProf(RPos(OrgPosT()), startMisPar, fixSel, "profStart" + nDpr);
1033
1034 const TString nDpz(this->Unique(Form("misParZ%d", iPar)));
1035 TH2 *hMisParZ = this->CreateHist2D(OrgPosT() += ZPos(), misPar, fixSel, nDpz, "BOX");
1036 TProfile *hProfParZ = this->CreateHistProf(OrgPosT() += ZPos(), misPar, fixSel, "prof" + nDpz);
1037
1038 TProfile *hProfParStartZ =
1039 !addStartMis ? 0 : this->CreateHistProf(OrgPosT() += ZPos(), startMisPar, fixSel, "profStart" + nDpz);
1040
1041 const TString nDpp(this->Unique(Form("misParPhi%d", iPar)));
1042 TH2 *hMisParPhi = this->CreateHist2D(Phi(OrgPosT()), misPar, fixSel, nDpp, "BOX");
1043 TProfile *hProfParPhi = this->CreateHistProf(Phi(OrgPosT()), misPar, fixSel, "prof" + nDpp);
1044
1045 TProfile *hProfParStartPhi =
1046 !addStartMis ? 0 : this->CreateHistProf(Phi(OrgPosT()), startMisPar, fixSel, "profStart" + nDpp);
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056 const TString nDpa0(this->Unique(Form("misParAlpha%d%d", vsEuler, iPar)));
1057 TH2 *hMisParAl0 = 0;
1058 if (vsEuler >= 0)
1059 hMisParAl0 = this->CreateHist2D(Alpha(OrgPosT(), (vsEuler == 0)), misPar, fixSel, nDpa0, "BOX");
1060 TProfile *hProfParAl0 = 0;
1061 if (vsEuler >= 0)
1062 hProfParAl0 = this->CreateHistProf(Alpha(OrgPosT(), (vsEuler == 0)), misPar, fixSel, "prof" + nDpa0);
1063
1064
1065 const TString nDpb0(this->Unique(Form("misParBeta%d%d", vsEuler, iPar)));
1066 TH2 *hMisParBet0 = 0;
1067 if (vsEuler >= 0)
1068 hMisParBet0 = this->CreateHist2D(Beta(OrgPosT(), (vsEuler == 0)), misPar, fixSel, nDpb0, "BOX");
1069 TProfile *hProfParBet0 = 0;
1070 if (vsEuler >= 0)
1071 hProfParBet0 = this->CreateHistProf(Beta(OrgPosT(), (vsEuler == 0)), misPar, fixSel, "prof" + nDpb0);
1072
1073
1074 const TString nDpg0(this->Unique(Form("misParGamma%d%d", vsEuler, iPar)));
1075 TH2 *hMisParGam0 = 0;
1076 if (vsEuler >= 0)
1077 hMisParGam0 = this->CreateHist2D(Gamma(OrgPosT(), (vsEuler == 0)), misPar, fixSel, nDpg0, "BOX");
1078 TProfile *hProfParGam0 = 0;
1079 if (vsEuler >= 0)
1080 hProfParGam0 = this->CreateHistProf(Gamma(OrgPosT(), (vsEuler == 0)), misPar, fixSel, "prof" + nDpg0);
1081
1082
1083 hMisParR->SetTitle(DelName(iPar) += " vs. r" + titleAdd + ";r[cm];" + DelNameU(iPar));
1084
1085 hMisParZ->SetTitle(DelName(iPar) += " vs. z" + titleAdd + ";z[cm];" + DelNameU(iPar));
1086 hMisParPhi->SetTitle(DelName(iPar) += +" vs. #phi" + titleAdd + ";#phi;" + DelNameU(iPar));
1087
1088 if (hMisParAl0)
1089 hMisParAl0->SetTitle(DelName(iPar) +=
1090 Form(" vs. euler #alpha^{%d};#alpha^{%d};", vsEuler, vsEuler) + DelNameU(iPar));
1091 if (hMisParBet0)
1092 hMisParBet0->SetTitle(DelName(iPar) +=
1093 Form(" vs. euler #beta^{%d};#beta^{%d};", vsEuler, vsEuler) + DelNameU(iPar));
1094 if (hMisParGam0)
1095 hMisParGam0->SetTitle(DelName(iPar) +=
1096 Form(" vs. euler #gamma^{%d};#gamma^{%d};", vsEuler, vsEuler) + DelNameU(iPar));
1097 hProfParR->SetTitle("#LT" + DelName(iPar) += "#GT vs. r" + titleAdd + ";r[cm];" + DelNameU(iPar));
1098 hProfParZ->SetTitle("#LT" + DelName(iPar) += "#GT vs. z" + titleAdd + ";z[cm];" + DelNameU(iPar));
1099 hProfParPhi->SetTitle("#LT" + DelName(iPar) += "#GT vs. #phi" + titleAdd + ";#phi;" + DelNameU(iPar));
1100
1101 if (hProfParAl0)
1102 hProfParAl0->SetTitle("#LT" + DelName(iPar) +=
1103 Form("#GT vs. euler #alpha^{%d};#alpha^{%d};", vsEuler, vsEuler) + DelNameU(iPar));
1104 if (hProfParBet0)
1105 hProfParBet0->SetTitle("#LT" + DelName(iPar) +=
1106 Form("#GT vs. euler #beta^{%d};#beta^{%d};", vsEuler, vsEuler) + DelNameU(iPar));
1107 if (hProfParGam0)
1108 hProfParGam0->SetTitle("#LT" + DelName(iPar) +=
1109 Form("#GT vs. euler #gamma^{%d};#gamma^{%d};", vsEuler, vsEuler) + DelNameU(iPar));
1110 if (addStartMis) {
1111 hProfParStartR->SetTitle("#LT" + DelName(iPar) += "#GT vs. r (start);r[cm];" + DelNameU(iPar));
1112 hProfParStartZ->SetTitle("#LT" + DelName(iPar) += "#GT vs. z (start);z[cm];" + DelNameU(iPar));
1113 hProfParStartPhi->SetTitle("#LT" + DelName(iPar) += "#GT vs. #phi (start);#phi;" + DelNameU(iPar));
1114
1115 }
1116
1117 fHistManager->AddHist(hMisParR, layer + nPlot);
1118 fHistManager->AddHist(hProfParR, layer + nPlot);
1119 if (addStartMis)
1120 fHistManager->AddHist(hProfParStartR, layer + nPlot);
1121
1122 fHistManager->AddHist(hMisParZ, layer + nPlot);
1123 fHistManager->AddHist(hProfParZ, layer + nPlot);
1124 if (addStartMis)
1125 fHistManager->AddHist(hProfParStartZ, layer + nPlot);
1126
1127 fHistManager->AddHist(hMisParPhi, layer + nPlot);
1128 fHistManager->AddHist(hProfParPhi, layer + nPlot);
1129 if (addStartMis)
1130 fHistManager->AddHist(hProfParStartPhi, layer + nPlot);
1131
1132
1133
1134
1135
1136 if (hMisParAl0)
1137 fHistManager->AddHist(hMisParAl0, layer + nPlot + 1);
1138 if (hProfParAl0)
1139 fHistManager->AddHist(hProfParAl0, layer + nPlot + 1);
1140
1141 if (hMisParBet0)
1142 fHistManager->AddHist(hMisParBet0, layer + nPlot + 1);
1143 if (hProfParBet0)
1144 fHistManager->AddHist(hProfParBet0, layer + nPlot + 1);
1145
1146 if (hMisParGam0)
1147 fHistManager->AddHist(hMisParGam0, layer + nPlot + 1);
1148 if (hProfParGam0)
1149 fHistManager->AddHist(hProfParGam0, layer + nPlot + 1);
1150
1151 fHistManager->SetNumHistsXY((addStartMis ? 3 : 2), 4, layer + nPlot);
1152 if (vsEuler >= 0) {
1153 fHistManager->SetNumHistsXY(2, 3, layer + nPlot + 1);
1154 ++nPlot;
1155 }
1156 ++nPlot;
1157 }
1158
1159 fHistManager->Draw();
1160 }
1161
1162
1163 void PlotMillePede::DrawPosMisVsLocation(bool addPlots, const TString &selection, Option_t *option) {
1164 const Int_t layer = this->PrepareAdd(addPlots);
1165
1166 TString opt(option);
1167 opt.ToLower();
1168 const bool addStart = (opt.Contains("start") ? true : false);
1169
1170 TString sel(selection);
1171 this->AddBasicSelection(sel);
1172
1173 const TString posNames[] = {"rphi", "r", "z", "x", "y", "phi"};
1174
1175 const TString titleAdd = this->TitleAdd();
1176 UInt_t nPlot = 0;
1177 for (UInt_t iPos = 0; iPos < sizeof(posNames) / sizeof(posNames[0]); ++iPos) {
1178 const TString &posName = posNames[iPos];
1179
1180 const TString toMu(this->ToMumMuRad(posName));
1181 const TString startPos(Parenth(DeltaPos(posName, MisPosT())) += toMu);
1182 const TString endPos(Parenth(DeltaPos(posName, PosT())) += toMu);
1183
1184 const TString nDpr(this->Unique(posName + "MisPosR"));
1185 TH2 *hEndPosR = this->CreateHist2D(RPos(OrgPosT()), endPos, sel, nDpr, "BOX");
1186 if (0. == hEndPosR->GetEntries()) {
1187 delete hEndPosR;
1188 continue;
1189 }
1190 TProfile *hProfPosR = this->CreateHistProf(RPos(OrgPosT()), endPos, sel, "prof" + nDpr);
1191 TProfile *hProfPosStartR = !addStart ? 0 : this->CreateHistProf(RPos(OrgPosT()), startPos, sel, "profStart" + nDpr);
1192
1193 const TString nDpz(this->Unique(posName + "MisPosZ"));
1194 TH2 *hEndPosZ = this->CreateHist2D(OrgPosT() += ZPos(), endPos, sel, nDpz, "BOX");
1195 TProfile *hProfPosZ = this->CreateHistProf(OrgPosT() += ZPos(), endPos, sel, "prof" + nDpz);
1196 TProfile *hProfPosStartZ =
1197 !addStart ? 0 : this->CreateHistProf(OrgPosT() += ZPos(), startPos, sel, "profStart" + nDpz);
1198
1199 const TString nDpp(this->Unique(posName + "MisPosPhi"));
1200 TH2 *hEndPosPhi = this->CreateHist2D(Phi(OrgPosT()), endPos, sel, nDpp, "BOX");
1201 TProfile *hProfPosPhi = this->CreateHistProf(Phi(OrgPosT()), endPos, sel, "prof" + nDpp);
1202 TProfile *hProfPosStartPhi =
1203 !addStart ? 0 : this->CreateHistProf(Phi(OrgPosT()), startPos, sel, "profStart" + nDpp);
1204
1205 const TString nDpy(this->Unique(posName + "MisPosY"));
1206 TH2 *hEndPosY = this->CreateHist2D(OrgPosT() += YPos(), endPos, sel, nDpy, "BOX");
1207 TProfile *hProfPosY = this->CreateHistProf(OrgPosT() += YPos(), endPos, sel, "prof" + nDpy);
1208 TProfile *hProfPosStartY =
1209 !addStart ? 0 : this->CreateHistProf(OrgPosT() += YPos(), startPos, sel, "profStart" + nDpy);
1210
1211 hEndPosR->SetTitle(DelName(posName) += " vs. r" + titleAdd + ";r[cm];" + DelNameU(posName));
1212 hEndPosZ->SetTitle(DelName(posName) += " vs. z" + titleAdd + ";z[cm];" + DelNameU(posName));
1213 hEndPosPhi->SetTitle(DelName(posName) += " vs. #phi" + titleAdd + ";#phi;" + DelNameU(posName));
1214 hEndPosY->SetTitle(DelName(posName) += " vs. y" + titleAdd + ";y[cm];" + DelNameU(posName));
1215 hProfPosR->SetTitle("#LT" + DelName(posName) += "#GT vs. r" + titleAdd + ";r[cm];" + DelNameU(posName));
1216 hProfPosZ->SetTitle("#LT" + DelName(posName) += "#GT vs. z" + titleAdd + ";z[cm];" + DelNameU(posName));
1217 hProfPosPhi->SetTitle("#LT" + DelName(posName) += "#GT vs. #phi" + titleAdd + ";#phi;" + DelNameU(posName));
1218 hProfPosY->SetTitle("#LT" + DelName(posName) += "#GT vs. y" + titleAdd + ";y[cm];" + DelNameU(posName));
1219 if (addStart) {
1220 hProfPosStartR->SetTitle("#LT" + DelName(posName) +=
1221 "#GT vs. r (start)" + titleAdd + ";r[cm];" + DelNameU(posName));
1222 hProfPosStartZ->SetTitle("#LT" + DelName(posName) += "#GT vs. z" + titleAdd + ";z[cm];" + DelNameU(posName));
1223 hProfPosStartPhi->SetTitle("#LT" + DelName(posName) += "#GT vs. #phi" + titleAdd + ";#phi;" + DelNameU(posName));
1224 hProfPosStartY->SetTitle("#LT" + DelName(posName) += "#GT vs. y" + titleAdd + ";y[cm];" + DelNameU(posName));
1225 }
1226
1227 fHistManager->AddHist(hEndPosR, layer + nPlot);
1228 fHistManager->AddHist(hProfPosR, layer + nPlot);
1229 if (addStart)
1230 fHistManager->AddHist(hProfPosStartR, layer + nPlot);
1231
1232 fHistManager->AddHist(hEndPosZ, layer + nPlot);
1233 fHistManager->AddHist(hProfPosZ, layer + nPlot);
1234 if (addStart)
1235 fHistManager->AddHist(hProfPosStartZ, layer + nPlot);
1236
1237 fHistManager->AddHist(hEndPosPhi, layer + nPlot);
1238 fHistManager->AddHist(hProfPosPhi, layer + nPlot);
1239 if (addStart)
1240 fHistManager->AddHist(hProfPosStartPhi, layer + nPlot);
1241
1242 fHistManager->AddHist(hEndPosY, layer + nPlot);
1243 fHistManager->AddHist(hProfPosY, layer + nPlot);
1244 if (addStart)
1245 fHistManager->AddHist(hProfPosStartY, layer + nPlot);
1246
1247 fHistManager->SetNumHistsXY((addStart ? 3 : 2), 4, layer + nPlot);
1248
1249 ++nPlot;
1250 }
1251
1252 fHistManager->Draw();
1253 }
1254
1255
1256 void PlotMillePede::DrawLabelDiffAbove(UInt_t iPar, float minDiff, bool addPlots) {
1257 const Int_t layer = this->PrepareAdd(addPlots);
1258 const TString titleAdd = this->TitleAdd();
1259 const TString misPar(this->FinalMisAlignment(iPar) += ToMumMuRad(iPar));
1260
1261 TString fixSel(Fixed(iPar, false) += AndL() += Abs(misPar)
1262 += Form(">%f", minDiff));
1263 this->AddBasicSelection(fixSel);
1264
1265 const TString name(this->Unique(Form("labelBigMis%d(100000,0,100000)", iPar)));
1266
1267 TH1 *hLabel = this->CreateHist(MpT() += "Label", fixSel, name);
1268
1269
1270 if (0. == hLabel->GetEntries())
1271 return;
1272
1273 hLabel->SetTitle("Label, " + Abs(DelNameU(iPar)) += Form(">%f", minDiff) + titleAdd + ";label");
1274
1275 fHistManager->AddHist(hLabel, layer);
1276
1277 fHistManager->Draw();
1278 }
1279
1280
1281 void PlotMillePede::DrawGlobCorr(bool addPlots, const TString &sel, Option_t *option, Float_t min, Float_t max) {
1282
1283
1284 const TString opt(option);
1285 const bool nal = opt.Contains("nal", TString::kIgnoreCase);
1286 const Int_t layer = this->PrepareAdd(addPlots);
1287 const TString titleAdd = this->TitleAdd();
1288
1289 UInt_t nPlot = 0;
1290 for (UInt_t iPar = 0; iPar < kNpar; ++iPar) {
1291 TString aSel(sel);
1292 if (aSel.Length())
1293 aSel.Prepend(Fixed(iPar, false) += AndL());
1294 else
1295 aSel = Fixed(iPar, false);
1296 if (opt.Contains("valid", TString::kIgnoreCase)) {
1297 if (aSel.Length())
1298 aSel.Prepend(Valid(iPar) += AndL());
1299 else
1300 aSel = Valid(iPar);
1301 }
1302
1303 this->AddBasicSelection(aSel);
1304
1305 const TString limits(nal ? "" : Form("(100,%f,%f)", min, max));
1306 const TString hName(this->Unique(Form("globCor%d", iPar)) += limits);
1307 TH1 *h = this->CreateHist(Cor(iPar), aSel, hName);
1308 if (0. == h->GetEntries())
1309 continue;
1310
1311 const TString limits2dR(nal ? "" : Form("(110,0,110, 100,%f,%f)", min, max));
1312 const TString hNameR(this->Unique(Form("rGlobCor%d", iPar)) += limits2dR);
1313 TH1 *hR = this->CreateHist2D(RPos(OrgPosT()), Cor(iPar), aSel, hNameR, "BOX");
1314
1315 const TString limits2dZ(nal ? "" : Form("(110,-275,275,100,%f,%f)", min, max));
1316 const TString hNameZ(this->Unique(Form("zGlobCor%d", iPar)) += limits2dZ);
1317 TH1 *hZ = this->CreateHist2D(OrgPosT() += ZPos(), Cor(iPar), aSel, hNameZ, "BOX");
1318
1319 const TString limits2dPhi(nal ? "" : Form("(100,-3.15,3.15,100,%f,%f)", min, max));
1320 const TString hNamePhi(this->Unique(Form("phiGlobCor%d", iPar)) += limits2dPhi);
1321 TH1 *hPhi = this->CreateHist2D(Phi(OrgPosT()), Cor(iPar), aSel, hNamePhi, "BOX");
1322
1323 h->SetTitle(this->DelName(iPar) += titleAdd + ";Global Correlation;#parameters");
1324 hR->SetTitle(this->DelName(iPar) += titleAdd + ";r[cm];Global Correlation");
1325 hZ->SetTitle(this->DelName(iPar) += titleAdd + ";z[cm];Global Correlation");
1326 hPhi->SetTitle(this->DelName(iPar) += titleAdd + ";#phi;Global Correlation");
1327 fHistManager->AddHist(h, layer);
1328 fHistManager->AddHist(hR, layer + nPlot + 1);
1329 fHistManager->AddHist(hZ, layer + nPlot + 1);
1330 fHistManager->AddHist(hPhi, layer + nPlot + 1);
1331 ++nPlot;
1332 }
1333
1334 fHistManager->Draw();
1335 }
1336
1337
1338 void PlotMillePede::DrawPull(Option_t *option) {
1339 const TString opt(option);
1340 const Int_t layer = this->PrepareAdd(opt.Contains("add", TString::kIgnoreCase));
1341
1342 const TString titleAdd = this->TitleAdd();
1343 UInt_t nPlot = 0;
1344 for (UInt_t iPar = 0; iPar < kNpar; ++iPar) {
1345 const TString misPar(this->FinalMisAlignment(iPar));
1346
1347 TString sel(Fixed(iPar, false));
1348 if (opt.Contains("valid", TString::kIgnoreCase)) {
1349 sel += AndL() += Valid(iPar);
1350 }
1351 this->AddBasicSelection(sel);
1352
1353 const TString hNameS(this->Unique(Form("sigma%d", iPar)));
1354 TH1 *hSi = this->CreateHist(ParSi(iPar) += ToMumMuRad(iPar), sel, hNameS);
1355 const TString hNameP(this->Unique(Form("pull%d", iPar)) +=
1356 (opt.Contains("nolimit", TString::kIgnoreCase) ? "" : "(100,-6,6)"));
1357 TH1 *hPull = this->CreateHist(misPar + Div() += ParSi(iPar), sel, hNameP);
1358 if (0. == hPull->GetEntries())
1359 continue;
1360
1361 hPull->SetTitle("pull " + DelName(iPar) += titleAdd + ";#Delta/" + Fun("#sigma", DelName(iPar)) += ";#parameters");
1362 hPull->Fit("gaus", "Q0L");
1363 hPull->GetFunction("gaus")->ResetBit(TF1::kNotDraw);
1364 hSi->SetTitle(DelName(iPar) += titleAdd + ";" + Fun("#sigma", DelName(iPar)) += Unit(iPar) += ";#parameters");
1365
1366 fHistManager->AddHist(hPull, layer);
1367 fHistManager->AddHist(hSi, layer + 1);
1368 ++nPlot;
1369 }
1370
1371 fHistManager->Draw();
1372 }
1373
1374
1375 void PlotMillePede::DrawErrorVsHit(bool addPlots, const TString &sel) {
1376 const Int_t layer = this->PrepareAdd(addPlots);
1377 const TString titleAdd = this->TitleAdd();
1378
1379 TString andSel(sel);
1380 if (andSel.Length())
1381 andSel.Prepend(AndL());
1382
1383 UInt_t nPlot = 0;
1384 for (UInt_t iPar = 0; iPar < kNpar; ++iPar) {
1385 TString aSel(sel);
1386 if (aSel.Length())
1387 aSel.Prepend(Fixed(iPar, false) += AndL());
1388 else
1389 aSel = Fixed(iPar, false);
1390 this->AddBasicSelection(aSel);
1391
1392 const TString toMum(this->ToMumMuRad(iPar));
1393 const TString nMiPaHit(this->Unique(Form("mis%dVsHit", iPar)));
1394 const TString misPar(this->FinalMisAlignment(iPar) += toMum);
1395
1396 TH2 *hDParVsHit = this->CreateHist2D(misPar, HitsX(), aSel, nMiPaHit, "BOX");
1397 if (0. == hDParVsHit->GetEntries())
1398 continue;
1399
1400 const TString nSiHit(this->Unique(Form("sigmaVsHit%d", iPar)));
1401 TH2 *hSiVsHit = this->CreateHist2D(ParSi(iPar) += toMum, HitsX(), ParSiOk(iPar) += AndL() += aSel, nSiHit, "BOX");
1402
1403 const TString nSiDiHit(this->Unique(Form("sigmaByHit%d", iPar)));
1404 TH1 *hSiByHit =
1405 this->CreateHist((ParSi(iPar) += toMum) += Div() += Sqrt(HitsX()), ParSiOk(iPar) += AndL() += aSel, nSiDiHit);
1406
1407 hDParVsHit->SetTitle(DelName(iPar) += ": remaining misalign vs. N_{hit,x}" + titleAdd + ";" +
1408 Fun("#Delta", DelName(iPar)) += Unit(iPar) += ";N_{hit,x}");
1409 hSiVsHit->SetTitle(Fun("#sigma", DelName(iPar)) +=
1410 " vs. N_{hit,x}" + titleAdd + ";" + Fun("#sigma", DelName(iPar)) += Unit(iPar) += ";N_{hit,x}");
1411 hSiByHit->SetTitle(Fun("#sigma", DelName(iPar)) += Div() += Sqrt("N_{hit,x}") +=
1412 titleAdd + ";" + Fun("#sigma", DelName(iPar)) += Div() += Sqrt("N_{hit,x}") += Unit(iPar) +=
1413 ";#parameters");
1414
1415 fHistManager->AddHist(hDParVsHit, layer);
1416
1417 if (hSiVsHit->GetEntries())
1418 fHistManager->AddHist(hSiVsHit, layer + 1);
1419 if (hSiByHit->GetEntries())
1420 fHistManager->AddHist(hSiByHit, layer + 2);
1421
1422 ++nPlot;
1423 }
1424
1425 fHistManager->Draw();
1426 }
1427
1428
1429 void PlotMillePede::DrawHitMaps(bool addPlots, bool inclFullFixed) {
1430 const Int_t layer = this->PrepareAdd(addPlots);
1431 TString sel(inclFullFixed ? "" : AnyFreePar().Data());
1432 this->AddBasicSelection(sel);
1433 const TString titleAdd = this->TitleAdd();
1434
1435 Option_t *drawOpt = "BOX";
1436
1437 TH2 *hRx = this->CreateHist2D(RPos(OrgPosT()), HitsX(), sel, this->Unique("xHitsR"), drawOpt);
1438 TProfile *hRxProf = this->CreateHistProf(RPos(OrgPosT()), HitsX(), sel, this->Unique("profxHitsR"));
1439
1440 TH2 *hRy = this->CreateHist2D(RPos(OrgPosT()), HitsY(), sel, this->Unique("yHitsR"), drawOpt);
1441 TProfile *hRyProf = this->CreateHistProf(RPos(OrgPosT()), HitsY(), sel, this->Unique("profyHitsR"));
1442
1443
1444 TH2 *hZx = this->CreateHist2D(OrgPosT() += ZPos(), HitsX(), sel, this->Unique("xHitsZ"), drawOpt);
1445 TProfile *hZxProf = this->CreateHistProf(OrgPosT() += ZPos(), HitsX(), sel, this->Unique("profxHitsZ"));
1446
1447 TH2 *hZy = this->CreateHist2D(OrgPosT() += ZPos(), HitsY(), sel, this->Unique("yHitsZ"), drawOpt);
1448 TProfile *hZyProf = this->CreateHistProf(OrgPosT() += ZPos(), HitsY(), sel, this->Unique("profyHitsZ"));
1449
1450
1451 TH2 *hPhiX = this->CreateHist2D(Phi(OrgPosT()), HitsX(), sel, this->Unique("xHitsPhi"), drawOpt);
1452 TProfile *hPhixProf = this->CreateHistProf(Phi(OrgPosT()), HitsX(), sel, this->Unique("profxHitsPhi"));
1453
1454 TH2 *hPhiY = this->CreateHist2D(Phi(OrgPosT()), HitsY(), sel, this->Unique("yHitsPhi"), drawOpt);
1455 TProfile *hPhiyProf = this->CreateHistProf(Phi(OrgPosT()), HitsY(), sel, this->Unique("profyHitsPhi"));
1456
1457
1458 TString selHitY(HitsY());
1459 if (!sel.IsNull())
1460 selHitY += AndL() += sel;
1461 TH1 *hNhitXlog = this->CreateHist(Fun("TMath::Log10", HitsX()),
1462 sel + AndL() += Parenth(HitsX() += ">0"),
1463 this->Unique("NhitXlog"));
1464 TH1 *hNhitX = this->CreateHist(HitsX(), sel, this->Unique("NhitX"));
1465 TH1 *hNhitY = this->CreateHist(HitsY(), sel, this->Unique("NhitY"));
1466 TH2 *hNhitXvsY = this->CreateHist2D(HitsX(), HitsY(), sel, this->Unique("NhitXvsY"), drawOpt);
1467 TH1 *hNhitDiffXy = this->CreateHist(HitsX() += Min() += HitsY(), selHitY, this->Unique("NhitDiffXy"));
1468 TH1 *hNhitDiffXyVsR =
1469 this->CreateHist2D(RPos(OrgPosT()), HitsX() += Min() += HitsY(), selHitY, this->Unique("NhitDiffXyVsR"), drawOpt);
1470
1471
1472 drawOpt = "COLZ";
1473 if (!sel.IsNull())
1474 sel = Parenth(sel).Prepend(Mal());
1475 TH2 *hTotPhiRx = this->CreateHist2D(OrgPosT() += XPos(),
1476 OrgPosT() += YPos(),
1477 HitsX() += sel,
1478 this->Unique("totXvsXY"),
1479 drawOpt);
1480 TH2 *hTotPhiRy =
1481 this->CreateHist2D(OrgPosT() += XPos(), OrgPosT() += YPos(), HitsY() += sel, this->Unique("totYvsXY"), drawOpt);
1482 TH2 *hTotRzX =
1483 this->CreateHist2D(OrgPosT() += ZPos(), RPos(OrgPosT()), HitsX() += sel, this->Unique("totXvsZR"), drawOpt);
1484 TH2 *hTotRzY =
1485 this->CreateHist2D(OrgPosT() += ZPos(), RPos(OrgPosT()), HitsY() += sel, this->Unique("totYvsZR"), drawOpt);
1486
1487 hRx->SetTitle("#hits_{x} vs. r" + titleAdd + ";r[cm];N_{hit,x}");
1488 hRxProf->SetTitle("#LT#hits_{x}> vs. r" + titleAdd + ";r[cm];N_{hit,x}");
1489 hRy->SetTitle("#hits_{y} vs. r" + titleAdd + ";r[cm];N_{hit,y}");
1490 hRyProf->SetTitle("#LT#hits_{y}> vs. r" + titleAdd + ";r[cm];N_{hit,y}");
1491 hZx->SetTitle("#hits_{x} vs. z" + titleAdd + ";z[cm];N_{hit,x}");
1492 hZxProf->SetTitle("#LT#hits_{x}> vs. z" + titleAdd + ";z[cm];N_{hit,x}");
1493 hZy->SetTitle("#hits_{y} vs. z" + titleAdd + ";z[cm];N_{hit,y}");
1494 hZyProf->SetTitle("#LT#hits_{y}> vs. z" + titleAdd + ";z[cm];N_{hit,y}");
1495 hPhiX->SetTitle("#hits_{x} vs. #phi" + titleAdd + ";#phi;N_{hit,x}");
1496 hPhixProf->SetTitle("#LT#hits_{x}> vs. #phi" + titleAdd + ";#phi;N_{hit,x}");
1497 hPhiY->SetTitle("#hits_{y} vs. #phi" + titleAdd + ";#phi;N_{hit,y}");
1498 hPhiyProf->SetTitle("#LT#hits_{y}> vs. #phi" + titleAdd + ";#phi;N_{hit,y}");
1499
1500 hNhitXlog->SetTitle("#hits_{x}: log_{10}" + titleAdd + ";log_{10}(N_{hit,x})");
1501 hNhitX->SetTitle("#hits_{x}" + titleAdd + ";N_{hit,x}");
1502 hNhitY->SetTitle("#hits_{y}" + titleAdd + ";N_{hit,y}");
1503 hNhitXvsY->SetTitle("#hits_{x} vs. #hits_{y}" + titleAdd + ";N_{hit,x};N_{hit,y}");
1504 hNhitDiffXy->SetTitle("#hits_{x} - #hits_{y} (#hits_{y}#neq0)" + titleAdd + ";N_{hit,x}-N_{hit,y};#alignables");
1505 hNhitDiffXyVsR->SetTitle("#hits_{x} - #hits_{y} vs. r (#hits_{y}#neq0)" + titleAdd + ";r[cm];N_{hit,x}-N_{hit,y}");
1506
1507 hTotPhiRx->SetTitle("#hits_{x}" + titleAdd + ";x [cm];y [cm]");
1508 hTotPhiRy->SetTitle("#hits_{y}" + titleAdd + ";x [cm];y [cm]");
1509 hTotRzX->SetTitle("#hits_{x}" + titleAdd + ";z [cm];r [cm]");
1510 hTotRzY->SetTitle("#hits_{y}" + titleAdd + ";z [cm];r [cm]");
1511
1512 const bool addYhists = (hNhitY->GetMean() || hNhitY->GetRMS());
1513
1514 fHistManager->AddHist(hRx, layer);
1515 fHistManager->AddHist(hRxProf, layer + 1);
1516 if (addYhists) {
1517 fHistManager->AddHist(hRy, layer);
1518 fHistManager->AddHist(hRyProf, layer + 1);
1519 }
1520 fHistManager->AddHist(hZx, layer);
1521 fHistManager->AddHist(hZxProf, layer + 1);
1522 if (addYhists) {
1523 fHistManager->AddHist(hZy, layer);
1524 fHistManager->AddHist(hZyProf, layer + 1);
1525 }
1526 fHistManager->AddHist(hPhiX, layer);
1527 fHistManager->AddHist(hPhixProf, layer + 1);
1528 if (addYhists) {
1529 fHistManager->AddHist(hPhiY, layer);
1530 fHistManager->AddHist(hPhiyProf, layer + 1);
1531 }
1532
1533 fHistManager->AddHist(hNhitXlog, layer + 2);
1534 fHistManager->AddHist(hNhitX, layer + 2);
1535 fHistManager->AddHist(hNhitY, layer + 2);
1536 if (addYhists) {
1537 fHistManager->AddHist(hNhitXvsY, layer + 2);
1538 fHistManager->AddHist(hNhitDiffXy, layer + 2);
1539 fHistManager->AddHist(hNhitDiffXyVsR, layer + 2);
1540 }
1541
1542 fHistManager->AddHist(hTotPhiRx, layer + 3);
1543 if (addYhists)
1544 fHistManager->AddHist(hTotPhiRy, layer + 3);
1545 fHistManager->AddHist(hTotRzX, layer + 3);
1546 if (addYhists)
1547 if (hNhitY->GetEntries())
1548 fHistManager->AddHist(hTotRzY, layer + 3);
1549
1550 fHistManager->Draw();
1551 }
1552
1553
1554 void PlotMillePede::DrawBigPullLabel(float minPull, bool addPlots) {
1555 const Int_t layer = this->PrepareAdd(addPlots);
1556
1557
1558 for (UInt_t iPar = 0; iPar < kNpar; ++iPar) {
1559 const TString misPar(this->FinalMisAlignment(iPar));
1560
1561 const TString hNameL(Form("bigPullLabel%d(100000,0,100000)", iPar));
1562 const TString isBigPull(Abs(misPar + Div() += ParSi(iPar)) += Form(" > %f", minPull));
1563 TString sel(isBigPull + AndL() += Fixed(iPar, false));
1564 this->AddBasicSelection(sel);
1565
1566 TH1 *hLabel = this->CreateHist(Label(iPar), sel, hNameL);
1567 if (0. == hLabel->GetEntries())
1568 continue;
1569
1570 hLabel->SetTitle((DelName(iPar) += Form(": |pull| > %f", minPull)) += ";Label");
1571 fHistManager->AddHist(hLabel, layer);
1572
1573 }
1574
1575 fHistManager->Draw();
1576 }
1577
1578
1579 void PlotMillePede::DrawBigPullPos(float minPull, bool addPlots) {
1580 const Int_t layer = this->PrepareAdd(addPlots);
1581 const TString titleAdd = this->TitleAdd();
1582
1583 UInt_t nPlot = 0;
1584 for (UInt_t iPar = 0; iPar < kNpar; ++iPar) {
1585 const TString misPar(this->FinalMisAlignment(iPar));
1586
1587 const TString isBigPull(Abs(misPar + Div() += ParSi(iPar)) += Form(" > %f", minPull));
1588 TString sel(isBigPull + AndL() += Fixed(iPar, false));
1589 this->AddBasicSelection(sel);
1590
1591 const TString hNameR(this->Unique(Form("bigPullR%d", iPar)));
1592
1593 TH1 *hR = this->CreateHist(RPos(OrgPosT()), sel, hNameR);
1594 if (0. == hR->GetEntries())
1595 continue;
1596 const TString hNameZ(this->Unique(Form("bigPullZ%d", iPar)));
1597 TH1 *hZ = this->CreateHist(OrgPosT() += ZPos(), sel, hNameZ);
1598 const TString hNameP(this->Unique(Form("bigPullP%d", iPar)));
1599 TH1 *hPhi = this->CreateHist(Phi(OrgPosT()), sel, hNameP);
1600
1601 hR->SetTitle((DelName(iPar) += Form(": |pull| > %f", minPull)) += titleAdd + ";r[cm]");
1602 hZ->SetTitle((DelName(iPar) += Form(": |pull| > %f", minPull)) += titleAdd + ";z[cm]");
1603 hPhi->SetTitle((DelName(iPar) += Form(": |pull| > %f", minPull)) += titleAdd + ";#phi");
1604 fHistManager->AddHist(hR, layer + nPlot);
1605 fHistManager->AddHist(hZ, layer + nPlot);
1606 fHistManager->AddHist(hPhi, layer + nPlot);
1607 ++nPlot;
1608 }
1609
1610 fHistManager->Draw();
1611 }
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645 void PlotMillePede::DrawSubDetId(bool addPlots) {
1646 const Int_t layer = this->PrepareAdd(addPlots);
1647 const TString titleAdd = this->TitleAdd();
1648
1649 const TString nameAll(this->Unique("subDetId"));
1650 const TString nameAct(this->Unique("subDetIdActive"));
1651
1652 TString sel;
1653 this->AddBasicSelection(sel);
1654 TH1 *hAll = this->CreateHist(SubDetId(), sel, nameAll);
1655
1656 if (!sel.IsNull())
1657 sel = Parenth(sel) += AndL();
1658 sel += AnyFreePar();
1659 TH1 *hAct = this->CreateHist(SubDetId(), sel, nameAct);
1660
1661 if (hAll->GetEntries()) {
1662 hAll->SetTitle("subDetId" + titleAdd + ";ID(subdet)");
1663 fHistManager->AddHist(hAll, layer);
1664 }
1665 if (hAct->GetEntries()) {
1666 hAct->SetTitle("subDetId (any free param)" + titleAdd + ";ID(subdet)");
1667 fHistManager->AddHist(hAct, layer);
1668 }
1669
1670 fHistManager->Draw();
1671 }
1672
1673
1674 void PlotMillePede::DrawXyArrow(Double_t factor, Option_t *option) {
1675 const TString opt(option);
1676 const Int_t layer = this->PrepareAdd(opt.Contains("add", TString::kIgnoreCase));
1677
1678
1679 TString sel;
1680 this->AddBasicSelection(sel);
1681
1682 delete this->Draw(OrgPosT() += XPos() += ":" + OrgPosT() += YPos(), sel, "", "goff");
1683
1684 TTree *tree = this->GetMainTree();
1685 const Long64_t size = tree->GetSelectedRows();
1686
1687 if (size == 0)
1688 return;
1689
1690
1691 Double_t maxX = tree->GetV1()[TMath::LocMax(size, tree->GetV1())];
1692 maxX *= (maxX < 0 ? 0.9 : 1.1);
1693 Double_t minX = tree->GetV1()[TMath::LocMin(size, tree->GetV1())];
1694 minX *= (minX > 0 ? 0.9 : 1.1);
1695 Double_t maxY = tree->GetV2()[TMath::LocMax(size, tree->GetV2())];
1696 maxY *= (maxY < 0 ? 0.9 : 1.1);
1697 Double_t minY = tree->GetV2()[TMath::LocMin(size, tree->GetV2())];
1698 minY *= (minY > 0 ? 0.9 : 1.1);
1699 TH1 *hFrame = new TH2F(this->Unique("frame"),
1700 Form("scale %g%s;x [cm];y [cm]", factor, this->TitleAdd().Data()),
1701 10,
1702 minX,
1703 maxX,
1704 10,
1705 minY,
1706 maxY);
1707 hFrame->SetOption("AXIS");
1708 hFrame->SetEntries(size);
1709 fHistManager->AddHist(hFrame, layer);
1710
1711
1712 const std::vector<double> xs(tree->GetV1(), tree->GetV1() + size);
1713 const std::vector<double> ys(tree->GetV2(), tree->GetV2() + size);
1714
1715
1716 delete this->Draw(
1717 DeltaPos("x", PosT()) += ":" + DeltaPos("y", PosT()) += ":" + DeltaPos("z", PosT()), sel, "", "goff");
1718
1719 const std::vector<double> deltaXs(tree->GetV1(), tree->GetV1() + size);
1720 const std::vector<double> deltaYs(tree->GetV2(), tree->GetV2() + size);
1721
1722 if (opt.Contains("zcirc", TString::kIgnoreCase)) {
1723
1724 const std::vector<double> deltaZs(tree->GetV3(), tree->GetV3() + size);
1725 const Double_t rootFactor = TMath::Sqrt(TMath::Abs(factor));
1726
1727 for (unsigned int i = 0; i < size; ++i) {
1728 if (deltaZs[i] == 0.)
1729 continue;
1730 TEllipse *circ = new TEllipse(xs[i], ys[i], TMath::Abs(rootFactor * deltaZs[i]));
1731
1732 if (deltaZs[i] < 0.)
1733 circ->SetFillColor(kRed);
1734 else
1735 circ->SetFillColor(kGreen);
1736 fHistManager->AddObject(circ, layer, 0);
1737 }
1738 }
1739
1740
1741 TArrow::SetDefaultAngle(30.);
1742 for (unsigned int i = 0; i < size; ++i) {
1743 TArrow *arr = new TArrow(xs[i], ys[i], xs[i] + factor * deltaXs[i], ys[i] + factor * deltaYs[i], 0.01, "|>");
1744 fHistManager->AddObject(arr, layer, 0);
1745 }
1746
1747 fHistManager->Draw();
1748 }
1749
1750
1751 void PlotMillePede::ScanSelection(const char *sel, const char *addColumns) {
1752 TString realSel;
1753 if (sel) {
1754 realSel = sel;
1755 } else {
1756 this->AddBasicSelection(realSel);
1757 const TString titleAdd(this->TitleAdd());
1758 if (titleAdd.Length()) {
1759 std::cout << "Active selection: " << titleAdd << std::endl;
1760 }
1761 }
1762
1763 const TString mpPar(MpT() += Par());
1764
1765 TString scan("Id:Pos:" + mpPar += ":HitsX:Sigma:Label");
1766 if (addColumns)
1767 scan += addColumns;
1768 this->GetMainTree()->Scan(scan, realSel);
1769 }
1770
1771
1772 void PlotMillePede::ScanPedeParAbove(UInt_t iPar, float absMinInMuMm) {
1773 const TString mpPar(MpT() += Par(iPar));
1774 TString sel(Fixed(iPar, false) += AndL() += Abs(mpPar) += Form(">%f", absMinInMuMm));
1775 this->AddBasicSelection(sel);
1776 const TString titleAdd(this->TitleAdd());
1777 if (titleAdd.Length())
1778 std::cout << titleAdd << std::endl;
1779
1780 this->ScanSelection(sel);
1781
1782 }
1783
1784
1785 void PlotMillePede::DrawCheck() {
1786 fHistManager->Clear();
1787
1788 const TString allTreeNames[] = {OrgPosT(), MisPosT(), PosT(), MisParT(), ParT(), MpT()};
1789
1790 const unsigned int nTree = sizeof(allTreeNames) / sizeof(allTreeNames[0]);
1791 for (unsigned int i = 0; i < nTree; ++i) {
1792 TH1 *hId = this->CreateHist((PosT() += "Id - ") + allTreeNames[i] + "Id", "");
1793 TH1 *hObjId = this->CreateHist((PosT() += "ObjId - ") + allTreeNames[i] + "ObjId", "");
1794 fHistManager->AddHist(hId, i);
1795 fHistManager->AddHist(hObjId, i);
1796 }
1797
1798 fHistManager->Draw();
1799 }
1800
1801
1802 Int_t PlotMillePede::PrepareAdd(bool addPlots) {
1803 if (addPlots) {
1804 return fHistManager->GetNumLayers();
1805 } else {
1806 fHistManager->Clear();
1807 return 0;
1808 }
1809 }
1810
1811
1812 TString PlotMillePede::Unique(const char *name) const {
1813 if (!gROOT->FindObject(name))
1814 return name;
1815
1816 UInt_t i = 1;
1817 while (gROOT->FindObject(Form("%s_%u", name, i)))
1818 ++i;
1819
1820 return Form("%s_%u", name, i);
1821 }
1822
1823
1824 void PlotMillePede::AddBasicSelection(TString &sel) const {
1825
1826 if (fHieraLevel >= 0) {
1827 if (sel.IsNull())
1828 sel = HieraLev(fHieraLevel);
1829 else
1830 sel = Parenth(sel) += AndL() += HieraLev(fHieraLevel);
1831 }
1832
1833
1834 TString subDetSel;
1835 for (Int_t iSub = 0; iSub < fSubDetIds.GetSize(); ++iSub) {
1836 const TString newSubDetSel(Parenth(SubDetId()) += Form("==%d", fSubDetIds[iSub]));
1837 if (subDetSel.IsNull())
1838 subDetSel = newSubDetSel;
1839 else
1840 subDetSel = Parenth(subDetSel) += OrL() += newSubDetSel;
1841 }
1842 if (!subDetSel.IsNull()) {
1843 if (sel.IsNull())
1844 sel = Parenth(subDetSel);
1845 else
1846 sel = Parenth(sel) += AndL() += Parenth(subDetSel);
1847 }
1848
1849
1850 if (fAlignableTypeId >= 0) {
1851 const TString alignableTypeSel(Parenth(AlignableTypeId() += Form("==%d", fAlignableTypeId)));
1852 if (sel.IsNull())
1853 sel = alignableTypeSel;
1854 else
1855 sel = Parenth(sel) += AndL() += alignableTypeSel;
1856 }
1857
1858 if (fAdditionalSel.Length()) {
1859 if (sel.IsNull())
1860 sel = fAdditionalSel;
1861 else
1862 sel = Parenth(sel) += AndL() += fAdditionalSel;
1863 }
1864 }
1865
1866
1867 void PlotMillePede::SetMaxDev(Float_t maxDev) {
1868
1869 fMaxDevUp = TMath::Abs(maxDev);
1870 fMaxDevDown = -TMath::Abs(maxDev);
1871 }
1872
1873
1874 void PlotMillePede::SetMaxDev(Float_t maxDevDown, Float_t maxDevUp) {
1875
1876 if (maxDevUp < maxDevDown) {
1877 ::Error(
1878 "PlotMillePede::SetMaxDev", "Upper limit %f smaller than lower limit %f => Swap them!", maxDevUp, maxDevDown);
1879 fMaxDevUp = maxDevDown;
1880 fMaxDevDown = maxDevUp;
1881 } else {
1882 fMaxDevUp = maxDevUp;
1883 fMaxDevDown = maxDevDown;
1884 }
1885 }
1886
1887
1888 void PlotMillePede::SetSubDetId(Int_t subDetId) {
1889
1890 if (subDetId == -1) {
1891 fSubDetIds.Set(0);
1892 } else {
1893 fSubDetIds.Set(1, &subDetId);
1894 }
1895 }
1896
1897
1898 void PlotMillePede::AddSubDetId(Int_t subDetId) {
1899
1900 const Int_t last = fSubDetIds.GetSize();
1901 fSubDetIds.Set(last + 1);
1902 fSubDetIds[last] = subDetId;
1903 }
1904
1905
1906 void PlotMillePede::SetSubDetIds(Int_t id1, Int_t id2, Int_t id3, Int_t id4, Int_t id5) {
1907 this->SetSubDetId(id1);
1908 const Int_t ids[] = {id2, id3, id4, id5};
1909 for (unsigned int i = 0; i < sizeof(ids) / sizeof(ids[0]); ++i) {
1910 if (ids[i] > 0)
1911 this->AddSubDetId(ids[i]);
1912 }
1913 }
1914
1915
1916 TString PlotMillePede::TitleAdd() const {
1917
1918 TString result;
1919 for (Int_t i = 0; i < fSubDetIds.GetSize(); ++i) {
1920 if (result.Length())
1921 result += " + ";
1922 switch (fSubDetIds[i]) {
1923 case 1:
1924 result += "BPIX";
1925 break;
1926 case 2:
1927 result += "FPIX";
1928 break;
1929 case 3:
1930 result += "TIB";
1931 break;
1932 case 4:
1933 result += "TID";
1934 break;
1935 case 5:
1936 result += "TOB";
1937 break;
1938 case 6:
1939 result += "TEC";
1940 break;
1941 default:
1942
1943 result += "unknown subDet";
1944 }
1945 }
1946
1947 if (fAlignableTypeId >= 0) {
1948 if (result.Length())
1949 result += ", ";
1950
1951 result += this->AlignableObjIdString(fAlignableTypeId);
1952 }
1953
1954 if (fHieraLevel != 0) {
1955 if (result.Length())
1956 result += ", ";
1957 if (fHieraLevel < 0)
1958 result += Form("all hierar. levels");
1959 else
1960 result += Form("hier. level %d", fHieraLevel);
1961 }
1962
1963 if (fAdditionalSelTitle.Length()) {
1964 if (result.Length())
1965 result += ", ";
1966 result += fAdditionalSelTitle;
1967 }
1968
1969 if (fTitle.Length()) {
1970 if (result.Length())
1971 result.Prepend(fTitle + ", ");
1972 else
1973 result.Prepend(fTitle);
1974 }
1975
1976 if (result.Length())
1977 result.Prepend(": ");
1978
1979 return result;
1980 }
1981
1982
1983 TString PlotMillePede::AlignableObjIdString(Int_t objId) const {
1984 switch (objId) {
1985 case 1:
1986 return "DetUnit";
1987 case 2:
1988 return "Det";
1989
1990 case 5:
1991 return "BPIXLayer";
1992 case 6:
1993 return "BPIXHalfBarrel";
1994
1995 case 11:
1996 return "FPIXHalfDisk";
1997
1998 case 15:
1999 return "TIBString";
2000
2001 case 19:
2002 return "TIBHalfBarrel";
2003 case 20:
2004 return "TIBBarrel";
2005
2006 case 25:
2007 return "TIDEndcap";
2008
2009 case 29:
2010 return "TOBHalfBarrel";
2011 case 30:
2012 return "TOBBarrel";
2013
2014 case 36:
2015 return "TECEndcap";
2016 default:
2017 ::Error("PlotMillePede::AlignableObjIdString",
2018 "Missing implementation for ObjId %d, see "
2019 "Alignment/CommonAlignment/interface/StructureType.h",
2020 objId);
2021 return Form("alignable obj id %d", objId);
2022 }
2023 }
2024
2025
2026 Int_t PlotMillePede::SetHieraLevel(Int_t hieraLevel) {
2027 const Int_t oldLevel = fHieraLevel;
2028 fHieraLevel = hieraLevel;
2029 return oldLevel;
2030 }
2031
2032
2033 Int_t PlotMillePede::SetAlignableTypeId(
2034 Int_t alignableTypeId) {
2035 const Int_t oldTypeId = fAlignableTypeId;
2036 fAlignableTypeId = alignableTypeId;
2037 return oldTypeId;
2038 }
2039
2040
2041 void PlotMillePede::AddAdditionalSel(const char *selection) {
2042 if (!selection || selection[0] == '\0') {
2043 ::Warning("PlotMillePede::AddAdditionalSel", "Ignore empty selection.");
2044 } else {
2045 if (fAdditionalSel.Length())
2046 fAdditionalSel = Parenth(fAdditionalSel) += AndL();
2047
2048 if (fAdditionalSelTitle.Length())
2049 fAdditionalSelTitle += ", ";
2050 const TString sel(selection);
2051
2052 if (sel == "StripDoubleOr1D") {
2053 fAdditionalSel += "(Id&3)==0";
2054 fAdditionalSelTitle += "Double sided or 1D layer/ring";
2055 } else if (sel == "StripRphi") {
2056 fAdditionalSel += "(Id&3)==2";
2057 fAdditionalSelTitle += "R#phi";
2058 } else if (sel == "StripStereo") {
2059 fAdditionalSel += "(Id&3)==1";
2060 fAdditionalSelTitle += "Stereo";
2061
2062 } else if (sel == "NotStripDoubleOr1D") {
2063 fAdditionalSel += "(Id&3)!=0";
2064 fAdditionalSelTitle += "!(Double sided or 1D layer/ring)";
2065 } else if (sel == "NotStripRphi") {
2066 fAdditionalSel += "(Id&3)!=2";
2067 fAdditionalSelTitle += "!R#phi";
2068 } else if (sel == "NotStripStereo") {
2069 fAdditionalSel += "(Id&3)!=1";
2070 fAdditionalSelTitle += "!Stereo";
2071
2072 } else {
2073 fAdditionalSel += selection;
2074 fAdditionalSelTitle += selection;
2075 }
2076 }
2077 }
2078
2079
2080 void PlotMillePede::AddAdditionalSel(const TString &xyzrPhiNhit, Float_t min, Float_t max) {
2081 const TString oldTitle = fAdditionalSelTitle;
2082 if (xyzrPhiNhit == "Nhit") {
2083 this->AddAdditionalSel(HitsX() += Form(">=%f && ", min) + HitsX() += Form("<%f", max));
2084 } else {
2085 this->AddAdditionalSel(OrgPos(xyzrPhiNhit) += Form(">=%f && ", min) + OrgPos(xyzrPhiNhit) += Form("<%f", max));
2086 }
2087
2088 fAdditionalSelTitle = oldTitle;
2089 if (fAdditionalSelTitle.Length())
2090 fAdditionalSelTitle += ", ";
2091 fAdditionalSelTitle += Form("%g #leq %s < %g", min, xyzrPhiNhit.Data(), max);
2092 }
2093
2094
2095 TString PlotMillePede::FinalMisAlignment(UInt_t iPar) const {
2096 return (fUseDiff ? DiffPar(ParT(), MisParT(), iPar) : Parenth(ParT() += Par(iPar)));
2097 }
2098
2099
2100 void PlotMillePede::CopyAddBinning(TString &name, const TH1 *h) const {
2101 if (!h)
2102 return;
2103
2104 name += Form("(%d,%f,%f", h->GetNbinsX(), h->GetXaxis()->GetXmin(), h->GetXaxis()->GetXmax());
2105 if (h->GetDimension() > 1) {
2106 name += Form(",%d,%f,%f", h->GetNbinsY(), h->GetYaxis()->GetXmin(), h->GetYaxis()->GetXmax());
2107 }
2108 if (h->GetDimension() > 2) {
2109 name += Form(",%d,%f,%f", h->GetNbinsZ(), h->GetZaxis()->GetXmin(), h->GetZaxis()->GetXmax());
2110 }
2111
2112 name += ')';
2113 }
2114
2115
2116 void PlotMillePede::SetOutName(const TString &name) { fHistManager->SetCanvasName(name); }