|
|
Go to the documentation of this file.
23 #include <unordered_set>
27 #include "TDecompLU.h"
41 return new TH1D(
name,
";residual (mm)", 2000, -3., +3.);
47 TGraph *
g =
new TGraph();
48 g->SetName(
name.c_str());
49 g->SetTitle(
title.c_str());
56 chisqn_lin_fitted =
new TH1D(
"chi^2 norm, lin, fitted",
";#chi^{2}/ndf;", 5000, 0., 500.);
57 chisqn_lin_selected =
new TH1D(
"chi^2 norm, lin, selected",
";#chi^{2}/ndf;", 5000, 0., 500.);
58 chisqn_log_fitted =
new TH1D(
"chi^2 norm, log, fitted",
";log_{10}(#chi^{2}/ndf);", 700, -1., 6.);
59 chisqn_log_selected =
new TH1D(
"chi^2 norm, log, selected",
";log_{10}(#chi^{2}/ndf);", 700, -1., 6.);
61 fitAxVsAyGraph_fitted =
newGraph(
"ax vs. ay, fitted",
";a_{x} (rad);a_{y} (rad)");
62 fitAxVsAyGraph_selected =
newGraph(
"ax vs. ay, selected",
";a_{x} (rad);a_{y} (rad)");
63 fitBxVsByGraph_fitted =
newGraph(
"bx vs. by, fitted",
";b_{x} (mm);b_{y} (mm)");
64 fitBxVsByGraph_selected =
newGraph(
"bx vs. by, selected",
";b_{x} (mm);b_{y} (mm)");
83 chisqn_lin_fitted->Write();
84 chisqn_lin_selected->Write();
85 chisqn_log_fitted->Write();
86 chisqn_log_selected->Write();
88 fitAxVsAyGraph_fitted->Write();
89 fitAxVsAyGraph_selected->Write();
90 fitBxVsByGraph_fitted->Write();
91 fitBxVsByGraph_selected->Write();
98 :
verbosity(ps.getUntrackedParameter<unsigned
int>(
"verbosity", 0)),
102 z0(ps.getParameter<double>(
"z0")),
112 maxTrackAx(ps.getParameter<double>(
"maxTrackAx")),
113 maxTrackAy(ps.getParameter<double>(
"maxTrackAy")),
150 vector<string> alNames(ps.
getParameter<vector<string> >(
"algorithms"));
151 for (
unsigned int i = 0;
i < alNames.size();
i++) {
154 if (alNames[
i] ==
"Ideal") {
157 }
else if (alNames[
i] ==
"Jan") {
165 throw cms::Exception(
"PPS") <<
"Unknown alignment algorithm `" << alNames[
i] <<
"'.";
170 if (ct.compare(
"fixedDetectors") == 0)
172 else if (ct.compare(
"standard") == 0)
175 throw cms::Exception(
"PPS") <<
"Unknown constraints type `" << ct <<
"'.";
178 string aars_str = ps.
getParameter<
string>(
"additionalAcceptedRPSets");
180 size_t idx_b = 0, idx_e = string::npos;
181 while (idx_b != string::npos) {
183 idx_e = aars_str.find(
';', idx_b);
184 size_t len = (idx_e == string::npos) ? string::npos : idx_e - idx_b;
185 string block = aars_str.substr(idx_b, len);
188 if (!
block.empty()) {
189 set<unsigned int> rpSet;
192 size_t bi_b = 0, bi_e = string::npos;
193 while (bi_b != string::npos) {
194 bi_e =
block.find(
',', bi_b);
195 size_t bit_len = (bi_e == string::npos) ? string::npos : bi_e - bi_b;
196 const string &bit =
block.substr(bi_b, bit_len);
198 unsigned int rp = atoi(bit.c_str());
201 bi_b = (bi_e == string::npos) ? string::npos : bi_e + 1;
208 idx_b = (idx_e == string::npos) ? string::npos : idx_e + 1;
218 for (vector<AlignmentAlgorithm *>::iterator it =
algorithms.begin(); it !=
algorithms.end(); ++it)
240 delete it->second.total_fitted;
241 delete it->second.total_selected;
242 delete it->second.selected_vs_chiSq;
243 for (
auto sit = it->second.perRPSet_fitted.begin(); sit != it->second.perRPSet_fitted.end(); ++sit)
245 for (
auto sit = it->second.perRPSet_selected.begin(); sit != it->second.perRPSet_selected.end(); ++sit)
297 printf(
"\n---------- StraightTrackAlignment::ProcessEvent (%u:%llu)\n", eventId.
run(), eventId.
event());
304 for (
const auto &
pv : uvPatternsStrip) {
306 const unsigned int rpDecId = detId.
arm() * 100 + detId.station() * 10 + detId.rp();
313 unsigned int n_U = 0, n_V = 0;
314 unsigned int idx_U = 0, idx_V = 0;
315 for (
unsigned int pi = 0;
pi <
pv.size();
pi++) {
318 switch (
pattern.projection()) {
334 if (n_U != 1 || n_V != 1)
338 if (!
pv[idx_U].fittable() || !
pv[idx_V].fittable())
342 for (
const auto &
pattern : {
pv[idx_U],
pv[idx_V]}) {
343 for (
const auto &hitsDetSet :
pattern.hits()) {
349 for (
auto &
hit : hitsDetSet)
356 for (
const auto &
pv : hitsDiamond) {
359 const unsigned int rpDecId = detId.
arm() * 100 + detId.station() * 10 + detId.rp();
371 for (
const auto &
h :
pv)
376 map<unsigned int, unsigned int> pixelPlaneMultiplicity;
377 for (
const auto &
pv : hitsPixel)
378 pixelPlaneMultiplicity[
pv.detId()] +=
pv.size();
380 for (
const auto &
pv : hitsPixel) {
383 const unsigned int rpDecId = detId.
arm() * 100 + detId.station() * 10 + detId.rp();
394 if (pixelPlaneMultiplicity[
pv.detId()] > 1)
397 for (
const auto &
h :
pv) {
400 const double dz2 = dg.getDirectionData(2).dz;
401 const double z = dg.z +
h.point().x() * dz1 +
h.point().y() * dz2;
403 double x_unc =
sqrt(
h.error().xx());
404 double y_unc =
sqrt(
h.error().yy());
419 for (
const auto &
pv : tracksPixel) {
421 const unsigned int rpDecId =
rpId.arm() * 100 +
rpId.station() * 10 +
rpId.rp();
428 unsigned int n_valid_tracks = 0;
429 for (
const auto &tr :
pv) {
434 if (n_valid_tracks > 1)
438 for (
const auto &tr :
pv) {
443 for (
const auto &hv : tr.hits()) {
450 for (
const auto &
h : hv) {
452 if (!
h.isUsedForFit())
457 const double dz2 = dg.getDirectionData(2).dz;
458 const double z = dg.z +
h.point().x() * dz1 +
h.point().y() * dz2;
460 double x_unc =
sqrt(
h.error().xx());
461 double y_unc =
sqrt(
h.error().yy());
484 set<unsigned int> selectedRPs;
487 const unsigned int decRPId = detId.
arm() * 100 + detId.
station() * 10 + detId.
rp();
488 selectedRPs.insert(decRPId);
496 bool top =
false, bottom =
false, horizontal =
false;
497 unordered_set<unsigned int>
units;
498 for (
const auto &rp : selectedRPs) {
499 unsigned int rpIdx = rp % 10;
500 unsigned int stId = rp / 10;
501 unsigned int unitId = stId * 10;
505 if (rpIdx == 0 || rpIdx == 4)
507 if (rpIdx == 1 || rpIdx == 5)
509 if (rpIdx == 2 || rpIdx == 3)
512 units.insert(unitId);
515 bool overlap = (top && horizontal) || (bottom && horizontal);
517 bool rp_set_accepted =
true;
521 rp_set_accepted =
false;
525 rp_set_accepted =
false;
528 rp_set_accepted =
false;
531 rp_set_accepted =
false;
536 rp_set_accepted =
true;
539 printf(
"* rp set accepted: %u\n", rp_set_accepted);
541 bool selected = rp_set_accepted;
554 printf(
"* event selected: %u\n", selected);
574 throw "StraightTrackAlignment: Number of tracks processed reached maximum";
580 const set<unsigned int> &selectedRPs,
582 bool trackSelected) {
616 it =
rpSetPlots.insert({selectedRPs, RPSetPlots(setToString(selectedRPs))}).first;
618 it->second.chisqn_lin_fitted->Fill(trackFit.
chiSqPerNdf());
619 it->second.chisqn_log_fitted->Fill(log10(trackFit.
chiSqPerNdf()));
620 it->second.fitAxVsAyGraph_fitted->SetPoint(it->second.fitAxVsAyGraph_fitted->GetN(), trackFit.
ax, trackFit.
ay);
621 it->second.fitBxVsByGraph_fitted->SetPoint(it->second.fitBxVsByGraph_fitted->GetN(), trackFit.
bx, trackFit.
by);
624 it->second.chisqn_lin_selected->Fill(trackFit.
chiSqPerNdf());
625 it->second.chisqn_log_selected->Fill(log10(trackFit.
chiSqPerNdf()));
626 it->second.fitAxVsAyGraph_selected->SetPoint(it->second.fitAxVsAyGraph_selected->GetN(), trackFit.
ax, trackFit.
ay);
627 it->second.fitBxVsByGraph_selected->SetPoint(it->second.fitBxVsByGraph_selected->GetN(), trackFit.
bx, trackFit.
by);
631 unsigned int id =
hit.
id;
634 const auto dirData =
geom.getDirectionData(
hit.dirIdx);
636 double m =
hit.position + dirData.s - (
hit.
z -
geom.z) * dirData.dz;
637 double x = trackFit.
ax *
hit.
z + trackFit.
bx;
638 double y = trackFit.
ay *
hit.
z + trackFit.
by;
639 double f =
x * dirData.dx +
y * dirData.dy;
646 sprintf(
buf,
"%u: total_fitted",
id);
648 sprintf(
buf,
"%u: total_selected",
id);
650 it->second.selected_vs_chiSq =
new TGraph();
651 sprintf(
buf,
"%u: selected_vs_chiSq",
id);
652 it->second.selected_vs_chiSq->SetName(
buf);
655 it->second.total_fitted->Fill(
R);
658 it->second.total_selected->Fill(
R);
659 it->second.selected_vs_chiSq->SetPoint(it->second.selected_vs_chiSq->GetN(), trackFit.
chiSqPerNdf(),
R);
662 auto sit = it->second.perRPSet_fitted.find(selectedRPs);
663 if (sit == it->second.perRPSet_fitted.end()) {
665 sprintf(
buf,
"%u: ",
id);
669 it->second.perRPSet_fitted.insert(pair<set<unsigned int>, TH1D *>(selectedRPs,
newResiduaHist(
label.c_str())))
673 sit->second->Fill(
R);
676 sit = it->second.perRPSet_selected.find(selectedRPs);
677 if (sit == it->second.perRPSet_selected.end()) {
679 sprintf(
buf,
"%u: ",
id);
682 sit = it->second.perRPSet_selected
687 sit->second->Fill(
R);
713 printf(
"----------------------------------------------------------------------------------------------------\n");
714 printf(
"\n>> StraightTrackAlignment::Finish\n");
719 printf(
"\n* events per RP set:\n");
720 printf(
"%30s %10s%10s\n",
"set of RPs",
"fitted",
"selected");
728 printf(
"%30s :%10lu%10lu\n",
label.c_str(), it->second,
sv);
732 printf(
"\n* hits per plane:\n");
736 printf(
" : %u\n", it.second);
757 printf(
"\n>> StraightTrackAlignment::Finish > %lu constraints built\n",
constraints.size());
763 vector<map<unsigned int, AlignmentResult> >
results;
765 TDirectory *
dir =
nullptr;
774 <<
"' algorithm has failed (return value " <<
rf <<
").";
779 printf(
"\n>> StraightTrackAlignment::Finish > Print\n");
785 signed int prevRPId = -1;
789 if (
rpId != prevRPId)
797 for (
unsigned int a = 0;
a <
results.size();
a++) {
798 const auto it =
results[
a].find(dit.first);
801 printf(
"%19s",
"----");
803 printf(
"%8s",
"----");
813 const auto &ac = it->second;
814 double v = 0.,
e = 0.;
835 printf(
"%+8.1f ± %7.1f",
v * 1E3,
e * 1E3);
837 printf(
"%+8.1f",
v * 1E3);
856 for (
unsigned int a = 0;
a <
results.size();
a++) {
861 const auto dir1 =
d.getDirectionData(1);
862 const auto dir2 =
d.getDirectionData(2);
865 const double sh_x = (
dir2.dy *
p.second.getShR1() -
dir1.dy *
p.second.getShR2()) / det;
866 const double sh_y = (-
dir2.dx *
p.second.getShR1() +
dir1.dx *
p.second.getShR2()) / det;
868 const double sh_x_e =
869 sqrt(
pow(
dir2.dy / det *
p.second.getShR1Unc(), 2) +
pow(
dir1.dy / det *
p.second.getShR2Unc(), 2));
870 const double sh_y_e =
871 sqrt(
pow(
dir2.dx / det *
p.second.getShR1Unc(), 2) +
pow(
dir1.dx / det *
p.second.getShR2Unc(), 2));
878 p.second.getShZUnc(),
884 p.second.getRotZUnc());
904 cumulativeAlignments.
addCorrections(convertedAlignments,
false,
true,
true);
914 const bool equalWeights =
false;
916 cumulativeAlignments, expandedAlignments, factoredAlignments,
task.
geometry, equalWeights,
verbosity);
950 unsigned int N =
s.size();
957 for (set<unsigned int>::iterator it =
s.begin(); it !=
s.end(); ++it, ++
i) {
958 sprintf(
buf,
"%u", *it);
970 for (
unsigned int i = 0;
i <
N;
i++)
977 printf(
"═════════════════════════╬");
979 for (
unsigned int a = 0;
a <
results.size();
a++) {
995 unsigned int size = 0;
996 for (
unsigned int a = 0;
a <
results.size();
a++)
1001 unsigned int space = (
size -
tag.size()) / 2;
1003 printf(
"%s",
tag.c_str());
1016 for (
unsigned int a = 0;
a <
results.size();
a++) {
1039 TDirectory *commonDir =
df->mkdir(
"common");
1040 gDirectory = commonDir;
1055 gDirectory = commonDir->mkdir(
"plots global");
1058 TDirectory *ppsDir = commonDir->mkdir(
"plots per RP set");
1060 gDirectory = ppsDir->mkdir(
setToString(it->first).c_str());
1065 TDirectory *resDir = commonDir->mkdir(
"residuals");
1066 for (map<unsigned int, ResiduaHistogramSet>::const_iterator it =
residuaHistograms.begin();
1070 sprintf(
buf,
"%u", it->first);
1071 gDirectory = resDir->mkdir(
buf);
1072 it->second.total_fitted->Write();
1073 it->second.total_selected->Write();
1074 it->second.selected_vs_chiSq->Write();
1084 gDirectory = gDirectory->mkdir(
"selected per RP set");
1085 TCanvas *
c =
new TCanvas;
1086 c->SetName(
"alltogether");
1087 unsigned int idx = 0;
1088 for (
map<set<unsigned int>, TH1D *>::const_iterator sit = it->second.perRPSet_selected.begin();
1089 sit != it->second.perRPSet_selected.end();
1091 sit->second->SetLineColor(
idx + 1);
1092 sit->second->Draw((
idx == 0) ?
"" :
"same");
1093 sit->second->Write();
1100 for (vector<AlignmentAlgorithm *>::const_iterator it =
algorithms.begin(); it !=
algorithms.end(); ++it) {
1101 TDirectory *algDir =
df->mkdir((*it)->getName().c_str());
1102 (*it)->saveDiagnostics(algDir);
signed int eventsFitted
counter of processed tracks
double dz
x, y and z components of the direction unit vector in global coordinates
T const * product() const
TH1D * chisqn_log_selected
signed int maxEvents
stops after this event number has been reached
bool requireAtLeast3PotsInOverlap
if a track goes through overlap, select it only if it leaves signal in at least 3 pots
detector shifts in first readout direction
const DirectionData & getDirectionData(unsigned int idx) const
std::map< std::set< unsigned int >, unsigned long > fittedTracksPerRPSet
counter of fitted tracks in a certain RP set
TGraph * fitAxVsAyGraph_fitted
plots ax vs. ay
Calculates the ideal result of the StraightTrackAlignment.
std::map< unsigned int, unsigned int > selectedHitsPerPlane
counter of selected hits per plane
signed int eventsTotal
counter of events
TFile * taskDataFile
the file with task data
signed int ndf
the number of degrees of freedom
double chiSqPerNdf() const
std::string quantityClassTag(QuantityClass) const
returns a string tag for the given quantity class
std::string diagnosticsFile
file name for some event selection statistics
bool saveIntermediateResults
whether itermediate results (S, CS matrices) of alignments shall be saved
std::map< unsigned int, ResiduaHistogramSet > residuaHistograms
residua histograms
std::string factoredFileNamePrefix
file name prefix for cumulative factored result files
detector rotations around z
TGraph * fitBxVsByGraph_selected
TH1D * chisqn_lin_fitted
normalised chi^2 histograms for all/selected tracks, in linear/logarithmic scale
TH1D * fitBxHist_selected
fit bx histograms for all/selected tracks
void printAlgorithmsLine(const std::vector< std::map< unsigned int, AlignmentResult > > &)
signed int eventsSelected
counter of processed tracks
LocalTrackFitter fitter
track fitter
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
void addCorrections(const CTPPSRPAlignmentCorrectionsData &, bool sumErrors=true, bool addSh=true, bool addRot=true)
adds (merges) corrections on top of the current values
void buildIndexMaps()
builds "mapMatrixIndeces" from "geometry"
std::string taskDataFileName
the name task data file
virtual void finish()
performs analyses and fill results variable
RPSetPlots globalPlots
global (all RP sets) chi^2 histograms
void buildFixedDetectorsConstraints(std::vector< AlignmentConstraint > &) const
builds a set of fixed-detector constraints
TString units(TString variable, Char_t axis)
TH1D * fitNdfHist_selected
fit num. of degrees of freedom histograms for all/selected tracks
A structure to hold relevant geometrical information about one detector/sensor.
CTPPSRPAlignmentCorrectionsData initialAlignments
(real geometry) alignments before this alignment iteration
AlignmentTask task
the alignment task to be solved
bool fit(HitCollection &, const AlignmentGeometry &, LocalTrackFit &) const
unsigned int verbosity
verbosity level
Abstract parent for all (track-based) alignment algorithms.
TString getName(TString structure, int layer, TString geometry)
void updateDiagnosticHistograms(const HitCollection &selection, const std::set< unsigned int > &selectedRPs, const LocalTrackFit &trackFit, bool trackSelected)
fills diagnostic (chi^2, residua, ...) histograms
Jan's alignment algorithm.
bool requireOverlap
if true, only track through vertical-horizontal overlap are seleceted
static void writeToXML(const CTPPSRPAlignmentCorrectionsDataSequence &seq, const std::string &fileName, bool precise=false, bool wrErrors=true, bool wrSh_xy=true, bool wrSh_z=false, bool wrRot_xy=false, bool wrRot_z=true)
writes sequence of alignment corrections into a single XML file
TH1D * fitPHist_selected
fit p-value histograms for all/selected tracks
bool cutOnChiSqPerNdf
whether to cut on chi^2/ndf
std::vector< std::set< unsigned int > > additionalAcceptedRPSets
list of RP sets accepted irrespective of the other "require" settings
const std::map< unsigned int, DetGeometry > & getSensorMap() const
A linear pattern in U or V projection. The intercept b is taken at the middle of a RP: (geometry->Get...
StraightTrackAlignment(const edm::ParameterSet &)
std::vector< unsigned int > excludePlanes
list of planes to be excluded from processing
bool removeImpossible
remove events with impossible signatures (i.e. simultaneously top and bottom)
map: detector id --> residua histogram
bool buildDiagnosticPlots
whether to build and save diagnostic plots
TH1D * newResiduaHist(const char *name)
creates a new residua histogram
virtual void processEvent(const edm::EventID &eventId, const edm::DetSetVector< TotemRPUVPattern > &uvPatternsStrip, const edm::DetSetVector< CTPPSDiamondRecHit > &hitsDiamond, const edm::DetSetVector< CTPPSPixelRecHit > &hitsPixel, const edm::DetSetVector< CTPPSPixelLocalTrack > &tracksPixel)
TH1D * chisqn_lin_selected
Base class for CTPPS detector IDs.
std::string fileNamePrefix
file name prefix for result files
bool preciseXMLFormat
whether to use long format (many decimal digits) when saving XML files
EventNumber_t event() const
std::map< std::set< unsigned int >, RPSetPlots > rpSetPlots
chi^2 histograms per RP set
void printLineSeparator(const std::vector< std::map< unsigned int, AlignmentResult > > &)
void buildConstraints(std::vector< AlignmentConstraint > &)
builds a selected set of constraints
const DetGeometry & get(unsigned int id) const
retrieves sensor geometry
for(int i=first, nt=offsets[nh];i< nt;i+=gridDim.x *blockDim.x)
void saveDiagnostics() const
saves a ROOT file with diagnostic plots
TGraph * newGraph(const string &name, const string &title)
TGraph * fitBxVsByGraph_fitted
plots bx vs. by
ConstraintsType constraintsType
the chosen type of constraints
std::vector< unsigned int > rpIds
list of RPs for which the alignment parameters shall be optimized
unsigned int requireNumberOfUnits
select only tracks with activity in minimal number of units
void buildStandardConstraints(std::vector< AlignmentConstraint > &) const
builds the standard constraints
void factorRPFromSensorCorrections(const CTPPSRPAlignmentCorrectionsData &input, CTPPSRPAlignmentCorrectionsData &expanded, CTPPSRPAlignmentCorrectionsData &factored, const AlignmentGeometry &, bool equalWeights=false, unsigned int verbosity=0)
static std::string setToString(const std::set< unsigned int > &)
converts a set to string
double maxTrackAx
cuts on absolute values of the track angle
std::map< std::set< unsigned int >, unsigned long > selectedTracksPerRPSet
counter of selected tracks in a certain RP set
double bx
intercepts in mm
void printId(unsigned int id)
std::vector< Hit > HitCollection
void printN(const char *str, unsigned int N)
result pretty printing routines
static void buildGeometry(const std::vector< unsigned int > &rpDecIds, const std::vector< unsigned int > &excludedSensors, const CTPPSGeometry *, double z0, AlignmentGeometry &geometry)
builds the alignment geometry
double chiSqPerNdfCut
the value of chi^2/ndf cut threshold
Local (linear) track description (or a fit result). Uses global reference system.
Container for CTPPS RP alignment corrections. The corrections are stored on two levels - RP and senso...
TGraph * fitAxVsAyGraph_selected
T getParameter(std::string const &) const
std::vector< QuantityClass > quantityClasses
list of quantity classes to be optimized
AlignmentGeometry geometry
the geometry for this task
void printQuantitiesLine(const std::vector< std::map< unsigned int, AlignmentResult > > &)
bool isValidSensorId(unsigned int id) const
check whether the sensor Id is valid (present in the map)
void print() const
Prints the geometry.
bool saveXMLUncertainties
whether to save uncertainties in the result XML files
Power< A, B >::type pow(const A &a, const B &b)
Alignment correction for an element of the CT-PPS detector. Within the geometry description,...
TH1D * fitAxHist_selected
fit ax histograms for all/selected tracks
std::vector< AlignmentAlgorithm * > algorithms
the collection of the alignment algorithms
TH1D * fitAyHist_selected
fit ay histograms for all/selected tracks
virtual void begin(edm::ESHandle< CTPPSRPAlignmentCorrectionsData > hRealAlignment, edm::ESHandle< CTPPSGeometry > hRealGeometry, edm::ESHandle< CTPPSGeometry > hMisalignedGeometry)
virtual ~StraightTrackAlignment()
detector shifts in second readout direction
void setSensorCorrection(unsigned int id, const CTPPSRPAlignmentCorrectionData &ac)
sets the alignment correction for the given sensor
TH1D * fitByHist_selected
fit by histograms for all/selected tracks
SeedingHitSet::ConstRecHitPointer Hit
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
std::string expandedFileNamePrefix
file name prefix for cumulative expanded result files