Skip to content
Snippets Groups Projects
Commit e0efdf97 authored by Simon Spannagel's avatar Simon Spannagel
Browse files

Merge branch 'cluster_error_plots' into 'master'

Add cluster error plots

See merge request corryvreckan/corryvreckan!617
parents 6bcb08d7 e9c217c3
No related branches found
No related tags found
No related merge requests found
Pipeline #5342337 failed
......@@ -85,6 +85,31 @@ void Clustering4D::initialize() {
clusterTimes = new TH1F("clusterTimes", title.c_str(), 3e6, 0, 3e9);
title = m_detector->getName() + " Cluster multiplicity;clusters;events";
clusterMultiplicity = new TH1F("clusterMultiplicity", title.c_str(), 50, -0.5, 49.5);
title = m_detector->getName() + " Cluster Uncertainty x;cluster uncertainty x [um];events";
clusterUncertaintyX = new TH1F(
"clusterUncertaintyX", title.c_str(), 100, 0, static_cast<double>(Units::convert(m_detector->getPitch().X(), "um")));
title = m_detector->getName() + " Cluster Uncertainty y;cluster uncertainty y [um];events";
clusterUncertaintyY = new TH1F(
"clusterUncertaintyY", title.c_str(), 100, 0, static_cast<double>(Units::convert(m_detector->getPitch().Y(), "um")));
title = m_detector->getName() + " Cluster Uncertainty x vs position;column [px];row [px];<cluster uncertainty x> [um]";
clusterUncertaintyXvsXY = new TProfile2D("clusterUncertaintyXvsXY",
title.c_str(),
m_detector->nPixels().X(),
-0.5,
m_detector->nPixels().X() - 0.5,
m_detector->nPixels().Y(),
-0.5,
m_detector->nPixels().Y() - 0.5);
title = m_detector->getName() + " Cluster Uncertainty Y vs position;column [px];row [px];<cluster uncertainty y> [um]";
clusterUncertaintyYvsXY = new TProfile2D("clusterUncertaintyYvsXY",
title.c_str(),
m_detector->nPixels().X(),
-0.5,
m_detector->nPixels().X() - 0.5,
m_detector->nPixels().Y(),
-0.5,
m_detector->nPixels().Y() - 0.5);
title =
m_detector->getName() + " pixel - seed pixel timestamp (all pixels w/o seed);ts_{pixel} - ts_ {seed} [ns];events";
pxTimeMinusSeedTime = new TH1F("pxTimeMinusSeedTime", title.c_str(), 1000, -99.5 * 1.5625, 900.5 * 1.5625);
......@@ -207,6 +232,12 @@ StatusCode Clustering4D::run(const std::shared_ptr<Clipboard>& clipboard) {
clusterPositionGlobal->Fill(cluster->global().x(), cluster->global().y());
clusterPositionLocal->Fill(cluster->column(), cluster->row());
clusterTimes->Fill(static_cast<double>(Units::convert(cluster->timestamp(), "ns")));
clusterUncertaintyX->Fill(static_cast<double>(Units::convert(cluster->errorX(), "um")));
clusterUncertaintyY->Fill(static_cast<double>(Units::convert(cluster->errorY(), "um")));
clusterUncertaintyXvsXY->Fill(
cluster->column(), cluster->row(), static_cast<double>(Units::convert(cluster->errorX(), "um")));
clusterUncertaintyYvsXY->Fill(
cluster->column(), cluster->row(), static_cast<double>(Units::convert(cluster->errorY(), "um")));
// to check that cluster timestamp = earliest pixel timestamp
if(cluster->size() > 1) {
......
......@@ -15,6 +15,7 @@
#include <TCanvas.h>
#include <TH1F.h>
#include <TH2F.h>
#include <TProfile2D.h>
#include <iostream>
#include "core/module/Module.hpp"
#include "objects/Cluster.hpp"
......@@ -53,6 +54,11 @@ namespace corryvreckan {
TH2F* clusterPositionLocal;
TH1F* clusterTimes;
TH1F* clusterMultiplicity;
TH1F* clusterUncertaintyX;
TH1F* clusterUncertaintyY;
TProfile2D* clusterUncertaintyXvsXY;
TProfile2D* clusterUncertaintyYvsXY;
TH1F* pxTimeMinusSeedTime;
TH2F* pxTimeMinusSeedTime_vs_pxCharge;
TH2F* pxTimeMinusSeedTime_vs_pxCharge_2px;
......
......@@ -35,6 +35,7 @@ For each detector the following plots are produced:
* 2D cluster positions in global coordinates
* Cluster times
* Cluster multiplicity
* Cluster uncertainty
* Time difference between seed pixel and other pixels in a cluster for different cluster sizes and vs. the pixel charge
### Usage
......
......@@ -63,6 +63,12 @@ void ClusteringSpatial::initialize() {
clusterTimes = new TH1F("clusterTimes", title.c_str(), 3e6, 0, 3e9);
title = m_detector->getName() + " Cluster multiplicity;clusters;events";
clusterMultiplicity = new TH1F("clusterMultiplicity", title.c_str(), 50, -0.5, 49.5);
title = m_detector->getName() + " Cluster Uncertainty x;cluster uncertainty x [um];events";
clusterUncertaintyX = new TH1F(
"clusterUncertaintyX", title.c_str(), 100, 0, static_cast<double>(Units::convert(m_detector->getPitch().X(), "um")));
title = m_detector->getName() + " Cluster Uncertainty y;cluster uncertainty y [um];events";
clusterUncertaintyY = new TH1F(
"clusterUncertaintyY", title.c_str(), 100, 0, static_cast<double>(Units::convert(m_detector->getPitch().Y(), "um")));
}
StatusCode ClusteringSpatial::run(const std::shared_ptr<Clipboard>& clipboard) {
......@@ -180,6 +186,8 @@ StatusCode ClusteringSpatial::run(const std::shared_ptr<Clipboard>& clipboard) {
clusterPositionGlobal->Fill(cluster->global().x(), cluster->global().y());
clusterPositionLocal->Fill(cluster->column(), cluster->row());
clusterTimes->Fill(static_cast<double>(Units::convert(cluster->timestamp(), "ns")));
clusterUncertaintyX->Fill(static_cast<double>(Units::convert(cluster->errorX(), "um")));
clusterUncertaintyY->Fill(static_cast<double>(Units::convert(cluster->errorY(), "um")));
LOG(DEBUG) << "cluster local: " << cluster->local();
deviceClusters.push_back(cluster);
......
......@@ -48,6 +48,8 @@ namespace corryvreckan {
TH2F* clusterPositionGlobal;
TH2F* clusterPositionLocal;
TH1F* clusterTimes;
TH1F* clusterUncertaintyX;
TH1F* clusterUncertaintyY;
bool useTriggerTimestamp;
bool chargeWeighting;
......
......@@ -28,6 +28,7 @@ For each detector the following plots are produced:
* 2D cluster positions in global and local coordinates
* Cluster times
* Cluster multiplicity
* Cluster uncertainty
### Usage
```toml
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment