REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
Data Structures | Public Types | Public Member Functions | Static Public Member Functions | Protected Attributes
TRestHits Class Reference

Detailed Description

It saves a 3-coordinate position and an energy for each punctual deposition.

TRestHits is a generic data holder that defines an arbitrary physical quantity (usually the energy) in a 3-dimensional space (x,y,z). However, it also defines an optional time member that might be used to add additional information to the spatial event, such as it could be the drift time in a TPC. However, the final meaning of that member must be interpreted in the context of the event data processing algorithms, and/or the data type using the hits, such as: TRestGeant4Hits, TRestDetectorHits, ...

On top of that, we may also define the hit type in particular scenarios where one of the spatial coordinates remains unknown, and we may have a REST_HitType that defines a XY, XZ, XYZ, etc, hit type.

This class defines typical transformations required by spatially defined physical quantities such as rotation or translation, basic hit distance calculation, and hit operations such as merging, adding/removing or swaping hits.

It contains also more sophisticated methods to perform physical calculations and parameterize the properties of a group of hits or cluster such as performing a gaussian fit to the hit distribution, such as TRestHits::GetGaussSigmaX, where two hits are added, one to each side of the event, and a Gaussian is fitted. The hits are added so that the fit works even for small events as shown in the figure below. The parameter sigma is extracted from the fit and its absolute value is returned.

Other methods determine the number of hits or the total energy contained in a particular geometrical shape, see for example TRestHits::GetEnergyInCylinder, and different physical quantities on such fiducialization, i.e. TRestHits::GetMeanPositionInPrism.


RESTsoft - Software for Rare Event Searches with TPCs

History of developments:

2015-Sep: First concept. Created as part of the conceptualization of existing REST software.

Author
Javier Galan (javie.nosp@m.r.ga.nosp@m.lan@u.nosp@m.niza.nosp@m.r.es)

2015-Nov: Changed vectors fX fY fZ and fEnergy from <Int_t> to < Float_t>.

Author
JuanAn Garcia (juana.nosp@m.n318.nosp@m.@gmai.nosp@m.l.co.nosp@m.m)

2022-July: Introducing gausian hits fitting

Author
Cristina Margalejo (cmarg.nosp@m.alej.nosp@m.o@uni.nosp@m.zar..nosp@m.es)

Definition at line 39 of file TRestHits.h.

#include <TRestHits.h>

Inheritance diagram for TRestHits:
TRestGeant4Hits TRestVolumeHits

Data Structures

class  TRestHits_Iterator
 

Public Types

typedef TRestHits_Iterator iterator
 

Public Member Functions

void AddHit (const TVector3 &position, Double_t energy, Double_t time=0, REST_HitType type=XYZ)
 Adds a new hit to the list of hits using a TVector3. More...
 
void AddHit (TRestHits &hits, Int_t n)
 Adds a new hit to the list of hits using the hit n inside another TRestHits object. More...
 
virtual Bool_t areXY () const
 It will return true only if all the hits inside are of type XY. More...
 
virtual Bool_t areXYZ () const
 It will return true only if all the hits inside are of type XYZ. More...
 
virtual Bool_t areXZ () const
 It will return true only if all the hits inside are of type XZ. More...
 
virtual Bool_t areYZ () const
 It will return true only if all the hits inside are of type YZ. More...
 
TRestHits_Iterator begin ()
 
TRestHits_Iterator end ()
 
Int_t GetClosestHit (const TVector3 &position) const
 It returns the closest hit to a given position. More...
 
Double_t GetDistance (int N, int M) const
 
Double_t GetDistance2 (int n, int m) const
 It returns the euclidian distance between hits n and m. More...
 
Double_t GetDistanceToNode (Int_t n) const
 It determines the distance required to travel from the first hit to the hit n adding all the distances of the hits that are found till the hit n. More...
 
Double_t GetEnergy (int n) const
 
Double_t GetEnergyInCylinder (const TVector3 &x0, const TVector3 &x1, Double_t radius) const
 It determines the total energy contained inside a cylinder with a given radius and delimited between x0 and y0 vertex. More...
 
Double_t GetEnergyInCylinder (Int_t i, Int_t j, Double_t radius) const
 It determines the total energy contained inside a cylinder with a given radius and delimited between the hit number i and the hit number j. More...
 
Double_t GetEnergyInPrism (const TVector3 &x0, const TVector3 &x1, Double_t sizeX, Double_t sizeY, Double_t theta) const
 It determines the total hit energy contained inside a prisma delimited between x0 and y0 vertex, and with face dimensions sizeX by sizeY. The angle theta should serve to rotate the prism along its axis to give full freedom. More...
 
Double_t GetEnergyInSphere (const TVector3 &pos0, Double_t radius) const
 It determines the total energy contained in a sphere with position pos0 for a given spherical radius. More...
 
Double_t GetEnergyInSphere (Double_t x, Double_t y, Double_t z, Double_t radius) const
 It determines the total energy contained in a sphere with position x0,y0 and y0 for a given radius. More...
 
const std::vector< Float_t > & GetEnergyVector () const
 
Double_t GetEnergyX () const
 It calculates the total energy of hits with a valid X coordinate. More...
 
Double_t GetEnergyY () const
 It calculates the total energy of hits with a valid Y coordinate. More...
 
Double_t GetGaussSigmaX (Double_t error=150.0, Int_t nHitsMin=100000)
 It computes the gaussian sigma in the X-coordinate. It adds a hit to the right and a hit to the left, with energy = 0 +/- user defined error in ADC. Then it fits a gaussian to the hits and extracts the sigma. The hits are just added for fitting purposes and do not go into any further processing. More...
 
Double_t GetGaussSigmaY (Double_t error=150.0, Int_t nHitsMin=100000)
 It computes the gaussian sigma in the Y-coordinate. It adds a hit to the right and a hit to the left, with energy = 0 +/- user defined error in ADC. Then it fits a gaussian to the hits and extracts the sigma. The hits are just added for fitting purposes and do not go into any further processing. More...
 
Double_t GetGaussSigmaZ (Double_t error=150.0, Int_t nHitsMin=100000)
 It computes the gaussian sigma in the Z-coordinate. It adds a hit to the right and a hit to the left, with energy = 0 +/- user defined error in ADC. Then it fits a gaussian to the hits and extracts the sigma. The hits are just added for fitting purposes and do not go into any further processing. More...
 
Double_t GetHitsPathLength (Int_t n=0, Int_t m=0) const
 It determines the distance required to travel from hit n to hit m adding all the distances of the hits that are found between both. More...
 
Double_t GetHitsTwist (Int_t n, Int_t m) const
 A parameter measuring how straight is a given sequence of hits. If the value is close to zero, the hits follow a straight path in average. I believe the value should be then -1 to 1 depending where the track is twisting. Or may be just a positive value giving the measurement of twist. Not 100% sure now. More...
 
Double_t GetHitsTwistWeighted (Int_t n, Int_t m) const
 Same as GetHitsTwist but weighting with the energy. More...
 
Double_t GetMaximumHitDistance () const
 It returns the maximum distance between 2-hits. More...
 
Double_t GetMaximumHitDistance2 () const
 It returns the maximum squared distance between 2-hits. More...
 
Double_t GetMaximumHitEnergy () const
 It returns the maximum hit energy. More...
 
Double_t GetMeanHitEnergy () const
 It returns the mean hits energy. More...
 
TVector3 GetMeanPosition () const
 It calculates the mean position weighting with the energy of the hits. Each coordinate is calculated considering the valid coordinates of each hit component. More...
 
TVector3 GetMeanPositionInCylinder (const TVector3 &x0, const TVector3 &x1, Double_t radius) const
 It determines the mean position using the hits contained inside a cylinder with a given radius and delimited between x0 and x1 vertex. More...
 
TVector3 GetMeanPositionInPrism (const TVector3 &x0, const TVector3 &x1, Double_t sizeX, Double_t sizeY, Double_t theta) const
 It determines the mean position of hits contained inside a prisma delimited between x0 and x1 vertex, and with face dimensions sizeX by sizeY. The angle theta should serve to rotate the prism along its axis to give full freedom. More...
 
Double_t GetMeanPositionX () const
 It calculates the mean X position weighting with the energy of the hits with a valid X coordinate. More...
 
Double_t GetMeanPositionXInCylinder (const TVector3 &x0, const TVector3 &x1, Double_t radius) const
 It determines the mean position X using the hits contained inside a cylinder with a given radius and delimited between x0 and x1 vertex. More...
 
Double_t GetMeanPositionXInPrism (const TVector3 &x0, const TVector3 &x1, Double_t sizeX, Double_t sizeY, Double_t theta) const
 It determines the mean X position of hits contained inside a prisma delimited between x0 and x1 vertex, and with face dimensions sizeX by sizeY. The angle theta should serve to rotate the prism along its axis to give full freedom. More...
 
Double_t GetMeanPositionY () const
 It calculates the mean Y position weighting with the energy of the hits with a valid Y coordinate. More...
 
Double_t GetMeanPositionYInCylinder (const TVector3 &x0, const TVector3 &x1, Double_t radius) const
 It determines the mean position Y using the hits contained inside a cylinder with a given radius and delimited between x0 and x1 vertex. More...
 
Double_t GetMeanPositionYInPrism (const TVector3 &x0, const TVector3 &x1, Double_t sizeX, Double_t sizeY, Double_t theta) const
 It determines the mean Y position of hits contained inside a prisma delimited between x0 and x1 vertex, and with face dimensions sizeX by sizeY. The angle theta should serve to rotate the prism along its axis to give full freedom. More...
 
Double_t GetMeanPositionZ () const
 It calculates the mean Z position weighting with the energy of the hits with a valid Z coordinate. More...
 
Double_t GetMeanPositionZInCylinder (const TVector3 &x0, const TVector3 &x1, Double_t radius) const
 It determines the mean position Z using the hits contained inside a cylinder with a given radius and delimited between x0 and x1 vertex. More...
 
Double_t GetMeanPositionZInPrism (const TVector3 &x0, const TVector3 &x1, Double_t sizeX, Double_t sizeY, Double_t theta) const
 It determines the mean Z position of hits contained inside a prisma delimited between x0 and x1 vertex, and with face dimensions sizeX by sizeY. The angle theta should serve to rotate the prism along its axis to give full freedom. More...
 
Double_t GetMinimumHitEnergy () const
 It returns the minimum hit energy. More...
 
Int_t GetMostEnergeticHitInRange (Int_t n, Int_t m) const
 It returns the most energetic hit found between hits n and m. More...
 
size_t GetNumberOfHits () const
 
Int_t GetNumberOfHitsInsideCylinder (const TVector3 &x0, const TVector3 &x1, Double_t radius) const
 It determines the total number of hits contained inside a cylinder with a given radius and delimited between x0 and y0 vertex. More...
 
Int_t GetNumberOfHitsInsideCylinder (Int_t i, Int_t j, Double_t radius) const
 It determines the total number of hits contained inside a cylinder with a given radius and delimited between the hit number i and the hit number j. More...
 
Int_t GetNumberOfHitsInsidePrism (const TVector3 &x0, const TVector3 &x1, Double_t sizeX, Double_t sizeY, Double_t theta) const
 It determines the total number of hits contained inside a prisma delimited between x0 and y0 vertex, and with face dimensions sizeX by sizeY. The angle theta should serve to rotate the prism along its axis to give full freedom. More...
 
Int_t GetNumberOfHitsX () const
 It returns the number of hits with a valid X coordinate. More...
 
Int_t GetNumberOfHitsY () const
 It returns the number of hits with a valid Y coordinate. More...
 
TVector3 GetPosition (int n) const
 It returns the position of hit number n. More...
 
TVector2 GetProjection (Int_t n, Int_t m, const TVector3 &position) const
 It returns the longitudinal and transversal projections of position to the axis defined by the hits n and m. More...
 
Double_t GetSigmaX () const
 It calculates the hits standard deviation in the X-coordinate. More...
 
Double_t GetSigmaXY2 () const
 It calculates the 2-dimensional hits variance. More...
 
Double_t GetSigmaY () const
 It calculates the hits standard deviation in the Y-coordinate. More...
 
Double_t GetSigmaZ2 () const
 It returns the hits distribution variance on the Z-axis. More...
 
Double_t GetSkewXY () const
 It returns the 2-dimensional skewness on the XY-plane which is a measure of the hits distribution asymmetry. More...
 
Double_t GetSkewZ () const
 It returns the hits distribution skewness, or asymmetry on the Z-axis. More...
 
const std::vector< Float_t > & GetTime () const
 
Double_t GetTime (int n) const
 
Double_t GetTotalDistance () const
 It determines the distance required to travel from the first to the last hit adding all the distances of the hits that are found between both. More...
 
Double_t GetTotalEnergy () const
 
Double_t GetTransversalProjection (const TVector3 &p0, const TVector3 &direction, const TVector3 &position) const
 It returns the transversal projection of position to the line defined by position and direction. More...
 
REST_HitType GetType (int n) const
 
TVector3 GetVector (int i, int j) const
 It returns the vector that goes from hit j to hit i. More...
 
const std::vector< Float_t > & GetX () const
 
Double_t GetX (int n) const
 
const std::vector< Float_t > & GetY () const
 
Double_t GetY (int n) const
 
const std::vector< Float_t > & GetZ () const
 
Double_t GetZ (int n) const
 
Bool_t isHitNInsideCylinder (Int_t n, const TVector3 &x0, const TVector3 &x1, Double_t radius) const
 It determines if hit n is contained inside a cylinder with a given radius and delimited between x0 and x1 vertex. More...
 
Bool_t isHitNInsidePrism (Int_t n, const TVector3 &x0, const TVector3 &x1, Double_t sizeX, Double_t sizeY, Double_t theta) const
 It determines if hit n is contained inside a prisma delimited between x0 and y0 vertex, and with face dimensions sizeX by sizeY. The angle theta should serve to rotate the prism along its axis to give full freedom. More...
 
Bool_t isHitNInsideSphere (Int_t n, const TVector3 &pos0, Double_t radius) const
 It determines if the hit n is contained in a sphere with position pos0 for a given sphereical radius. More...
 
Bool_t isHitNInsideSphere (Int_t n, Double_t x0, Double_t y0, Double_t z0, Double_t radius) const
 It determines the total energy contained in a sphere with position x0,y0 and y0 for a given radius. More...
 
Bool_t isNaN (Int_t n) const
 It will return true only if all the 3-coordinates of hit number n are not a number,. More...
 
Bool_t isSortedByEnergy () const
 It returns true if the hits are ordered in increasing energies. More...
 
virtual void MergeHits (int n, int m)
 It merges hits n and m being the resulting hit placed at the weighted center and being its final energy the addition of the energies of the hits n and m. More...
 
virtual void PrintHits (Int_t nHits=-1) const
 It prints on screen the first nHits from the list. More...
 
virtual void RemoveHit (int n)
 It removes the hit at position n from the list. More...
 
virtual void RemoveHits ()
 It removes all hits inside the class. More...
 
void Rotate (Int_t n, Double_t alpha, const TVector3 &vAxis, const TVector3 &vMean)
 It rotates hit n by an angle akpha along the vAxis with center at vMean. More...
 
void RotateIn3D (Int_t n, Double_t alpha, Double_t beta, Double_t gamma, const TVector3 &center)
 It rotates hit n following rotations in Z, Y and X by angles gamma, beta and alpha. The rotation is performed with center at vMean. More...
 
virtual void SwapHits (Int_t i, Int_t j)
 It exchanges hits n and m affecting to the ordering of the hits inside the list of hits. More...
 
void Translate (Int_t n, Double_t x, Double_t y, Double_t z)
 It moves hit n by a given amount (x,y,z). More...
 
 TRestHits ()
 Default constructor.
 
void WriteHitsToTextFile (TString filename)
 It writes the hits to a plain text file. More...
 
 ~TRestHits ()
 Default destructor.
 

Static Public Member Functions

static void GetBoundaries (std::vector< double > &val, double &max, double &min, int &nBins, double offset=10)
 TODO This method is not using any TRestHits member. This probably means that it should be placed somewhere else. More...
 

Protected Attributes

std::vector< Float_t > fEnergy
 Energy deposited at each 3-coordinate position (units keV) More...
 
size_t fNHits = 0
 Number of punctual energy depositions, it is the length for all the arrays. More...
 
std::vector< Float_t > fTime
 Absolute time information for each punctual deposition (units us, 0 is time of decay) More...
 
Double_t fTotalEnergy = 0
 Event total energy. More...
 
std::vector< REST_HitType > fType
 The type of hit X,Y,XY,XYZ, ... More...
 
std::vector< Float_t > fX
 Position on X axis for each punctual deposition (units mm) More...
 
std::vector< Float_t > fY
 Position on Y axis for each punctual deposition (units mm) More...
 
std::vector< Float_t > fZ
 Position on Z axis for each punctual deposition (units mm) More...
 

Member Typedef Documentation

◆ iterator

Definition at line 275 of file TRestHits.h.

Member Function Documentation

◆ AddHit() [1/2]

void TRestHits::AddHit ( const TVector3 &  position,
Double_t  energy,
Double_t  time = 0,
REST_HitType  type = XYZ 
)

Adds a new hit to the list of hits using a TVector3.

Definition at line 345 of file TRestHits.cxx.

◆ AddHit() [2/2]

void TRestHits::AddHit ( TRestHits hits,
Int_t  n 
)

Adds a new hit to the list of hits using the hit n inside another TRestHits object.

Definition at line 357 of file TRestHits.cxx.

◆ areXY()

Bool_t TRestHits::areXY ( ) const
virtual

It will return true only if all the hits inside are of type XY.

Reimplemented in TRestVolumeHits.

Definition at line 101 of file TRestHits.cxx.

◆ areXYZ()

Bool_t TRestHits::areXYZ ( ) const
virtual

It will return true only if all the hits inside are of type XYZ.

Reimplemented in TRestVolumeHits.

Definition at line 146 of file TRestHits.cxx.

◆ areXZ()

Bool_t TRestHits::areXZ ( ) const
virtual

It will return true only if all the hits inside are of type XZ.

Reimplemented in TRestVolumeHits.

Definition at line 116 of file TRestHits.cxx.

◆ areYZ()

Bool_t TRestHits::areYZ ( ) const
virtual

It will return true only if all the hits inside are of type YZ.

Reimplemented in TRestVolumeHits.

Definition at line 131 of file TRestHits.cxx.

◆ begin()

TRestHits_Iterator TRestHits::begin ( )
inline

Definition at line 273 of file TRestHits.h.

◆ end()

TRestHits_Iterator TRestHits::end ( )
inline

Definition at line 274 of file TRestHits.h.

◆ GetBoundaries()

void TRestHits::GetBoundaries ( std::vector< double > &  val,
double &  max,
double &  min,
int &  nBins,
double  offset = 10 
)
static

TODO This method is not using any TRestHits member. This probably means that it should be placed somewhere else.

Definition at line 738 of file TRestHits.cxx.

◆ GetClosestHit()

Int_t TRestHits::GetClosestHit ( const TVector3 &  position) const

It returns the closest hit to a given position.

Definition at line 1244 of file TRestHits.cxx.

◆ GetDistance()

Double_t TRestHits::GetDistance ( int  N,
int  M 
) const
inline

Definition at line 173 of file TRestHits.h.

◆ GetDistance2()

Double_t TRestHits::GetDistance2 ( int  n,
int  m 
) const

It returns the euclidian distance between hits n and m.

Definition at line 1201 of file TRestHits.cxx.

◆ GetDistanceToNode()

Double_t TRestHits::GetDistanceToNode ( Int_t  n) const

It determines the distance required to travel from the first hit to the hit n adding all the distances of the hits that are found till the hit n.

Definition at line 1217 of file TRestHits.cxx.

◆ GetEnergy()

Double_t TRestHits::GetEnergy ( int  n) const
inline

Definition at line 105 of file TRestHits.h.

◆ GetEnergyInCylinder() [1/2]

Double_t TRestHits::GetEnergyInCylinder ( const TVector3 &  x0,
const TVector3 &  x1,
Double_t  radius 
) const

It determines the total energy contained inside a cylinder with a given radius and delimited between x0 and y0 vertex.

Definition at line 262 of file TRestHits.cxx.

◆ GetEnergyInCylinder() [2/2]

Double_t TRestHits::GetEnergyInCylinder ( Int_t  i,
Int_t  j,
Double_t  radius 
) const

It determines the total energy contained inside a cylinder with a given radius and delimited between the hit number i and the hit number j.

Definition at line 254 of file TRestHits.cxx.

◆ GetEnergyInPrism()

Double_t TRestHits::GetEnergyInPrism ( const TVector3 &  x0,
const TVector3 &  x1,
Double_t  sizeX,
Double_t  sizeY,
Double_t  theta 
) const

It determines the total hit energy contained inside a prisma delimited between x0 and y0 vertex, and with face dimensions sizeX by sizeY. The angle theta should serve to rotate the prism along its axis to give full freedom.

Definition at line 200 of file TRestHits.cxx.

◆ GetEnergyInSphere() [1/2]

Double_t TRestHits::GetEnergyInSphere ( const TVector3 &  pos0,
Double_t  radius 
) const

It determines the total energy contained in a sphere with position pos0 for a given spherical radius.

Definition at line 296 of file TRestHits.cxx.

◆ GetEnergyInSphere() [2/2]

Double_t TRestHits::GetEnergyInSphere ( Double_t  x,
Double_t  y,
Double_t  z,
Double_t  radius 
) const

It determines the total energy contained in a sphere with position x0,y0 and y0 for a given radius.

Definition at line 304 of file TRestHits.cxx.

◆ GetEnergyVector()

const std::vector< Float_t > & TRestHits::GetEnergyVector ( ) const
inline

Definition at line 99 of file TRestHits.h.

◆ GetEnergyX()

Double_t TRestHits::GetEnergyX ( ) const

It calculates the total energy of hits with a valid X coordinate.

Definition at line 564 of file TRestHits.cxx.

◆ GetEnergyY()

Double_t TRestHits::GetEnergyY ( ) const

It calculates the total energy of hits with a valid Y coordinate.

Definition at line 578 of file TRestHits.cxx.

◆ GetGaussSigmaX()

Double_t TRestHits::GetGaussSigmaX ( Double_t  error = 150.0,
Int_t  nHitsMin = 100000 
)

It computes the gaussian sigma in the X-coordinate. It adds a hit to the right and a hit to the left, with energy = 0 +/- user defined error in ADC. Then it fits a gaussian to the hits and extracts the sigma. The hits are just added for fitting purposes and do not go into any further processing.

Definition at line 762 of file TRestHits.cxx.

◆ GetGaussSigmaY()

Double_t TRestHits::GetGaussSigmaY ( Double_t  error = 150.0,
Int_t  nHitsMin = 100000 
)

It computes the gaussian sigma in the Y-coordinate. It adds a hit to the right and a hit to the left, with energy = 0 +/- user defined error in ADC. Then it fits a gaussian to the hits and extracts the sigma. The hits are just added for fitting purposes and do not go into any further processing.

Definition at line 836 of file TRestHits.cxx.

◆ GetGaussSigmaZ()

Double_t TRestHits::GetGaussSigmaZ ( Double_t  error = 150.0,
Int_t  nHitsMin = 100000 
)

It computes the gaussian sigma in the Z-coordinate. It adds a hit to the right and a hit to the left, with energy = 0 +/- user defined error in ADC. Then it fits a gaussian to the hits and extracts the sigma. The hits are just added for fitting purposes and do not go into any further processing.

Definition at line 907 of file TRestHits.cxx.

◆ GetHitsPathLength()

Double_t TRestHits::GetHitsPathLength ( Int_t  n = 0,
Int_t  m = 0 
) const

It determines the distance required to travel from hit n to hit m adding all the distances of the hits that are found between both.

Definition at line 1179 of file TRestHits.cxx.

◆ GetHitsTwist()

Double_t TRestHits::GetHitsTwist ( Int_t  n,
Int_t  m 
) const

A parameter measuring how straight is a given sequence of hits. If the value is close to zero, the hits follow a straight path in average. I believe the value should be then -1 to 1 depending where the track is twisting. Or may be just a positive value giving the measurement of twist. Not 100% sure now.

Definition at line 1301 of file TRestHits.cxx.

◆ GetHitsTwistWeighted()

Double_t TRestHits::GetHitsTwistWeighted ( Int_t  n,
Int_t  m 
) const

Same as GetHitsTwist but weighting with the energy.

Definition at line 1325 of file TRestHits.cxx.

◆ GetMaximumHitDistance()

Double_t TRestHits::GetMaximumHitDistance ( ) const

It returns the maximum distance between 2-hits.

Definition at line 1352 of file TRestHits.cxx.

◆ GetMaximumHitDistance2()

Double_t TRestHits::GetMaximumHitDistance2 ( ) const

It returns the maximum squared distance between 2-hits.

Definition at line 1366 of file TRestHits.cxx.

◆ GetMaximumHitEnergy()

Double_t TRestHits::GetMaximumHitEnergy ( ) const

It returns the maximum hit energy.

Definition at line 426 of file TRestHits.cxx.

◆ GetMeanHitEnergy()

Double_t TRestHits::GetMeanHitEnergy ( ) const

It returns the mean hits energy.

Definition at line 452 of file TRestHits.cxx.

◆ GetMeanPosition()

TVector3 TRestHits::GetMeanPosition ( ) const

It calculates the mean position weighting with the energy of the hits. Each coordinate is calculated considering the valid coordinates of each hit component.

Definition at line 658 of file TRestHits.cxx.

◆ GetMeanPositionInCylinder()

TVector3 TRestHits::GetMeanPositionInCylinder ( const TVector3 &  x0,
const TVector3 &  x1,
Double_t  radius 
) const

It determines the mean position using the hits contained inside a cylinder with a given radius and delimited between x0 and x1 vertex.

Definition at line 1169 of file TRestHits.cxx.

◆ GetMeanPositionInPrism()

TVector3 TRestHits::GetMeanPositionInPrism ( const TVector3 &  x0,
const TVector3 &  x1,
Double_t  sizeX,
Double_t  sizeY,
Double_t  theta 
) const

It determines the mean position of hits contained inside a prisma delimited between x0 and x1 vertex, and with face dimensions sizeX by sizeY. The angle theta should serve to rotate the prism along its axis to give full freedom.

Definition at line 1095 of file TRestHits.cxx.

◆ GetMeanPositionX()

Double_t TRestHits::GetMeanPositionX ( ) const

It calculates the mean X position weighting with the energy of the hits with a valid X coordinate.

Definition at line 593 of file TRestHits.cxx.

◆ GetMeanPositionXInCylinder()

Double_t TRestHits::GetMeanPositionXInCylinder ( const TVector3 &  x0,
const TVector3 &  x1,
Double_t  radius 
) const

It determines the mean position X using the hits contained inside a cylinder with a given radius and delimited between x0 and x1 vertex.

Definition at line 1107 of file TRestHits.cxx.

◆ GetMeanPositionXInPrism()

Double_t TRestHits::GetMeanPositionXInPrism ( const TVector3 &  x0,
const TVector3 &  x1,
Double_t  sizeX,
Double_t  sizeY,
Double_t  theta 
) const

It determines the mean X position of hits contained inside a prisma delimited between x0 and x1 vertex, and with face dimensions sizeX by sizeY. The angle theta should serve to rotate the prism along its axis to give full freedom.

Definition at line 1030 of file TRestHits.cxx.

◆ GetMeanPositionY()

Double_t TRestHits::GetMeanPositionY ( ) const

It calculates the mean Y position weighting with the energy of the hits with a valid Y coordinate.

Definition at line 615 of file TRestHits.cxx.

◆ GetMeanPositionYInCylinder()

Double_t TRestHits::GetMeanPositionYInCylinder ( const TVector3 &  x0,
const TVector3 &  x1,
Double_t  radius 
) const

It determines the mean position Y using the hits contained inside a cylinder with a given radius and delimited between x0 and x1 vertex.

Definition at line 1128 of file TRestHits.cxx.

◆ GetMeanPositionYInPrism()

Double_t TRestHits::GetMeanPositionYInPrism ( const TVector3 &  x0,
const TVector3 &  x1,
Double_t  sizeX,
Double_t  sizeY,
Double_t  theta 
) const

It determines the mean Y position of hits contained inside a prisma delimited between x0 and x1 vertex, and with face dimensions sizeX by sizeY. The angle theta should serve to rotate the prism along its axis to give full freedom.

Definition at line 1052 of file TRestHits.cxx.

◆ GetMeanPositionZ()

Double_t TRestHits::GetMeanPositionZ ( ) const

It calculates the mean Z position weighting with the energy of the hits with a valid Z coordinate.

Definition at line 637 of file TRestHits.cxx.

◆ GetMeanPositionZInCylinder()

Double_t TRestHits::GetMeanPositionZInCylinder ( const TVector3 &  x0,
const TVector3 &  x1,
Double_t  radius 
) const

It determines the mean position Z using the hits contained inside a cylinder with a given radius and delimited between x0 and x1 vertex.

Definition at line 1149 of file TRestHits.cxx.

◆ GetMeanPositionZInPrism()

Double_t TRestHits::GetMeanPositionZInPrism ( const TVector3 &  x0,
const TVector3 &  x1,
Double_t  sizeX,
Double_t  sizeY,
Double_t  theta 
) const

It determines the mean Z position of hits contained inside a prisma delimited between x0 and x1 vertex, and with face dimensions sizeX by sizeY. The angle theta should serve to rotate the prism along its axis to give full freedom.

Definition at line 1074 of file TRestHits.cxx.

◆ GetMinimumHitEnergy()

Double_t TRestHits::GetMinimumHitEnergy ( ) const

It returns the minimum hit energy.

Definition at line 439 of file TRestHits.cxx.

◆ GetMostEnergeticHitInRange()

Int_t TRestHits::GetMostEnergeticHitInRange ( Int_t  n,
Int_t  m 
) const

It returns the most energetic hit found between hits n and m.

Definition at line 1229 of file TRestHits.cxx.

◆ GetNumberOfHits()

size_t TRestHits::GetNumberOfHits ( ) const
inline

Definition at line 93 of file TRestHits.h.

◆ GetNumberOfHitsInsideCylinder() [1/2]

Int_t TRestHits::GetNumberOfHitsInsideCylinder ( const TVector3 &  x0,
const TVector3 &  x1,
Double_t  radius 
) const

It determines the total number of hits contained inside a cylinder with a given radius and delimited between x0 and y0 vertex.

Definition at line 283 of file TRestHits.cxx.

◆ GetNumberOfHitsInsideCylinder() [2/2]

Int_t TRestHits::GetNumberOfHitsInsideCylinder ( Int_t  i,
Int_t  j,
Double_t  radius 
) const

It determines the total number of hits contained inside a cylinder with a given radius and delimited between the hit number i and the hit number j.

Definition at line 275 of file TRestHits.cxx.

◆ GetNumberOfHitsInsidePrism()

Int_t TRestHits::GetNumberOfHitsInsidePrism ( const TVector3 &  x0,
const TVector3 &  x1,
Double_t  sizeX,
Double_t  sizeY,
Double_t  theta 
) const

It determines the total number of hits contained inside a prisma delimited between x0 and y0 vertex, and with face dimensions sizeX by sizeY. The angle theta should serve to rotate the prism along its axis to give full freedom.

Definition at line 216 of file TRestHits.cxx.

◆ GetNumberOfHitsX()

Int_t TRestHits::GetNumberOfHitsX ( ) const

It returns the number of hits with a valid X coordinate.

Definition at line 536 of file TRestHits.cxx.

◆ GetNumberOfHitsY()

Int_t TRestHits::GetNumberOfHitsY ( ) const

It returns the number of hits with a valid Y coordinate.

Definition at line 550 of file TRestHits.cxx.

◆ GetPosition()

TVector3 TRestHits::GetPosition ( int  n) const

It returns the position of hit number n.

Definition at line 515 of file TRestHits.cxx.

◆ GetProjection()

TVector2 TRestHits::GetProjection ( Int_t  n,
Int_t  m,
const TVector3 &  position 
) const

It returns the longitudinal and transversal projections of position to the axis defined by the hits n and m.

Definition at line 1265 of file TRestHits.cxx.

◆ GetSigmaX()

Double_t TRestHits::GetSigmaX ( ) const

It calculates the hits standard deviation in the X-coordinate.

Definition at line 685 of file TRestHits.cxx.

◆ GetSigmaXY2()

Double_t TRestHits::GetSigmaXY2 ( ) const

It calculates the 2-dimensional hits variance.

Definition at line 666 of file TRestHits.cxx.

◆ GetSigmaY()

Double_t TRestHits::GetSigmaY ( ) const

It calculates the hits standard deviation in the Y-coordinate.

Definition at line 703 of file TRestHits.cxx.

◆ GetSigmaZ2()

Double_t TRestHits::GetSigmaZ2 ( ) const

It returns the hits distribution variance on the Z-axis.

Definition at line 997 of file TRestHits.cxx.

◆ GetSkewXY()

Double_t TRestHits::GetSkewXY ( ) const

It returns the 2-dimensional skewness on the XY-plane which is a measure of the hits distribution asymmetry.

Definition at line 979 of file TRestHits.cxx.

◆ GetSkewZ()

Double_t TRestHits::GetSkewZ ( ) const

It returns the hits distribution skewness, or asymmetry on the Z-axis.

Definition at line 1013 of file TRestHits.cxx.

◆ GetTime() [1/2]

const std::vector< Float_t > & TRestHits::GetTime ( ) const
inline

Definition at line 98 of file TRestHits.h.

◆ GetTime() [2/2]

Double_t TRestHits::GetTime ( int  n) const
inline

Definition at line 104 of file TRestHits.h.

◆ GetTotalDistance()

Double_t TRestHits::GetTotalDistance ( ) const

It determines the distance required to travel from the first to the last hit adding all the distances of the hits that are found between both.

Definition at line 1192 of file TRestHits.cxx.

◆ GetTotalEnergy()

Double_t TRestHits::GetTotalEnergy ( ) const

Definition at line 1393 of file TRestHits.cxx.

◆ GetTransversalProjection()

Double_t TRestHits::GetTransversalProjection ( const TVector3 &  p0,
const TVector3 &  direction,
const TVector3 &  position 
) const

It returns the transversal projection of position to the line defined by position and direction.

Definition at line 1284 of file TRestHits.cxx.

◆ GetType()

REST_HitType TRestHits::GetType ( int  n) const
inline

Definition at line 107 of file TRestHits.h.

◆ GetVector()

TVector3 TRestHits::GetVector ( int  i,
int  j 
) const

It returns the vector that goes from hit j to hit i.

Definition at line 531 of file TRestHits.cxx.

◆ GetX() [1/2]

const std::vector< Float_t > & TRestHits::GetX ( ) const
inline

Definition at line 95 of file TRestHits.h.

◆ GetX() [2/2]

Double_t TRestHits::GetX ( int  n) const
inline

Definition at line 101 of file TRestHits.h.

◆ GetY() [1/2]

const std::vector< Float_t > & TRestHits::GetY ( ) const
inline

Definition at line 96 of file TRestHits.h.

◆ GetY() [2/2]

Double_t TRestHits::GetY ( int  n) const
inline

Definition at line 102 of file TRestHits.h.

◆ GetZ() [1/2]

const std::vector< Float_t > & TRestHits::GetZ ( ) const
inline

Definition at line 97 of file TRestHits.h.

◆ GetZ() [2/2]

Double_t TRestHits::GetZ ( int  n) const
inline

Definition at line 103 of file TRestHits.h.

◆ isHitNInsideCylinder()

Bool_t TRestHits::isHitNInsideCylinder ( Int_t  n,
const TVector3 &  x0,
const TVector3 &  x1,
Double_t  radius 
) const

It determines if hit n is contained inside a cylinder with a given radius and delimited between x0 and x1 vertex.

Definition at line 230 of file TRestHits.cxx.

◆ isHitNInsidePrism()

Bool_t TRestHits::isHitNInsidePrism ( Int_t  n,
const TVector3 &  x0,
const TVector3 &  x1,
Double_t  sizeX,
Double_t  sizeY,
Double_t  theta 
) const

It determines if hit n is contained inside a prisma delimited between x0 and y0 vertex, and with face dimensions sizeX by sizeY. The angle theta should serve to rotate the prism along its axis to give full freedom.

TODO: It seems to me there is a problem with the rotation of the hits, which are rotated along Z axis, and not along the prism axis = x1-x0. As soon as the prism is aligned with Z no problem though.

Definition at line 176 of file TRestHits.cxx.

◆ isHitNInsideSphere() [1/2]

Bool_t TRestHits::isHitNInsideSphere ( Int_t  n,
const TVector3 &  pos0,
Double_t  radius 
) const

It determines if the hit n is contained in a sphere with position pos0 for a given sphereical radius.

Definition at line 322 of file TRestHits.cxx.

◆ isHitNInsideSphere() [2/2]

Bool_t TRestHits::isHitNInsideSphere ( Int_t  n,
Double_t  x0,
Double_t  y0,
Double_t  z0,
Double_t  radius 
) const

It determines the total energy contained in a sphere with position x0,y0 and y0 for a given radius.

Definition at line 330 of file TRestHits.cxx.

◆ isNaN()

Bool_t TRestHits::isNaN ( Int_t  n) const

It will return true only if all the 3-coordinates of hit number n are not a number,.

Definition at line 162 of file TRestHits.cxx.

◆ isSortedByEnergy()

Bool_t TRestHits::isSortedByEnergy ( ) const

It returns true if the hits are ordered in increasing energies.

Definition at line 490 of file TRestHits.cxx.

◆ MergeHits()

void TRestHits::MergeHits ( int  n,
int  m 
)
virtual

It merges hits n and m being the resulting hit placed at the weighted center and being its final energy the addition of the energies of the hits n and m.

Definition at line 458 of file TRestHits.cxx.

◆ PrintHits()

void TRestHits::PrintHits ( Int_t  nHits = -1) const
virtual

It prints on screen the first nHits from the list.

Definition at line 1380 of file TRestHits.cxx.

◆ RemoveHit()

void TRestHits::RemoveHit ( int  n)
virtual

It removes the hit at position n from the list.

Reimplemented in TRestVolumeHits.

Definition at line 503 of file TRestHits.cxx.

◆ RemoveHits()

void TRestHits::RemoveHits ( )
virtual

It removes all hits inside the class.

Reimplemented in TRestVolumeHits.

Definition at line 371 of file TRestHits.cxx.

◆ Rotate()

void TRestHits::Rotate ( Int_t  n,
Double_t  alpha,
const TVector3 &  vAxis,
const TVector3 &  vMean 
)

It rotates hit n by an angle akpha along the vAxis with center at vMean.

Definition at line 409 of file TRestHits.cxx.

◆ RotateIn3D()

void TRestHits::RotateIn3D ( Int_t  n,
Double_t  alpha,
Double_t  beta,
Double_t  gamma,
const TVector3 &  center 
)

It rotates hit n following rotations in Z, Y and X by angles gamma, beta and alpha. The rotation is performed with center at vMean.

Definition at line 393 of file TRestHits.cxx.

◆ SwapHits()

void TRestHits::SwapHits ( Int_t  i,
Int_t  j 
)
virtual

It exchanges hits n and m affecting to the ordering of the hits inside the list of hits.

Reimplemented in TRestVolumeHits.

Definition at line 478 of file TRestHits.cxx.

◆ Translate()

void TRestHits::Translate ( Int_t  n,
Double_t  x,
Double_t  y,
Double_t  z 
)

It moves hit n by a given amount (x,y,z).

Definition at line 383 of file TRestHits.cxx.

◆ WriteHitsToTextFile()

void TRestHits::WriteHitsToTextFile ( TString  filename)

It writes the hits to a plain text file.

Definition at line 721 of file TRestHits.cxx.

Field Documentation

◆ fEnergy

std::vector<Float_t> TRestHits::fEnergy
protected

Energy deposited at each 3-coordinate position (units keV)

Definition at line 59 of file TRestHits.h.

◆ fNHits

size_t TRestHits::fNHits = 0
protected

Number of punctual energy depositions, it is the length for all the arrays.

Definition at line 42 of file TRestHits.h.

◆ fTime

std::vector<Float_t> TRestHits::fTime
protected

Absolute time information for each punctual deposition (units us, 0 is time of decay)

Definition at line 56 of file TRestHits.h.

◆ fTotalEnergy

Double_t TRestHits::fTotalEnergy = 0
protected

Event total energy.

Definition at line 44 of file TRestHits.h.

◆ fType

std::vector<REST_HitType> TRestHits::fType
protected

The type of hit X,Y,XY,XYZ, ...

Definition at line 62 of file TRestHits.h.

◆ fX

std::vector<Float_t> TRestHits::fX
protected

Position on X axis for each punctual deposition (units mm)

Definition at line 47 of file TRestHits.h.

◆ fY

std::vector<Float_t> TRestHits::fY
protected

Position on Y axis for each punctual deposition (units mm)

Definition at line 50 of file TRestHits.h.

◆ fZ

std::vector<Float_t> TRestHits::fZ
protected

Position on Z axis for each punctual deposition (units mm)

Definition at line 53 of file TRestHits.h.


The documentation for this class was generated from the following files: