|
|
Go to the documentation of this file.
13 #include "TMatrixDSymEigen.h"
14 #include "TDecompSVD.h"
63 unsigned int id = it.first;
67 sprintf(
buf,
"%u: m distribution",
id);
68 s.m_dist =
new TH1D(
buf,
";u or v (mm)", 100, -25., 25.);
70 sprintf(
buf,
"%u: R distribution",
id);
71 s.R_dist =
new TH1D(
buf,
";R (mm)", 500, -0.5, 0.5);
75 s.coefHist.push_back(
new TH1D(
buf,
";coefficient", 100, -2., +2.));
78 TGraph *
g =
new TGraph();
80 g->SetTitle(
";coefficient;residual (mm)");
81 s.resVsCoef.push_back(
g);
95 printf(
"\n>> JanAlignmentAlgorithm::Feed\n");
100 double hax = trackFit.
ax;
101 double hay = trackFit.
ay;
116 set<unsigned int> rpSet;
120 const unsigned int rpDecId = 100 * detId.
arm() + 10 * detId.
station() + detId.
rp();
121 rpSet.insert(rpDecId);
129 unsigned int id =
hit->
id;
132 const auto &dirData =
d.getDirectionData(
hit->dirIdx);
134 A(
j, 0) =
hit->
z * dirData.dx;
135 A(
j, 1) = dirData.dx;
136 A(
j, 2) =
hit->
z * dirData.dy;
137 A(
j, 3) = dirData.dy;
139 m(
j) =
hit->position + dirData.s - (
hit->
z -
d.z) * dirData.dz;
141 Vi(
j,
j) = 1. /
hit->sigma /
hit->sigma;
143 double C = dirData.dx,
S = dirData.dy;
145 double hx = hax *
hit->
z + hbx;
146 double hy = hay *
hit->
z + hby;
147 double R =
m(
j) - (hx *
C + hy *
S);
164 Ga[
i][
j][matrixIndex] = -1.;
168 Ga[
i][
j][matrixIndex] = -1.;
172 Ga[
i][
j][matrixIndex] = hax *
C + hay *
S;
176 Ga[
i][
j][matrixIndex] = (hax *
hit->
z + hbx -
d.sx) * (-
S) + (hay *
hit->
z + hby -
d.sy) *
C;
181 double c = Ga[
i][
j][matrixIndex];
183 s.coefHist[
i]->Fill(
c);
184 s.resVsCoef[
i]->SetPoint(
s.resVsCoef[
i]->GetN(),
c,
R);
187 map<set<unsigned int>,
ScatterPlot>::iterator
hit =
s.resVsCoefRot_perRPSet.find(rpSet);
188 if (
hit ==
s.resVsCoefRot_perRPSet.end()) {
191 sp.
h =
new TH2D(
"",
"", 40, -20., +20., 60, -0.15, +0.15);
192 hit =
s.resVsCoefRot_perRPSet.insert(pair<set<unsigned int>,
ScatterPlot>(rpSet, sp)).first;
194 hit->second.g->SetPoint(
hit->second.g->GetN(),
c,
R);
195 hit->second.h->Fill(
c,
R);
202 TMatrixD AT(TMatrixD::kTransposed,
A);
203 TMatrixD ATViA(4, 4);
205 TMatrixD ATViAI(ATViA);
206 ATViAI = ATViA.Invert();
208 sigma -= Vi *
A * ATViAI * AT * Vi;
214 GaT[
i].Transpose(Ga[
i]);
223 if (
Mc[
i].GetNrows() < 1)
232 if (
Sc[
i][
j].GetNrows() < 1 ||
Sc[
i][
j].GetNcols() < 1)
235 Sc[
i][
j] += GaT[
i] * sigma * Ga[
j];
240 printf(
"* checking normalized residuals, selection.size = %u\n",
selection.size());
244 printf(
"- class %u\n",
i);
245 TVectorD
t(
Mc[
i].GetNrows());
246 for (
int j = 0;
j <
t.GetNrows();
j++)
253 tt = sigma * Ga[
i] *
t;
255 double ttn =
sqrt(
tt.Norm2Sqr());
256 printf(
"|tt| = %E\n", ttn);
270 printf(
"\n>> JanAlignmentAlgorithm::Analyze\n");
273 unsigned int dim = 0;
275 dim +=
Mc[
i].GetNrows();
279 printf(
"\tfull dimension: %u\n", dim);
292 S.ResizeTo(dim, dim);
293 unsigned int r_offset = 0, c_offset = 0;
296 unsigned int r_size = 0, c_size = 0;
298 r_size =
Sc[
i][
j].GetNrows();
299 c_size =
Sc[
i][
j].GetNcols();
301 if (r_size < 1 || c_size < 1)
304 TMatrixDSub(
S, r_offset, r_offset + r_size - 1, c_offset, c_offset + c_size - 1) =
Sc[
i][
j];
312 double maxDiff = 0., maxElem = 0.;
313 for (
unsigned int i = 0;
i < dim;
i++) {
314 for (
unsigned int j = 0;
j < dim;
j++) {
318 if (
S[
i][
j] > maxElem)
323 printf(
"\n* S matrix:\n\tdimension = %i\n\tmaximum asymmetry: %E\t(ratio to maximum element %E)\n",
330 TMatrixDSym S_sym(dim);
331 for (
unsigned int j = 0;
j < dim;
j++) {
332 for (
unsigned int i = 0;
i < dim;
i++) {
333 S_sym[
i][
j] =
S[
i][
j];
338 TMatrixDSymEigen S_eig(S_sym);
339 const TVectorD &S_eigVal_temp = S_eig.GetEigenValues();
340 S_eigVal.ResizeTo(S_eigVal_temp.GetNrows());
342 const TMatrixD &S_eigVec_temp = S_eig.GetEigenVectors();
358 printf(
"\n* S singular modes\n | ");
364 printf(
"---------- ");
367 for (
unsigned int j = 0;
j < dim;
j++) {
375 printf(
"\n* S has no singular modes\n");
382 map<unsigned int, AlignmentResult> &
results,
385 printf(
">> JanAlignmentAlgorithm::Solve\n");
390 unsigned int dim =
S.GetNrows();
396 const TVectorD &
cv =
constraints[
i].coef.find(quantityClass)->second;
397 for (
int k = 0;
k <
cv.GetNrows();
k++) {
406 printf(
"\n* constraint matrix\n");
413 for (
int j = 0;
j <
S.GetNrows();
j++)
421 for (
unsigned int j = 0;
j < dim;
j++) {
422 for (
unsigned int i = 0;
i < dim;
i++) {
429 for (
unsigned int j = 0;
j < dim;
j++) {
430 CS[
j][dim +
i] = CS[dim +
i][
j] =
C(
j,
i);
431 CS2[
j][dim +
i] = CS2[dim +
i][
j] = C2(
j,
i);
436 TMatrixDSymEigen CS_eig(CS);
437 TVectorD CS_eigVal = CS_eig.GetEigenValues();
438 TMatrixD CS_eigVec = CS_eig.GetEigenVectors();
442 printf(
"\n* eigen values of CS and S matrices (events = %u)\n",
events);
443 printf(
" # CS norm. CS S norm. S\n");
446 unsigned int singularModeCount = 0;
447 vector<unsigned int> weakModeIdx;
448 for (
int i = 0;
i < CS_eigVal.GetNrows();
i++) {
449 const double CS_nev = CS_eigVal[
i] /
events;
455 printf(
"%4i%+12.2E%+12.2E",
i, CS_eigVal[
i], CS_nev);
461 weakModeIdx.push_back(
i);
470 printf(
"%+12.2E%+12.2E",
S_eigVal[
i], S_nev);
483 if (!weakModeIdx.empty()) {
485 unsigned int first = 0;
487 while (
first < weakModeIdx.size()) {
489 if (
last >= weakModeIdx.size())
490 last = weakModeIdx.size();
492 printf(
"\n* CS weak modes\n | ");
494 printf(
"%+10.3E ", CS_eigVal[weakModeIdx[
i]]);
498 printf(
"---------- ");
506 double v = fabs(CS_eigVec(weakModeIdx[
i],
j));
516 double v = CS_eigVec(weakModeIdx[
i],
j);
517 if (fabs(
v) / maxs[
i -
first] > 1E-3)
518 printf(
"%+10.3E ",
v);
528 printf(
"\n* CS has no weak modes\n");
533 if (E.GetNcols() ==
C.GetNcols()) {
534 TMatrixD CTE(
C, TMatrixD::kTransposeMult, E);
535 print(CTE,
"* CTE matrix:");
536 const double &det = CTE.Determinant();
538 "\n* det(CTE) = %E, max(CTE) = %E, det(CTE)/max(CTE) = %E\n\tmax(C) = %E, max(E) = %E, "
539 "det(CTE)/max(C)/max(E) = %E\n",
545 det /
C.Max() / E.Max());
547 printf(
">> JanAlignmentAlgorithm::Solve > WARNING: C matrix has %u, while E matrix %u columns.\n",
555 LogError(
"PPS") <<
"\n>> JanAlignmentAlgorithm::Solve > ERROR: There are " << singularModeCount
556 <<
" singular modes in CS matrix. Stopping.";
562 for (
unsigned int i = 0;
i < dim;
i++)
568 TMatrixD CSI(TMatrixD::kInverted, CS);
569 TMatrixD CS2I(TMatrixD::kInverted, CS2);
580 EM2 = CS2I * S0 * CS2I;
582 TMatrixD EMdiff(EM2 - EM);
585 double max1 = -1., max2 = -1.,
maxDiff = -1.;
586 for (
int i = 0;
i < EMdiff.GetNrows();
i++) {
587 for (
int j = 0;
j < EMdiff.GetNcols();
j++) {
591 if (max1 < fabs(EM(
i,
j)))
592 max1 = fabs(EM(
i,
j));
594 if (max2 < fabs(EM2(
i,
j)))
595 max2 = fabs(EM2(
i,
j));
599 printf(
"EM max = %E, EM2 max = %E, EM2 - EM max = %E\n", max1, max2,
maxDiff);
603 TMatrixD &
U = CS_eigVec;
604 TMatrixD UT(TMatrixD::kTransposed,
U);
615 for (
int i = 0;
i < EMEi.GetNrows();
i++) {
616 for (
int j = 0;
j < EMEi.GetNcols();
j++) {
617 if (
max < EMEi(
i,
j))
622 printf(
"max = %E\n",
max);
627 printf(
"\n* Lambda (from the contribution of singular modes to MV)\n");
629 printf(
"\t%u (%25s)\t%+10.1E +- %10.1E\n",
633 sqrt(EM[
i + dim][
i + dim]) * 1E3);
655 double e =
sqrt(EM[fi][fi]);
691 CS_eigVal.Write(
"CS_eigen_values");
692 CS_eigVec.Write(
"CS_eigen_vectors");
722 for (map<unsigned int, DetStat>::iterator it =
statistics.begin(); it !=
statistics.end(); ++it) {
724 sprintf(
buf,
"%u", it->first);
725 gDirectory =
dir->mkdir(
buf);
727 it->second.m_dist->Write();
728 it->second.R_dist->Write();
731 it->second.coefHist[
c]->Write();
732 it->second.resVsCoef[
c]->Write();
735 gDirectory = gDirectory->mkdir(
"R vs. rot. coef, per RP set");
736 TCanvas *
c =
new TCanvas;
737 c->SetName(
"R vs. rot. coef, overlapped");
738 TH2D *
frame =
new TH2D(
"frame",
"frame", 100, -20., +20., 100, -0.15, +0.15);
740 unsigned int idx = 0;
741 for (
map<set<unsigned int>,
ScatterPlot>::iterator iit = it->second.resVsCoefRot_perRPSet.begin();
742 iit != it->second.resVsCoefRot_perRPSet.end();
746 for (set<unsigned int>::iterator sit = iit->first.begin(); sit != iit->first.end(); ++sit) {
748 sprintf(
buf,
"%u", *sit);
758 iit->second.g->SetTitle(
";rotation coefficient (mm);residual (mm)");
759 iit->second.g->SetMarkerColor(
idx + 1);
760 iit->second.g->SetName(
label.c_str());
761 iit->second.g->Draw((
idx == 0) ?
"p" :
"p");
762 iit->second.g->Write();
764 iit->second.h->SetName((
label +
" (hist)").c_str());
765 iit->second.h->SetTitle(
";rotation coefficient (mm);residual (mm)");
766 iit->second.h->Write();
769 gDirectory->cd(
"..");
JanAlignmentAlgorithm()
dummy constructor (not to be used)
double z0
a characteristic z in mm
signed int getQuantityIndex(QuantityClass cl, unsigned int detId) const
returns measurement index (if non-existent, returns -1)
The manager class for TOTEM RP geometry.
Represents an alignment task.
detector shifts in first readout direction
void feed(const HitCollection &, const LocalTrackFit &) override
process one track
void analyze() override
analyzes the data collected
std::string quantityClassTag(QuantityClass) const
returns a string tag for the given quantity class
~JanAlignmentAlgorithm() override
float maxDiff(float one, float two, float three, float four)
detector rotations around z
void begin(const CTPPSGeometry *geometryReal, const CTPPSGeometry *geometryMisaligned) override
prepare for processing
unsigned int events
event count
TMatrixD S_eigVec
matrix of S eigenvectors
A structure to hold relevant geometrical information about one detector/sensor.
unsigned int solve(const std::vector< AlignmentConstraint > &, std::map< unsigned int, AlignmentResult > &results, TDirectory *dir) override
double weakLimit
normalized eigen value below which the (CS) eigen vectors are considered as weak
Abstract parent for all (track-based) alignment algorithms.
a scatter plot, with graph and histogram representations
bool stopOnSingularModes
whether to stop when singular modes are identified
signed int getMeasurementIndex(QuantityClass cl, unsigned int detId, unsigned int dirIdx) const
returns measurement index (if non-existent, returns -1)
const std::map< unsigned int, DetGeometry > & getSensorMap() const
structure holding statistical information for one detector
void print(TMatrixD &m, const char *label=nullptr, bool mathematicaFormat=false)
AlignmentTask * task
the tasked to be completed
Base class for CTPPS detector IDs.
unsigned int quantitiesOfClass(QuantityClass) const
returns the number of quantities of the given class
Log< level::Error, false > LogError
const DetGeometry & get(unsigned int id) const
retrieves sensor geometry
Result of CTPPS track-based alignment.
TVectorD S_eigVal
eigen values of the S matrix
double bx
intercepts in mm
std::map< unsigned int, DetStat > statistics
statistical data collection
std::vector< Hit > HitCollection
Local (linear) track description (or a fit result). Uses global reference system.
T getParameter(std::string const &) const
std::vector< QuantityClass > quantityClasses
list of quantity classes to be optimized
double singularLimit
eigenvalues in (-singularLimit, singularLimit) are treated as singular
AlignmentGeometry geometry
the geometry for this task
void end() override
cleans up after processing
bool buildDiagnosticPlots
flag whether to build statistical plots
double z0
the point where intercepts are measured, in mm
void saveDiagnostics(TDirectory *) override
saves diagnostic histograms/plots
std::vector< SingularMode > singularModes
a list of the singular modes of the S matrix
detector shifts in second readout direction
ParameterSet const & getParameterSet(std::string const &) const
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
unsigned int getNumberOfDetectors() const
returns the number of detectors in the collection