calorimeter
Contents of reconstruction/calorimeter
-
class AnomalyGenerator : public Gaugi::AlgTool
- #include <AnomalyGenerator.h>
Tool to inject anomalies and defects into the calorimeter simulation.
This tool allows the user to simulate dead modules (by zeroing out the pulse) or noisy cells (by adding extra gaussian noise) based on configuration properties and event numbers. This is useful for testing the robustness of reconstruction algorithms against detector defects.
Properties:
DeadModules: List of booleans indicating if a block of events/cells is dead.
Cells: List of cell hashes to apply the anomaly to.
NoiseMean/Std: Parameters for the gaussian noise injection.
EventNumberRange: Range of event numbers where the anomaly is active.
This tool can simulate dead modules or add extra noise to specific cells to mimic realistic detector conditions.
Public Functions
-
AnomalyGenerator(std::string name)
Constructor
-
virtual ~AnomalyGenerator()
=====================================================================
-
virtual StatusCode initialize() override
=====================================================================
-
virtual StatusCode finalize() override
=====================================================================
-
virtual StatusCode execute(SG::EventContext &ctx, Gaugi::EDM *edm) const override
Apply anomalies to the current cell.
=====================================================================
Executes the anomaly injection.
Checks if the current event number and cell hash match any of the configured anomaly rules. If a match is found, it either zeros the pulse (dead module) or adds Gaussian noise.
- Parameters:
ctx – Event context.
edm – Pointer to the CaloDetDescriptor (cell).
ctx – EventContext.
edm – Pointer to the xAOD::CaloDetDescriptor (cell).
Private Functions
-
void AddGaussianNoise(std::vector<float> &pulse, float noiseMean, float noiseStddev) const
=====================================================================
Private Members
-
float m_noiseMean
-
float m_noiseStd
-
std::vector<bool> m_deadModules
-
std::vector<std::vector<int>> m_cells
-
std::vector<std::vector<int>> m_eventNumberRange
-
std::vector<float> m_noiseStdFactor
-
std::string m_inputEventKey
-
int m_outputLevel
Output level message
-
mutable TRandom3 m_rng
Random generator
-
class CaloCellMaker : public Gaugi::Algorithm
- #include <CaloCellMaker.h>
Constructs calorimeter cells from simulated hits.
Algorithm to build calorimeter cells from hits.
This algorithm is responsible for the digitization process. It reads energy deposits (hits), applies pulse shape simulation, noise, and other electronic effects to produce digitized cells. It works on a per-sampling basis (e.g., separate instance for ECAL Barrel Layer 1).
Properties:
InputHitsKey: StoreGate key for input energy deposits.
OutputCollectionKey: StoreGate key for output cells.
Eta/PhiBins: Binning definition for the readout segmentation.
ZMin/Max: Longitudinal extent of the calorimeter module.
Sampling: Calorimeter sampling ID (e.g., EMB1, Tile1).
BunchIdStart/End: Readout window in bunch crossings.
This algorithm is responsible for converting simulated hits (energy deposits) into calorimeter cells. It defines the readout geometry and applies the digitization steps (pulse shaping, noise, etc) by scheduling sub-tools.
Public Functions
-
CaloCellMaker(std::string name)
Contructor
-
~CaloCellMaker() = default
Destructor
-
virtual StatusCode initialize() override
=====================================================================
initialize the algorithm
-
virtual StatusCode bookHistograms(SG::EventContext &ctx) const override
=====================================================================
Book all histograms into the current storegate
-
virtual StatusCode execute(SG::EventContext &ctx, const G4Step *step) const override
=====================================================================
Execute in step action step from geant core
-
virtual StatusCode execute(SG::EventContext &ctx, int) const override
=====================================================================
Execute in ComponentAccumulator
-
virtual StatusCode pre_execute(SG::EventContext &ctx) const override
=====================================================================
execute before start the step action
-
virtual StatusCode post_execute(SG::EventContext &ctx) const override
=====================================================================
execute after the step action
Post-execution step: Digitization loop.
Iterates over all hits associated with this sampling.
Accumulates energy deposits per bunch crossing.
Runs the PulseGenerator to simulate the electronic signal shape.
Runs other tools (e.g., OptimalFilter, CrossTalk) to produce the final cell object.
-
virtual StatusCode fillHistograms(SG::EventContext &ctx) const override
=====================================================================
fill histogram in the end
-
virtual StatusCode finalize() override
=====================================================================
finalize the algorithm
-
void push_back(Gaugi::AlgTool *tool)
Add a tool to the execution list.
=====================================================================
- Parameters:
tool – Pointer to the AlgTool to be added.
-
void setPulseGenerator(Gaugi::AlgTool *tool)
Set the pulse generator tool.
- Parameters:
tool – Pointer to the PulseGenerator tool.
Private Functions
-
int find(const std::vector<float> &vec, float value) const
=====================================================================
-
unsigned long int hash(unsigned bin) const
=====================================================================
Private Members
-
std::string m_collectionKey
output collection key
-
std::string m_hitsKey
input hits key
-
std::string m_histPath
Base histogram path
-
int m_sampling
Sampling id for this reconstruction
-
int m_segment
Segment index for this sample calorimeter
-
int m_detector
Detector id for this sampling
-
int m_bcid_start
The start bunch crossing id for energy estimation
-
int m_bcid_end
The end bunch crossing id for energy estimation
-
int m_bc_nsamples
The number of samples per bunch crossing
-
float m_bc_duration
The time space (in ns) between two bunch crossings
-
std::vector<float> m_etaBins
-
std::vector<float> m_phiBins
-
float m_zMin
-
float m_zMax
-
float m_z
-
bool m_detailedHistograms
-
unsigned int m_nEtaBins
-
unsigned int m_nPhiBins
-
class CaloCellMerge : public Gaugi::Algorithm
- #include <CaloCellMerge.h>
Algorithm to merge multiple cell collections into a single container.
Algorithm to merge multiple cell collections.
This algorithm gathers partial cell collections (e.g. from different samplings like EMB1, EMB2, Tile1, etc.) produced by separate CaloCellMaker instances and merges them into a single global CaloCellContainer. It handles both reconstructed cells and truth cells.
Properties:
InputCollectionKeys: List of keys for the partial collections.
OutputCellsKey: Key for the final merged reconstructed cell container.
OutputTruthCellsKey: Key for the final merged truth cell container.
Collects cells from different calorimeter samplings (created by separate CaloCellMaker instances) and merges them into a single global CaloCellContainer.
Public Functions
-
CaloCellMerge(std::string name)
Contructor
-
~CaloCellMerge()
=====================================================================
Destructor
-
virtual StatusCode initialize() override
=====================================================================
initialize the algorithm
-
virtual StatusCode bookHistograms(SG::EventContext &ctx) const override
=====================================================================
Book all histograms into the current storegate
-
virtual StatusCode execute(SG::EventContext &ctx, const G4Step *step) const override
=====================================================================
Execute in step action step from geant core
-
virtual StatusCode execute(SG::EventContext &ctx, int) const override
=====================================================================
Execute in ComponentAccumulator
-
virtual StatusCode pre_execute(SG::EventContext &ctx) const override
=====================================================================
execute before start the step action
-
virtual StatusCode post_execute(SG::EventContext &ctx) const override
=====================================================================
execute after the step action
Core merge logic.
Iterates through all input keys, retrieves the corresponding
CaloDetDescriptorCollection, creates newxAOD::CaloCellobjects (copying info from descriptors), and pushes them into the output containers.
-
virtual StatusCode fillHistograms(SG::EventContext &ctx) const override
=====================================================================
fill hisogram in the end
-
virtual StatusCode finalize() override
=====================================================================
finalize the algorithm
-
class CaloClusterMaker : public Gaugi::Algorithm
- #include <CaloClusterMaker.h>
Reconstruction algorithm to perform Seeded-Clustering.
Algorithm to reconstruct calorimeter clusters.
This algorithm groups calorimeter cells into clusters initiated by seed particles (typically truth particles or reconstructed candidates). It scans for the “hottest” cell near a seed and aggregates energy within a fixed rectangular window in eta-phi.
Properties:
InputCellsKey: Collection of calorimeter cells.
InputSeedsKey: Collection of seeds (positions).
OutputClusterKey: Output collection of reconstructed clusters.
Eta/PhiWindow: Size of the clustering window.
MinCenterEnergy: Minimum energy required for the central cell to seed a cluster.
Uses a seed-based approach to group calorimeter cells into clusters. It identifies the highest energy cell (seed) and aggregates surrounding cells within a defined window.
Public Functions
-
CaloClusterMaker(std::string)
Constructor
-
virtual ~CaloClusterMaker()
=====================================================================
-
virtual StatusCode initialize() override
=====================================================================
-
virtual StatusCode bookHistograms(SG::EventContext &ctx) const override
=====================================================================
Book all histograms into the current storegate
-
virtual StatusCode execute(SG::EventContext &ctx, const G4Step *step) const override
=====================================================================
Execute in step action step from geant core
-
virtual StatusCode execute(SG::EventContext &ctx, int) const override
=====================================================================
Execute in ComponentAccumulator
-
virtual StatusCode pre_execute(SG::EventContext &ctx) const override
=====================================================================
execute before start the step action
-
virtual StatusCode post_execute(SG::EventContext &ctx) const override
=====================================================================
execute after the step action
Core clustering logic executed for each event.
Retrieves seeds and cell containers.
Iterates over seeds to find the corresponding “hot cell” (maximum energy) in the Second Layer (EMB2/EMEC2).
Verifies if the energy in a 0.1x0.1 core is above threshold (MinCenterEnergy).
If qualified, creates a CaloCluster, collects all cells within the Eta/Phi window, and calculates shower shapes.
-
virtual StatusCode fillHistograms(SG::EventContext &ctx) const override
=====================================================================
fill histogram in the end
-
virtual StatusCode finalize() override
=====================================================================
Private Functions
-
void fillCluster(SG::EventContext &ctx, xAOD::CaloCluster *clus, std::string key) const
=====================================================================
-
float dR(float eta1, float phi1, float eta2, float phi2) const
=====================================================================
-
class CaloHitMaker : public Gaugi::Algorithm
- #include <CaloHitMaker.h>
Collects energy deposits (Hits) during Geant4 simulation.
Algorithm to build calorimeter hits during simulation.
This algorithm runs at every Geant4 step. It checks if the energy deposit occurred within the defined volume (Sampling) and integrates the energy into a CaloHit object. It maps geometric position (x,y,z) to readout identifiers (eta, phi bins).
Properties:
OutputCollectionKey: StoreGate key for the collection of hits.
Eta/PhiBins: Readout segmentation.
RMin/Max, ZMin/Max: Spatial boundaries of the sensitive volume.
Sampling: Identifier for the calorimeter layer.
Collects energy deposits from Geant4 steps and integrates them into hits corresponding to readout cells.
Public Functions
-
CaloHitMaker(std::string name)
Contructor
-
~CaloHitMaker() = default
Destructor
-
virtual StatusCode initialize() override
=====================================================================
initialize the algorithm
-
virtual StatusCode bookHistograms(SG::EventContext &ctx) const override
=====================================================================
Book all histograms into the current storegate
-
virtual StatusCode execute(SG::EventContext &ctx, const G4Step *step) const override
=====================================================================
Execute in step action step from geant core
Execution per Geant4 Step.
Invoked by the framework for every particle step.
Checks if the step is within the geometric bounds of this calorimeter volume.
Calculates the corresponding eta/phi bin.
Retrieves the specific Hit object for that bin.
Adds the energy deposit to the Hit.
-
virtual StatusCode execute(SG::EventContext &ctx, int) const override
=====================================================================
Execute in ComponentAccumulator
-
virtual StatusCode pre_execute(SG::EventContext &ctx) const override
=====================================================================
execute before start the step action
Prepare the hit collection before processing steps.
Initializes the CaloHitCollection in the event store. Pre-populates it with empty hits for every defined readout cell (bin) to be ready for energy accumulation. This ensures that every cell has a corresponding object, even if empty.
-
virtual StatusCode post_execute(SG::EventContext &ctx) const override
=====================================================================
execute after the step action
-
virtual StatusCode fillHistograms(SG::EventContext &ctx) const override
=====================================================================
fill histogram in the end
-
virtual StatusCode finalize() override
=====================================================================
Private Functions
-
int find(const std::vector<float> &vec, float value) const
=====================================================================
-
unsigned long int hash(unsigned bin) const
=====================================================================
Private Members
-
std::string m_collectionKey
collection key
-
std::vector<float> m_etaBins
-
std::vector<float> m_phiBins
-
float m_rMin
-
float m_rMax
-
float m_zMin
-
float m_zMax
-
float m_noiseStd
-
int m_sampling
Sampling id for this reconstruction
-
int m_segment
Sampling id segment
-
int m_detector
Detector id for this sampling
-
int m_bcid_start
The start bunch crossing id for energy estimation
-
int m_bcid_end
The end bunch crossing id for energy estimation
-
float m_bc_duration
The time space (in ns) between two bunch crossings
-
std::string m_histPath
Base histogram path
-
bool m_detailedHistograms
detailed histogram flags
-
unsigned int m_nEtaBins
-
unsigned int m_nPhiBins
-
class CaloHitMerge : public Gaugi::Algorithm
- #include <CaloHitMerge.h>
Merges multiple CaloHit collections.
Algorithm to merge hit collections.
Collects
xAOD::CaloHitCollectionobjects (pointers to hits) from different sources and merges them into a singlexAOD::CaloHitContainer(owning the hits). This is typically the final step of hit generation before digitization.Properties:
InputCollectionKeys: Keys of collections to merge.
OutputHitsKey: Key of the output container.
Aggregates separate hit collections (usually from different detector regions) into a single container.
Public Functions
-
CaloHitMerge(std::string name)
Contructor
-
~CaloHitMerge()
=====================================================================
Destructor
-
virtual StatusCode initialize() override
=====================================================================
initialize the algorithm
-
virtual StatusCode bookHistograms(SG::EventContext &ctx) const override
=====================================================================
Book all histograms into the current storegate
-
virtual StatusCode execute(SG::EventContext &ctx, const G4Step *step) const override
=====================================================================
Execute in step action step from geant core
-
virtual StatusCode execute(SG::EventContext &ctx, int) const override
=====================================================================
Execute in ComponentAccumulator
-
virtual StatusCode pre_execute(SG::EventContext &ctx) const override
=====================================================================
execute before start the step action
-
virtual StatusCode post_execute(SG::EventContext &ctx) const override
=====================================================================
execute after the step action
Merges the collections.
Iterates over all input collections, creates new specific Hit objects (copying data), and adds them to the output container.
-
virtual StatusCode fillHistograms(SG::EventContext &ctx) const override
=====================================================================
fill hisogram in the end
-
virtual StatusCode finalize() override
=====================================================================
finalize the algorithm
-
class RingSet
- #include <CaloRingsMaker.h>
Helper class to manage a set of rings for a specific layer.
accumulating energy in concentric rings centered around a specific point.
Public Functions
-
RingSet(std::vector<CaloSampling> &samplings, unsigned nrings, float deta, float dphi)
Constructor
=====================================================================
=====================================================================
-
RingSet() = default
Destructor
-
void push_back(const xAOD::CaloCell *cell, float eta_center, float phi_center)
Add the cell energy to the correct ring position.
=====================================================================
Checks validity of the cell (sampling) and calculates its distance to the center to assign energy to the appropriate ring.
Adds a cell’s energy to the appropriate ring.
calculates the distance (deta, dphi) from the center. If the cell falls within the ring acceptance, adds its energy normalized by cosh(eta) (Transverse Energy) to the ring sum.
- Parameters:
cell – Pointer to the calorimeter cell.
eta_center – Eta of the ring center.
phi_center – Phi of the ring center.
cell – The calorimeter cell.
eta_center – Center eta of the ring system.
phi_center – Center phi of the ring system.
-
const std::vector<float> &rings() const
=====================================================================
Get the ringer shaper pattern for this RingSet
-
size_t size() const
=====================================================================
The number of rings in this RingSet
-
bool isValid(const xAOD::CaloCell*) const
Check if cell belongs to this RingSet’s samplings.
=====================================================================
- Returns:
true if valid, false otherwise.
-
void clear()
=====================================================================
Zeroize all energy values
-
RingSet(std::vector<CaloSampling> &samplings, unsigned nrings, float deta, float dphi)
-
class CaloRingsMaker : public Gaugi::Algorithm
- #include <CaloRingsMaker.h>
Algorithm to build RingSets (Ringer variable) from CaloClusters.
Algorithm to build RingSets from CaloClusters.
“Rings” are concentric energy sums around the cluster center, calculated per longitudinal layer. This compressed representation is often used for fast electron identification (neural networks).
Properties:
InputClusterKey: Input clusters to calculate rings for.
OutputRingerKey: Output ring container.
DeltaEta/PhiRings: Segmentation step size for rings.
NRings: Number of rings per layer.
LayerRings: Definitions of which sampling layers belong to which RingSet.
This algorithm creates the “Rings” feature set used for electron identification. It iterates over clusters and computes energy sums in concentric rings for each calorimeter layer.
Public Functions
-
CaloRingsMaker(std::string)
Constructor
-
virtual ~CaloRingsMaker()
=====================================================================
-
virtual StatusCode initialize() override
=====================================================================
-
virtual StatusCode bookHistograms(SG::EventContext &ctx) const override
=====================================================================
-
virtual StatusCode pre_execute(SG::EventContext &ctx) const override
=====================================================================
-
virtual StatusCode execute(SG::EventContext &ctx, const G4Step *step) const override
=====================================================================
-
virtual StatusCode execute(SG::EventContext &ctx, int) const override
=====================================================================
Execute in ComponentAccumulator
-
virtual StatusCode post_execute(SG::EventContext &ctx) const override
=====================================================================
Post-execution step to create rings
Core Rings calculation logic.
Reads the CaloClusters.
Configures the RingSet objects based on layer definitions.
For each cluster:
Identifies the max-energy cell to center the rings.
Iterates over all cells associated with the cluster.
Assigns cells to the appropriate RingSet and Ring index (radius).
Stores the resulting CaloRings object.
-
virtual StatusCode fillHistograms(SG::EventContext &ctx) const override
=====================================================================
-
virtual StatusCode finalize() override
=====================================================================
Private Functions
Private Members
-
int m_maxRingSets
-
int m_maxRingsAccumulated
-
std::string m_histPath
-
std::string m_clusterKey
-
std::string m_ringerKey
-
std::vector<float> m_detaRings
-
std::vector<float> m_dphiRings
-
std::vector<int> m_nRings
-
std::vector<float> m_etaRange
-
std::vector<std::vector<int>> m_layerRings
-
int m_outputLevel
-
bool m_doForward
-
bool m_DoSigmaCut
-
float m_SigmaCut
-
class CaloRingsMerge : public Gaugi::Algorithm
- #include <CaloRingsMerge.h>
Merges multiple CaloRings containers.
Algorithm to merge CaloRings containers.
Reads multiple
xAOD::CaloRingsContainerobjects from StoreGate and consolidates them into a single output container.Properties:
CollectionKeys: List of input keys.
OutputRingerKey: Output key.
Merges multiple CaloRings containers into a single one. This is useful when rings are built in parallel or in different steps.
Public Functions
-
CaloRingsMerge(std::string name)
Contructor
-
~CaloRingsMerge()
=====================================================================
Destructor
-
virtual StatusCode initialize() override
=====================================================================
initialize the algorithm
-
virtual StatusCode bookHistograms(SG::EventContext &ctx) const override
=====================================================================
Book all histograms into the current storegate
-
virtual StatusCode execute(SG::EventContext &ctx, const G4Step *step) const override
=====================================================================
Execute in step action step from geant core
-
virtual StatusCode execute(SG::EventContext &ctx, int) const override
=====================================================================
Execute in ComponentAccumulator
-
virtual StatusCode pre_execute(SG::EventContext &ctx) const override
=====================================================================
execute before start the step action
-
virtual StatusCode post_execute(SG::EventContext &ctx) const override
=====================================================================
execute after the step action
Performs the merge.
Iterates through input keys, retrieves containers, and deep-copies rings into the new master container.
-
virtual StatusCode fillHistograms(SG::EventContext &ctx) const override
=====================================================================
fill hisogram in the end
-
virtual StatusCode finalize() override
=====================================================================
finalize the algorithm
-
class ConstrainedOptimalFilter : public Gaugi::AlgTool
- #include <ConstrainedOptimalFilter.h>
Calculates energy and time using the Constrained Optimal Filtering method.
AlgTool for signal reconstruction using a Constrained Optimal Filter.
This tool reconstructs the cell energy from the digitized samples. Unlike the standard OF, it calculates the weights dynamically or uses specific constraints (e.g. baseline restoration) to minimize noise and pileup effects. It builds the autocorrelation matrix and solves the linear system to find the amplitude.
Properties:
PulsePath: Path to the reference pulse shape file.
Threshold: Minimum amplitude threshold for processing.
NSamples: Number of samples used in the filter.
Reconstructs the amplitude and time of the signal from the digitized samples using the Optimal Filtering technique with additional constraints (e.g., pedestal constraints).
Public Functions
-
ConstrainedOptimalFilter(std::string name)
Constructor
-
virtual ~ConstrainedOptimalFilter()
=====================================================================
-
virtual StatusCode initialize() override
=====================================================================
-
virtual StatusCode finalize() override
=====================================================================
-
void ReadShaper(std::string filepath)
-
void GeneratePulse(std::vector<float> &pulse) const
-
virtual StatusCode execute(SG::EventContext &ctx, Gaugi::EDM *edm) const override
Apply the filter to the cell.
=====================================================================
Executes the COF algorithm.
Generates the reference pulse.
Builds the covariance matrix (H).
Inverts the matrix to solve for the amplitude (a_hat).
Applies an iterative procedure to select valid samples (passing threshold).
Re-calculates the final energy using the selected samples.
- Parameters:
ctx – Event context.
edm – Pointer to the CaloDetDescriptor (cell).
-
class CrossTalkMaker : public Gaugi::Algorithm
- #include <CrossTalkMaker.h>
Simulates cell-to-cell cross-talk effects.
Algorithm to simulate cell-to-cell cross-talk.
This algorithm modifies the pulse shapes of calorimeter cells by mixing in contributions from their neighbors. It models both capacitive and inductive coupling. It creates a new collection of “cross-talked” cells.
Properties:
AmpCapacitive/Inductive: Coupling amplitudes.
MinEnergy: Process only cells above this energy threshold.
Simulates the leakage of signal from one cell to its neighbors due to capacitive or inductive coupling in the readout electronics.
Public Functions
-
CrossTalkMaker(std::string name)
Constructor
-
~CrossTalkMaker() = default
-
virtual StatusCode initialize() override
=====================================================================
-
virtual StatusCode bookHistograms(SG::EventContext &ctx) const override
=====================================================================
Book all histograms into the current storegate
-
virtual StatusCode pre_execute(SG::EventContext &ctx) const override
=====================================================================
execute before start the step action
-
virtual StatusCode execute(SG::EventContext &ctx, const G4Step *step) const override
=====================================================================
Execute in step action step from geant core
-
virtual StatusCode execute(SG::EventContext &ctx, int) const override
=====================================================================
Execute in ComponentAccumulator
Executes the cross-talk simulation.
Copies the original cells to a new container.
Iterates over valid central cells (above threshold).
Finds the 3x3 window of neighbor cells.
Calculates the distorted pulse for the central cell by summing contributions from neighbors (XTalkTF).
Updates the central cell’s pulse with the distorted version.
Runs downstream tools (e.g. OptimalFilter) on the modified cells.
-
virtual StatusCode post_execute(SG::EventContext &ctx) const override
=====================================================================
execute after the step action
-
virtual StatusCode finalize() override
=====================================================================
-
virtual StatusCode fillHistograms(SG::EventContext&) const override
=====================================================================
Fill all histograms into the current store gate
Private Functions
-
double XTalk(double x, bool type) const
-
double CellFunction(double x, bool type) const
-
float XTalkTF(float sample, int samp_index, bool diagonal, bool inductive) const
Private Members
-
std::vector<Gaugi::AlgTool*> m_toolHandles
The tool list that will be executed into the post execute step
-
float m_minEnergy
-
std::string m_collectionKey
-
std::string m_xtcollectionKey
-
std::string m_histPath
-
uint m_Nsamples = 5
-
uint m_tSamp = 25
-
double m_Cx = 47.000
-
double m_Rf = 0.078
-
double m_Rin = 1.200
-
double m_taud = 15.820
-
double m_taupa = 17.310
-
double m_td = 420.000
-
double m_tmax2 = 600.000
-
double m_C1 = 50.000
-
double m_ToNormXtC = 0.022206
-
double m_AmpXt_C
-
double m_AmpXt_L
-
double m_AmpXt_R
-
double m_AmpNoise = 50
-
double tau_0_mean = 0
-
double tau_std = 0.5
Warning
doxygenfile: Found multiple matches for file “LinkDef.h
Warning
doxygenfile: Found multiple matches for file “LinkDef.h
Warning
doxygenfile: Found multiple matches for file “LinkDef.h
Warning
doxygenfile: Found multiple matches for file “LinkDef.h
-
class OptimalFilter : public Gaugi::AlgTool
- #include <OptimalFilter.h>
Reconstructs cell energy and time using fixed Optimal Filter weights.
AlgTool for signal reconstruction using standard Optimal Filtering (OF).
Applies the standard Optimal Filtering technique (OF2 usually) where the energy and time are linear combinations of the digitized samples. The weights are pre-calculated to minimize noise and pileup.
Properties:
WeightsEnergy: Vector of weights for energy estimation.
WeightsTime: Vector of weights for time estimation.
Uses pre-calculated weights to estimate the energy and time from the digitized samples.
Public Functions
-
OptimalFilter(std::string name)
Constructor
-
virtual ~OptimalFilter()
=====================================================================
-
virtual StatusCode initialize() override
=====================================================================
-
virtual StatusCode finalize() override
=====================================================================
-
virtual StatusCode execute(SG::EventContext &ctx, Gaugi::EDM *edm) const override
Apply the optimal filter to the cell.
=====================================================================
Calculates Energy and Time.
Computes the dot product of the pulse samples with the energy/time weights. Sets the reconstructed Energy and Tau in the cell object.
- Parameters:
ctx – Event context.
edm – Pointer to the CaloDetDescriptor (cell).
-
class PileupMerge : public Gaugi::Algorithm
- #include <PileupMerge.h>
Merges pileup events (minbias) into the signal event.
Algorithm to overlay pileup events onto a hard-scatter event.
This algorithm simulates the effect of high-luminosity pileup by overlaying hits from pre-generated minimum bias events onto the current signal event. It handles the time window (bunch crossings) and uses Poisson statistics to determine the number of pileup interactions per bunch crossing.
Properties:
PileupAvg: Average number of pileup interactions (mu).
PileupSigma: Fluctuations in mu.
BunchIdStart/End: Time window in bunch crossings (-20 to 20 usually).
reads simulated hits from “minbias” events (pileup) and adds their energy contributions to the hits of the main signal event, simulating high-luminosity conditions.
Public Functions
-
PileupMerge(std::string)
Constructor
-
virtual ~PileupMerge()
=====================================================================
-
virtual StatusCode initialize() override
=====================================================================
-
virtual StatusCode bookHistograms(SG::EventContext &ctx) const override
=====================================================================
-
virtual StatusCode pre_execute(SG::EventContext &ctx) const override
=====================================================================
-
virtual StatusCode execute(SG::EventContext &ctx, const G4Step *step) const override
=====================================================================
-
virtual StatusCode execute(SG::EventContext&, int) const override
=====================================================================
-
virtual StatusCode post_execute(SG::EventContext &ctx) const override
=====================================================================
Executes the pileup merging process.
Reads the signal event hits.
Attempts to merge the pileup hits (with retries in case of I/O errors).
Saves the merged hit container and updated EventInfo (with new avgMu).
-
virtual StatusCode fillHistograms(SG::EventContext &ctx) const override
=====================================================================
-
virtual StatusCode finalize() override
=====================================================================
Private Functions
-
template<class T>
TBranch *InitBranch(TTree *fChain, std::string branch_name, T *param) const =====================================================================
-
int poisson(double nAvg) const
=====================================================================
-
void Read(SG::EventContext &ctx, const std::vector<std::string> &paths, std::string name) const
=====================================================================
-
void allocate(SG::ReadHandle<xAOD::CaloHitContainer> &container, std::vector<xAOD::CaloHit*> &vec_hits) const
=====================================================================
-
float merge(SG::EventContext &ctx, std::vector<xAOD::CaloHit*> &vec_hits) const
=====================================================================
Helper function to perform the actual merging.
Iterates through bunch crossings (BCID) from start to end. For each BCID, it determines the number of pileup interactions (Poisson) and randomly selects events from the low/high pileup input files/trees. It then reads the hits from those events and adds their energy to the corresponding signal hits.
- Returns:
The average number of pileup interactions added per bunch crossing.
Private Members
-
std::string m_inputHitsKey
-
std::string m_outputHitsKey
-
std::string m_inputEventKey
-
std::string m_outputEventKey
-
int m_outputLevel
-
mutable TRandom3 m_rng
-
float m_pileupAvg
-
float m_pileupSigma
-
float m_seed
-
int m_maxRetry
-
float m_trunc_mu
-
int m_bcid_start
The start bunch crossing id for energy estimation
-
int m_bcid_end
The end bunch crossing id for energy estimation
-
std::vector<std::string> m_lowPileupInputFiles
-
std::vector<std::string> m_highPileupInputFiles
-
std::string m_ntupleName
-
class PulseGenerator : public Gaugi::AlgTool
- #include <PulseGenerator.h>
Generates electronic pulse shapes for calorimeter cells.
Tool to simulate the electronic pulse shape.
Simulates the response of the readout electronics to an energy deposit. It uses a reference shaper function (read from a file) to convolve the energy deposit in time. It also adds electronic noise and random deformations.
Properties:
ShaperFile: File containing the reference pulse shape points.
Pedestal: Baseline voltage.
NoiseMean/Std: Electronic noise parameters.
SamplingRate: Readout sampling rate (usually 25ns).
This tool takes the energy deposit in a cell and generates a time-sampled electronic pulse (using a shaper function). It also adds electronic noise and can simulate defects.
Public Functions
-
PulseGenerator(std::string name)
Constructor
-
virtual ~PulseGenerator()
=====================================================================
-
virtual StatusCode initialize() override
=====================================================================
-
virtual StatusCode finalize() override
=====================================================================
-
virtual StatusCode execute(SG::EventContext &ctx, Gaugi::EDM *edm) const override
Execute the pulse generation for a specific cell.
=====================================================================
Generates the pulse for a cell.
Loops over bunch crossings and sums the contribution of energy deposits from each BCID, weighted by the shaper function at the appropriate time lag. Adds noise and sets the final pulse in the cell object.
- Parameters:
ctx – Event context.
edm – Pointer to the CaloDetDescriptor (the cell).
- Returns:
Status code indicating success or failure.
Private Functions
-
void ReadShaper(std::string)
=====================================================================
-
void GenerateDeterministicPulse(std::vector<float> &pulse, float amplitude, float phase, float lag) const
=====================================================================
-
void AddGaussianNoise(std::vector<float> &pulse, float noiseMean, float noiseStddev) const
=====================================================================
Private Members
-
int m_nsamples
Number of samples to be generated
-
int m_startSamplingBC
-
int m_shaperZeroIndex
-
float m_pedestal
-
float m_deformationMean
-
float m_deformationStd
-
float m_samplingRate
-
float m_shaperResolution
-
float m_noiseMean
-
float m_noiseStd
-
bool m_doDefects
-
bool m_deadModules
-
std::vector<std::vector<int>> m_cellHash
-
std::vector<float> m_noiseFactor
-
std::vector<std::vector<int>> m_noisyEvents
-
std::vector<float> m_shaper
-
std::vector<float> m_timeSeries
-
std::string m_shaperFile
The shaper configuration path
-
int m_outputLevel
Output level message
-
mutable TRandom3 m_rng
Random generator
-
class ShowerShapes : public Gaugi::AlgTool
- #include <ShowerShapes.h>
Utility tool for calculating calorimetric shower shape variables.
Tool to calculate shower shape variables.
Computes various quantities describing the longitudinal and lateral development of the shower, such as Reta, Rphi, weta2, f1, f3, etc. Use for electron identification/discrimination.
Calculates discriminating variables used for particle identification (e.g., e/gamma separation) based on the energy distribution within a cluster.
Public Functions
-
ShowerShapes(std::string)
Constructor
-
virtual ~ShowerShapes() = default
-
virtual StatusCode initialize() override
=====================================================================
-
virtual StatusCode finalize() override
=====================================================================
-
virtual StatusCode execute(SG::EventContext &ctx, Gaugi::EDM*) const override
=====================================================================
Execution entry point.
Calculates all defined shower shapes for a given CaloCluster. Populates the cluster object with the computed variables (setters).
The calculated variables include:
Pre-Sampler: e0, f0
EM1 (Strip): e1, f1, emaxs1, e2tsts1, eratio, weta1 (not implemented but reserved)
EM2 (Middle): e2, f2, reta, rphi, weta2, e233, e237, e277
EM3 (Back): e3, f3
Hadronic: ehad1, ehad2, ehad3, rhad, rhad1
Forward Moments: secondR, lambdaCenter, secondLambda, fracMax, lateralMom, longitudinalMom (for forward electrons)
- Parameters:
ctx – EventContext (unused here but required by interface).
edm – Pointer to the CaloCluster object.
-
inline void setForwardMoments(bool doForward)
Enable or disable calculation of forward moments.
- Parameters:
doForward – Boolean flag.
Private Functions
-
float calculateWeta2(xAOD::CaloCluster*, unsigned eta_ncell = 3, unsigned phi_ncell = 5) const
=====================================================================
Calculates Weta2 (Lateral Shower Width in EM2).
Formula: sqrt( (Sum(E_i * eta_i^2) / TotalE) - (Sum(E_i * eta_i) / TotalE)^2 ) Basically the standard deviation of the energy distribution in eta.
- Parameters:
clus – Cluster object.
eta_ncell – Eta window size.
phi_ncell – Phi window size.
- Returns:
Weta2 value.
-
float sumEnergyEM(xAOD::CaloCluster*, int sampling, unsigned eta_ncell = 1000, unsigned phi_ncell = 1000) const
=====================================================================
Sums energy in specified EM layer within a window.
- Parameters:
clus – The cluster.
sampling – 0=PS, 1=EM1, 2=EM2, 3=EM3.
eta_ncell – Window size in eta (units of cell size).
phi_ncell – Window size in phi (units of cell size).
- Returns:
Sum of energy of cells matching validity criteria.
-
float sumEnergyHAD(xAOD::CaloCluster*, int sampling, unsigned eta_ncell = 1000, unsigned phi_ncell = 1000) const
=====================================================================
Sums energy in specified Hadronic layer within a window.
- Parameters:
clus – The cluster.
sampling – 0=HAD1, 1=HAD2, 2=HAD3.
eta_ncell – Window size in eta.
phi_ncell – Window size in phi.
- Returns:
Sum of energy of cells matching validity criteria.
-
float calculateSecondR(xAOD::CaloCluster*, std::vector<TVector3>) const
Calculates SecondR (Lateral spread perpendicular to the shower axis).
- Parameters:
axis – Result from calculateShowerAxis.
- Returns:
Second moment lateral.
-
float calculateSecondLambda(xAOD::CaloCluster*, std::vector<TVector3>) const
Calculates SecondLambda (Longitudinal spread along the shower axis).
- Parameters:
axis – Result from calculateShowerAxis.
- Returns:
Second moment longitudinal.
-
float calculateLambdaCenter(xAOD::CaloCluster*, std::vector<TVector3>) const
Calculates LambdaCenter (Distance of shower center from calorimeter front face along axis).
- Parameters:
axis – Result from calculateShowerAxis.
- Returns:
Distance in mm (approximation).
-
float calculateFracMax(xAOD::CaloCluster*, std::vector<TVector3>) const
Calculates FracMax (Fraction of energy in the max-E cell).
- Parameters:
axis – Result from calculateShowerAxis.
- Returns:
MaxCelEnergy / TotalEnergy
-
float calculateLateralMom(xAOD::CaloCluster*, std::vector<TVector3>) const
Calculates Lateral Moment (describes lateral shape).
Code derived from standard Forward Electron ID implementation.
- Parameters:
axis – Result from calculateShowerAxis.
- Returns:
Lateral moment.
-
float calculateLongitudinalMom(xAOD::CaloCluster*, std::vector<TVector3>) const
Calculates Longitudinal Moment (describes longitudinal shape).
Code derived from standard Forward Electron ID implementation.
- Parameters:
axis – Result from calculateShowerAxis.
- Returns:
Longitudinal moment.
-
std::vector<TVector3> calculateShowerAxis(xAOD::CaloCluster*) const
Calculates the Shower Axis using Principal Component Analysis (Eigen analysis).
Constructs the covariance matrix of energy deposits in space (3D). The principal eigenvector (highest eigenvalue) corresponds to the main shower axis. Used for Forward forward electron identification where Calo pointing is crucial.
- Returns:
Vector containing: [ShowerAxis, ShowerCenter, Energies(Total, Max1, Max2), EMEC1_POS]
Private Members
-
bool m_doForwardMoments
-
ShowerShapes(std::string)