Commit f1338fc0 authored by Jens Kroeger's avatar Jens Kroeger
Browse files

AnalysisDUT: Add new TProfiles showing in-pixel distribution of cluster charge...

AnalysisDUT: Add new TProfiles showing in-pixel distribution of cluster charge for different cluster sizes
parent 373b846a
Pipeline #2731830 passed with stages
in 15 minutes and 21 seconds
...@@ -419,6 +419,54 @@ void AnalysisDUT::initialize() { ...@@ -419,6 +419,54 @@ void AnalysisDUT::initialize() {
-pitch_y / 2., -pitch_y / 2.,
pitch_y / 2.); pitch_y / 2.);
title = "Mean cluster charge map (1-pixel);" + mod_axes + "<cluster charge> [ke]";
qvsxmym_1px = new TProfile2D("qvsxmym_1px",
title.c_str(),
static_cast<int>(pitch_x),
-pitch_x / 2.,
pitch_x / 2.,
static_cast<int>(pitch_y),
-pitch_y / 2.,
pitch_y / 2.,
0,
250);
title = "Mean cluster charge map (2-pixel);" + mod_axes + "<cluster charge> [ke]";
qvsxmym_2px = new TProfile2D("qvsxmym_2px",
title.c_str(),
static_cast<int>(pitch_x),
-pitch_x / 2.,
pitch_x / 2.,
static_cast<int>(pitch_y),
-pitch_y / 2.,
pitch_y / 2.,
0,
250);
title = "Mean cluster charge map (3-pixel);" + mod_axes + "<cluster charge> [ke]";
qvsxmym_3px = new TProfile2D("qvsxmym_3px",
title.c_str(),
static_cast<int>(pitch_x),
-pitch_x / 2.,
pitch_x / 2.,
static_cast<int>(pitch_y),
-pitch_y / 2.,
pitch_y / 2.,
0,
250);
title = "Mean cluster charge map (4-pixel);" + mod_axes + "<cluster charge> [ke]";
qvsxmym_4px = new TProfile2D("qvsxmym_4px",
title.c_str(),
static_cast<int>(pitch_x),
-pitch_x / 2.,
pitch_x / 2.,
static_cast<int>(pitch_y),
-pitch_y / 2.,
pitch_y / 2.,
0,
250);
residualsTime = new TH1F("residualsTime", residualsTime = new TH1F("residualsTime",
"Time residual;time_{track}-time_{hit} [ns];#entries", "Time residual;time_{track}-time_{hit} [ns];#entries",
n_timebins_, n_timebins_,
...@@ -634,8 +682,8 @@ StatusCode AnalysisDUT::run(const std::shared_ptr<Clipboard>& clipboard) { ...@@ -634,8 +682,8 @@ StatusCode AnalysisDUT::run(const std::shared_ptr<Clipboard>& clipboard) {
// Calculate in-pixel position of track in microns // Calculate in-pixel position of track in microns
auto inpixel = m_detector->inPixel(localIntercept); auto inpixel = m_detector->inPixel(localIntercept);
auto xmod = inpixel.X() * 1000.; // convert mm -> um auto xmod_um = inpixel.X() * 1000.; // convert mm -> um
auto ymod = inpixel.Y() * 1000.; // convert mm -> um auto ymod_um = inpixel.Y() * 1000.; // convert mm -> um
// Loop over all associated DUT clusters: // Loop over all associated DUT clusters:
for(auto assoc_cluster : track->getAssociatedClusters(m_detector->getName())) { for(auto assoc_cluster : track->getAssociatedClusters(m_detector->getName())) {
...@@ -752,11 +800,11 @@ StatusCode AnalysisDUT::run(const std::shared_ptr<Clipboard>& clipboard) { ...@@ -752,11 +800,11 @@ StatusCode AnalysisDUT::run(const std::shared_ptr<Clipboard>& clipboard) {
clusterWidthColAssoc->Fill(static_cast<double>(assoc_cluster->columnWidth())); clusterWidthColAssoc->Fill(static_cast<double>(assoc_cluster->columnWidth()));
// Fill in-pixel plots: (all as function of track position within pixel cell) // Fill in-pixel plots: (all as function of track position within pixel cell)
qvsxmym->Fill(xmod, ymod, cluster_charge); // cluster charge profile qvsxmym->Fill(xmod_um, ymod_um, cluster_charge); // cluster charge profile
qMoyalvsxmym->Fill(xmod, ymod, exp(-normalized_charge / 3.5)); // norm. cluster charge profile qMoyalvsxmym->Fill(xmod_um, ymod_um, exp(-normalized_charge / 3.5)); // norm. cluster charge profile
// mean charge of cluster seed // mean charge of cluster seed
pxqvsxmym->Fill(xmod, ymod, assoc_cluster->getSeedPixel()->charge()); pxqvsxmym->Fill(xmod_um, ymod_um, assoc_cluster->getSeedPixel()->charge());
if(assoc_cluster->size() > 1) { if(assoc_cluster->size() > 1) {
for(auto& px : assoc_cluster->pixels()) { for(auto& px : assoc_cluster->pixels()) {
...@@ -790,20 +838,28 @@ StatusCode AnalysisDUT::run(const std::shared_ptr<Clipboard>& clipboard) { ...@@ -790,20 +838,28 @@ StatusCode AnalysisDUT::run(const std::shared_ptr<Clipboard>& clipboard) {
} }
// mean cluster size // mean cluster size
npxvsxmym->Fill(xmod, ymod, static_cast<double>(assoc_cluster->size())); npxvsxmym->Fill(xmod_um, ymod_um, static_cast<double>(assoc_cluster->size()));
if(assoc_cluster->size() == 1) if(assoc_cluster->size() == 1) {
npx1vsxmym->Fill(xmod, ymod); npx1vsxmym->Fill(xmod_um, ymod_um);
if(assoc_cluster->size() == 2) qvsxmym_1px->Fill(xmod_um, ymod_um, cluster_charge);
npx2vsxmym->Fill(xmod, ymod); }
if(assoc_cluster->size() == 3) if(assoc_cluster->size() == 2) {
npx3vsxmym->Fill(xmod, ymod); npx2vsxmym->Fill(xmod_um, ymod_um);
if(assoc_cluster->size() == 4) qvsxmym_2px->Fill(xmod_um, ymod_um, cluster_charge);
npx4vsxmym->Fill(xmod, ymod); }
if(assoc_cluster->size() == 3) {
npx3vsxmym->Fill(xmod_um, ymod_um);
qvsxmym_3px->Fill(xmod_um, ymod_um, cluster_charge);
}
if(assoc_cluster->size() == 4) {
npx4vsxmym->Fill(xmod_um, ymod_um);
qvsxmym_4px->Fill(xmod_um, ymod_um, cluster_charge);
}
// residual MAD x, y, combined (sqrt(x*x + y*y)) // residual MAD x, y, combined (sqrt(x*x + y*y))
rmsxvsxmym->Fill(xmod, ymod, xabsdistance); rmsxvsxmym->Fill(xmod_um, ymod_um, xabsdistance);
rmsyvsxmym->Fill(xmod, ymod, yabsdistance); rmsyvsxmym->Fill(xmod_um, ymod_um, yabsdistance);
rmsxyvsxmym->Fill(xmod, ymod, fabs(sqrt(xdistance * xdistance + ydistance * ydistance))); rmsxyvsxmym->Fill(xmod_um, ymod_um, fabs(sqrt(xdistance * xdistance + ydistance * ydistance)));
hAssociatedTracksGlobalPosition->Fill(globalIntercept.X(), globalIntercept.Y()); hAssociatedTracksGlobalPosition->Fill(globalIntercept.X(), globalIntercept.Y());
hAssociatedTracksLocalPosition->Fill(m_detector->getColumn(localIntercept), m_detector->getRow(localIntercept)); hAssociatedTracksLocalPosition->Fill(m_detector->getColumn(localIntercept), m_detector->getRow(localIntercept));
......
...@@ -74,6 +74,7 @@ namespace corryvreckan { ...@@ -74,6 +74,7 @@ namespace corryvreckan {
TProfile2D *rmsxvsxmym, *rmsyvsxmym, *rmsxyvsxmym; TProfile2D *rmsxvsxmym, *rmsyvsxmym, *rmsxyvsxmym;
TProfile2D *qvsxmym, *qMoyalvsxmym, *pxqvsxmym; TProfile2D *qvsxmym, *qMoyalvsxmym, *pxqvsxmym;
TProfile2D *qvsxmym_1px, *qvsxmym_2px, *qvsxmym_3px, *qvsxmym_4px;
TProfile2D* npxvsxmym; TProfile2D* npxvsxmym;
TH2F *npx1vsxmym, *npx2vsxmym, *npx3vsxmym, *npx4vsxmym; TH2F *npx1vsxmym, *npx2vsxmym, *npx3vsxmym, *npx4vsxmym;
......
Markdown is supported
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment