REST-for-Physics
v2.3
Rare Event Searches ToolKit for Physics
|
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.
2015-Nov: Changed vectors fX fY fZ and fEnergy from <Int_t> to < Float_t>.
2022-July: Introducing gausian hits fitting
Definition at line 39 of file TRestHits.h.
#include <TRestHits.h>
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. | |
void | AddHit (TRestHits &hits, Int_t n) |
Adds a new hit to the list of hits using the hit n inside another TRestHits object. | |
virtual Bool_t | areXY () const |
It will return true only if all the hits inside are of type XY. | |
virtual Bool_t | areXYZ () const |
It will return true only if all the hits inside are of type XYZ. | |
virtual Bool_t | areXZ () const |
It will return true only if all the hits inside are of type XZ. | |
virtual Bool_t | areYZ () const |
It will return true only if all the hits inside are of type YZ. | |
TRestHits_Iterator | begin () |
TRestHits_Iterator | end () |
Int_t | GetClosestHit (const TVector3 &position) const |
It returns the closest hit to a given position . | |
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 . | |
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 . | |
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. | |
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 . | |
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. | |
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 . | |
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 . | |
const std::vector< Float_t > & | GetEnergyVector () const |
Double_t | GetEnergyX () const |
It calculates the total energy of hits with a valid X coordinate. | |
Double_t | GetEnergyY () const |
It calculates the total energy of hits with a valid Y coordinate. | |
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. | |
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. | |
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. | |
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. | |
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. | |
Double_t | GetHitsTwistWeighted (Int_t n, Int_t m) const |
Same as GetHitsTwist but weighting with the energy. | |
Double_t | GetMaximumHitDistance () const |
It returns the maximum distance between 2-hits. | |
Double_t | GetMaximumHitDistance2 () const |
It returns the maximum squared distance between 2-hits. | |
Double_t | GetMaximumHitEnergy () const |
It returns the maximum hit energy. | |
Double_t | GetMeanHitEnergy () const |
It returns the mean hits energy. | |
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. | |
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. | |
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. | |
Double_t | GetMeanPositionX () const |
It calculates the mean X position weighting with the energy of the hits with a valid X coordinate. | |
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. | |
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. | |
Double_t | GetMeanPositionY () const |
It calculates the mean Y position weighting with the energy of the hits with a valid Y coordinate. | |
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. | |
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. | |
Double_t | GetMeanPositionZ () const |
It calculates the mean Z position weighting with the energy of the hits with a valid Z coordinate. | |
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. | |
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. | |
Double_t | GetMinimumHitEnergy () const |
It returns the minimum hit energy. | |
Int_t | GetMostEnergeticHitInRange (Int_t n, Int_t m) const |
It returns the most energetic hit found between hits n and m . | |
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. | |
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 . | |
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. | |
Int_t | GetNumberOfHitsX () const |
It returns the number of hits with a valid X coordinate. | |
Int_t | GetNumberOfHitsY () const |
It returns the number of hits with a valid Y coordinate. | |
TVector3 | GetPosition (int n) const |
It returns the position of hit number n . | |
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 . | |
Double_t | GetSigmaX () const |
It calculates the hits standard deviation in the X-coordinate. | |
Double_t | GetSigmaXY2 () const |
It calculates the 2-dimensional hits variance. | |
Double_t | GetSigmaY () const |
It calculates the hits standard deviation in the Y-coordinate. | |
Double_t | GetSigmaZ2 () const |
It returns the hits distribution variance on the Z-axis. | |
Double_t | GetSkewXY () const |
It returns the 2-dimensional skewness on the XY-plane which is a measure of the hits distribution asymmetry. | |
Double_t | GetSkewZ () const |
It returns the hits distribution skewness, or asymmetry on the Z-axis. | |
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. | |
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 . | |
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 . | |
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. | |
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 . | |
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 . | |
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,. | |
Bool_t | isSortedByEnergy () const |
It returns true if the hits are ordered in increasing energies. | |
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 . | |
virtual void | PrintHits (Int_t nHits=-1) const |
It prints on screen the first nHits from the list. | |
virtual void | RemoveHit (int n) |
It removes the hit at position n from the list. | |
virtual void | RemoveHits () |
It removes all hits inside the class. | |
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 . | |
void | RotateIn3D (Int_t n, Double_t alpha, Double_t beta, Double_t gamma, const TVector3 ¢er) |
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 . | |
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. | |
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). | |
TRestHits () | |
Default constructor. | |
void | WriteHitsToTextFile (TString filename) |
It writes the hits to a plain text file. | |
~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. | |
Protected Attributes | |
std::vector< Float_t > | fEnergy |
Energy deposited at each 3-coordinate position (units keV) | |
size_t | fNHits = 0 |
Number of punctual energy depositions, it is the length for all the arrays. | |
std::vector< Float_t > | fTime |
Absolute time information for each punctual deposition (units us, 0 is time of decay) | |
Double_t | fTotalEnergy = 0 |
Event total energy. | |
std::vector< REST_HitType > | fType |
The type of hit X,Y,XY,XYZ, ... | |
std::vector< Float_t > | fX |
Position on X axis for each punctual deposition (units mm) | |
std::vector< Float_t > | fY |
Position on Y axis for each punctual deposition (units mm) | |
std::vector< Float_t > | fZ |
Position on Z axis for each punctual deposition (units mm) | |
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.